Line data Source code
1 : ! This file is part of ELPA.
2 : !
3 : ! The ELPA library was originally created by the ELPA consortium,
4 : ! consisting of the following organizations:
5 : !
6 : ! - Max Planck Computing and Data Facility (MPCDF), formerly known as
7 : ! Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
8 : ! - Bergische Universität Wuppertal, Lehrstuhl für angewandte
9 : ! Informatik,
10 : ! - Technische Universität München, Lehrstuhl für Informatik mit
11 : ! Schwerpunkt Wissenschaftliches Rechnen ,
12 : ! - Fritz-Haber-Institut, Berlin, Abt. Theorie,
13 : ! - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
14 : ! Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
15 : ! and
16 : ! - IBM Deutschland GmbH
17 : !
18 : !
19 : ! More information can be found here:
20 : ! http://elpa.mpcdf.mpg.de/
21 : !
22 : ! ELPA is free software: you can redistribute it and/or modify
23 : ! it under the terms of the version 3 of the license of the
24 : ! GNU Lesser General Public License as published by the Free
25 : ! Software Foundation.
26 : !
27 : ! ELPA is distributed in the hope that it will be useful,
28 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : ! GNU Lesser General Public License for more details.
31 : !
32 : ! You should have received a copy of the GNU Lesser General Public License
33 : ! along with ELPA. If not, see <http://www.gnu.org/licenses/>
34 : !
35 : ! ELPA reflects a substantial effort on the part of the original
36 : ! ELPA consortium, and we ask you to respect the spirit of the
37 : ! license that we chose: i.e., please contribute any changes you
38 : ! may have back to the original ELPA library distribution, and keep
39 : ! any derivatives of ELPA under the same license that we chose for
40 : ! the original distribution, the GNU Lesser General Public License.
41 : !
42 : !
43 : #include "config-f90.h"
44 : #include "../assert.h"
45 : !>
46 : !> Fortran test programm to demonstrates the use of
47 : !> ELPA 2 real case library.
48 : !> If "HAVE_REDIRECT" was defined at build time
49 : !> the stdout and stderr output of each MPI task
50 : !> can be redirected to files if the environment
51 : !> variable "REDIRECT_ELPA_TEST_OUTPUT" is set
52 : !> to "true".
53 : !>
54 : !> By calling executable [arg1] [arg2] [arg3] [arg4]
55 : !> one can define the size (arg1), the number of
56 : !> Eigenvectors to compute (arg2), and the blocking (arg3).
57 : !> If these values are not set default values (500, 150, 16)
58 : !> are choosen.
59 : !> If these values are set the 4th argument can be
60 : !> "output", which specifies that the EV's are written to
61 : !> an ascii file.
62 : !>
63 : !> The real ELPA 2 kernel is set as the default kernel.
64 : !> However, this can be overriden by setting
65 : !> the environment variable "REAL_ELPA_KERNEL" to an
66 : !> appropiate value.
67 : !>
68 192 : program test_real2_double_banded
69 :
70 : !-------------------------------------------------------------------------------
71 : ! Standard eigenvalue problem - REAL version
72 : !
73 : ! This program demonstrates the use of the ELPA module
74 : ! together with standard scalapack routines
75 : !
76 : ! Copyright of the original code rests with the authors inside the ELPA
77 : ! consortium. The copyright of any additional modifications shall rest
78 : ! with their original authors, but shall adhere to the licensing terms
79 : ! distributed along with the original code in the file "COPYING".
80 : !
81 : !-------------------------------------------------------------------------------
82 192 : use elpa
83 :
84 : use test_util
85 : use test_read_input_parameters
86 : use test_check_correctness
87 : use test_setup_mpi
88 : use test_blacs_infrastructure
89 : use test_prepare_matrix
90 : #ifdef HAVE_REDIRECT
91 : use test_redirect
92 : #endif
93 : use test_output_type
94 : implicit none
95 :
96 : !-------------------------------------------------------------------------------
97 : ! Please set system size parameters below!
98 : ! na: System size
99 : ! nev: Number of eigenvectors to be calculated
100 : ! nblk: Blocking factor in block cyclic distribution
101 : !-------------------------------------------------------------------------------
102 :
103 : integer(kind=ik) :: nblk
104 : integer(kind=ik) :: na, nev
105 :
106 : integer(kind=ik) :: np_rows, np_cols, na_rows, na_cols
107 :
108 : integer(kind=ik) :: myid, nprocs, my_prow, my_pcol, mpi_comm_rows, mpi_comm_cols
109 : integer(kind=ik) :: i, mpierr, my_blacs_ctxt, sc_desc(9), info, nprow, npcol
110 :
111 : integer(kind=ik), external :: numroc
112 :
113 192 : real(kind=rk8), allocatable :: a(:,:), z(:,:), as(:,:), ev(:)
114 :
115 : integer(kind=ik) :: STATUS
116 : #ifdef WITH_OPENMP
117 : integer(kind=ik) :: omp_get_max_threads, required_mpi_thread_level, provided_mpi_thread_level
118 : #endif
119 : integer(kind=ik) :: success
120 : integer(kind=ik) :: numberOfDevices
121 : type(output_t) :: write_to_file
122 : character(len=8) :: task_suffix
123 : integer(kind=ik) :: j
124 : integer(kind=ik) :: global_row, global_col, local_row, local_col
125 : integer(kind=ik) :: bandwidth
126 : class(elpa_t), pointer :: e
127 : #define DOUBLE_PRECISION_REAL 1
128 :
129 :
130 192 : call read_input_parameters_traditional(na, nev, nblk, write_to_file)
131 :
132 : !-------------------------------------------------------------------------------
133 : ! MPI Initialization
134 192 : call setup_mpi(myid, nprocs)
135 :
136 192 : STATUS = 0
137 :
138 : #define REALCASE
139 :
140 : !-------------------------------------------------------------------------------
141 : ! Selection of number of processor rows/columns
142 : ! We try to set up the grid square-like, i.e. start the search for possible
143 : ! divisors of nprocs with a number next to the square root of nprocs
144 : ! and decrement it until a divisor is found.
145 :
146 192 : do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
147 0 : if(mod(nprocs,np_cols) == 0 ) exit
148 : enddo
149 : ! at the end of the above loop, nprocs is always divisible by np_cols
150 :
151 192 : np_rows = nprocs/np_cols
152 :
153 192 : if(myid==0) then
154 128 : print *
155 128 : print '(a)','Standard eigenvalue problem - REAL version'
156 128 : print *
157 128 : print '(3(a,i0))','Matrix size=',na,', Number of eigenvectors=',nev,', Block size=',nblk
158 128 : print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
159 128 : print *
160 : endif
161 :
162 : !-------------------------------------------------------------------------------
163 : ! Set up BLACS context and MPI communicators
164 : !
165 : ! The BLACS context is only necessary for using Scalapack.
166 : !
167 : ! For ELPA, the MPI communicators along rows/cols are sufficient,
168 : ! and the grid setup may be done in an arbitrary way as long as it is
169 : ! consistent (i.e. 0<=my_prow<np_rows, 0<=my_pcol<np_cols and every
170 : ! process has a unique (my_prow,my_pcol) pair).
171 :
172 : call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, 'C', &
173 192 : my_blacs_ctxt, my_prow, my_pcol)
174 :
175 192 : if (myid==0) then
176 128 : print '(a)','| Past BLACS_Gridinfo.'
177 : end if
178 :
179 : call set_up_blacs_descriptor(na ,nblk, my_prow, my_pcol, np_rows, np_cols, &
180 192 : na_rows, na_cols, sc_desc, my_blacs_ctxt, info)
181 :
182 192 : if (myid==0) then
183 128 : print '(a)','| Past scalapack descriptor setup.'
184 : end if
185 :
186 : !-------------------------------------------------------------------------------
187 : ! Allocate matrices and set up a test matrix for the eigenvalue problem
188 192 : allocate(a (na_rows,na_cols))
189 192 : allocate(z (na_rows,na_cols))
190 192 : allocate(as(na_rows,na_cols))
191 :
192 192 : allocate(ev(na))
193 :
194 192 : call prepare_matrix_random(na, myid, sc_desc, a, z, as)
195 :
196 : ! set values outside of the bandwidth to zero
197 192 : bandwidth = nblk
198 :
199 19392 : do local_row = 1, na_rows
200 19200 : global_row = index_l2g( local_row, nblk, my_prow, np_rows )
201 2899200 : do local_col = 1, na_cols
202 2880000 : global_col = index_l2g( local_col, nblk, my_pcol, np_cols )
203 :
204 2880000 : if (ABS(global_row-global_col) > bandwidth) then
205 2281216 : a(local_row, local_col) = 0.0
206 2281216 : as(local_row, local_col) = 0.0
207 : end if
208 : end do
209 : end do
210 :
211 192 : if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
212 0 : print *, "ELPA API version not supported"
213 0 : stop 1
214 : endif
215 192 : e => elpa_allocate()
216 :
217 192 : call e%set("na", na, success)
218 192 : assert_elpa_ok(success)
219 192 : call e%set("nev", nev, success)
220 192 : assert_elpa_ok(success)
221 192 : call e%set("local_nrows", na_rows, success)
222 192 : assert_elpa_ok(success)
223 192 : call e%set("local_ncols", na_cols, success)
224 192 : assert_elpa_ok(success)
225 192 : call e%set("nblk", nblk, success)
226 192 : assert_elpa_ok(success)
227 192 : call e%set("mpi_comm_parent", MPI_COMM_WORLD, success)
228 192 : assert_elpa_ok(success)
229 192 : call e%set("process_row", my_prow, success)
230 192 : assert_elpa_ok(success)
231 192 : call e%set("process_col", my_pcol, success)
232 192 : assert_elpa_ok(success)
233 :
234 192 : call e%set("bandwidth", bandwidth, success)
235 192 : assert_elpa_ok(success)
236 :
237 192 : assert(e%setup() .eq. ELPA_OK)
238 :
239 192 : call e%set("solver", ELPA_SOLVER_2STAGE, success)
240 192 : assert_elpa_ok(success)
241 :
242 192 : call e%eigenvectors(a, ev, z, success)
243 192 : assert_elpa_ok(success)
244 192 : call elpa_deallocate(e)
245 :
246 192 : call elpa_uninit()
247 :
248 :
249 : !-------------------------------------------------------------------------------
250 : ! Test correctness of result (using plain scalapack routines)
251 :
252 :
253 192 : status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
254 :
255 :
256 192 : deallocate(a)
257 192 : deallocate(as)
258 :
259 192 : deallocate(z)
260 192 : deallocate(ev)
261 :
262 : #ifdef WITH_MPI
263 128 : call blacs_gridexit(my_blacs_ctxt)
264 128 : call mpi_finalize(mpierr)
265 : #endif
266 192 : call EXIT(STATUS)
267 : end
268 :
269 : !-------------------------------------------------------------------------------
|