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 : !>
45 : !> Fortran test programm to demonstrates the use of
46 : !> ELPA 2 complex case library.
47 : !> If "HAVE_REDIRECT" was defined at build time
48 : !> the stdout and stderr output of each MPI task
49 : !> can be redirected to files if the environment
50 : !> variable "REDIRECT_ELPA_TEST_OUTPUT" is set
51 : !> to "true".
52 : !>
53 : !> By calling executable [arg1] [arg2] [arg3] [arg4]
54 : !> one can define the size (arg1), the number of
55 : !> Eigenvectors to compute (arg2), and the blocking (arg3).
56 : !> If these values are not set default values (500, 150, 16)
57 : !> are choosen.
58 : !> If these values are set the 4th argument can be
59 : !> "output", which specifies that the EV's are written to
60 : !> an ascii file.
61 : !>
62 : !> The complex ELPA 2 kernel is set as the default kernel.
63 : !> However, this can be overriden by setting
64 : !> the environment variable "COMPLEX_ELPA_KERNEL" to an
65 : !> appropiate value.
66 : !>
67 96 : program test_complex2
68 :
69 : !-------------------------------------------------------------------------------
70 : ! Standard eigenvalue problem - COMPLEX version
71 : !
72 : ! This program demonstrates the use of the ELPA module
73 : ! together with standard scalapack routines
74 : !
75 : ! Copyright of the original code rests with the authors inside the ELPA
76 : ! consortium. The copyright of any additional modifications shall rest
77 : ! with their original authors, but shall adhere to the licensing terms
78 : ! distributed along with the original code in the file "COPYING".
79 : !-------------------------------------------------------------------------------
80 96 : use elpa_driver
81 : use elpa_utilities, only : error_unit
82 : use test_util
83 : use test_read_input_parameters
84 : use test_check_correctness
85 : use test_setup_mpi
86 : use test_blacs_infrastructure
87 : use test_prepare_matrix
88 : #ifdef HAVE_REDIRECT
89 : use test_redirect
90 : #endif
91 :
92 : use test_output_type
93 : implicit none
94 :
95 : !-------------------------------------------------------------------------------
96 : ! Please set system size parameters below!
97 : ! na: System size
98 : ! nev: Number of eigenvectors to be calculated
99 : ! nblk: Blocking factor in block cyclic distribution
100 : !-------------------------------------------------------------------------------
101 :
102 : integer(kind=ik) :: nblk
103 : integer(kind=ik) :: na, nev
104 :
105 : integer(kind=ik) :: np_rows, np_cols, na_rows, na_cols
106 :
107 : integer(kind=ik) :: myid, nprocs, my_prow, my_pcol, mpi_comm_rows, mpi_comm_cols
108 : integer(kind=ik) :: i, mpierr, my_blacs_ctxt, sc_desc(9), info, nprow, npcol
109 :
110 : integer, external :: numroc
111 :
112 : complex(kind=ck4), parameter :: CZERO = (0.0,0.0), CONE = (1.0,0.0)
113 96 : real(kind=rk4), allocatable :: ev(:)
114 :
115 192 : complex(kind=ck4), allocatable :: a(:,:), z(:,:), as(:,:)
116 :
117 : integer(kind=ik) :: STATUS
118 : #ifdef WITH_OPENMP
119 : integer(kind=ik) :: omp_get_max_threads, required_mpi_thread_level, provided_mpi_thread_level
120 : #endif
121 : type(output_t) :: write_to_file
122 : logical :: success
123 : character(len=8) :: task_suffix
124 : integer(kind=ik) :: j
125 :
126 96 : success = .true.
127 :
128 96 : call read_input_parameters_traditional(na, nev, nblk, write_to_file)
129 : !-------------------------------------------------------------------------------
130 : ! MPI Initialization
131 96 : call setup_mpi(myid, nprocs)
132 96 : STATUS = 0
133 :
134 : #define DATATYPE COMPLEX
135 : #include "../../elpa_print_headers.F90"
136 :
137 : !-------------------------------------------------------------------------------
138 : ! Selection of number of processor rows/columns
139 : ! We try to set up the grid square-like, i.e. start the search for possible
140 : ! divisors of nprocs with a number next to the square root of nprocs
141 : ! and decrement it until a divisor is found.
142 :
143 96 : do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
144 0 : if(mod(nprocs,np_cols) == 0 ) exit
145 : enddo
146 : ! at the end of the above loop, nprocs is always divisible by np_cols
147 :
148 96 : np_rows = nprocs/np_cols
149 :
150 96 : if(myid==0) then
151 64 : print *
152 64 : print '(a)','Standard eigenvalue problem - COMPLEX version'
153 64 : print *
154 64 : print '(3(a,i0))','Matrix size=',na,', Number of eigenvectors=',nev,', Block size=',nblk
155 64 : print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
156 64 : print *
157 : endif
158 :
159 : !-------------------------------------------------------------------------------
160 : ! Set up BLACS context and MPI communicators
161 : !
162 : ! The BLACS context is only necessary for using Scalapack.
163 : !
164 : ! For ELPA, the MPI communicators along rows/cols are sufficient,
165 : ! and the grid setup may be done in an arbitrary way as long as it is
166 : ! consistent (i.e. 0<=my_prow<np_rows, 0<=my_pcol<np_cols and every
167 : ! process has a unique (my_prow,my_pcol) pair).
168 :
169 : call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, 'C', &
170 96 : my_blacs_ctxt, my_prow, my_pcol)
171 :
172 96 : if (myid==0) then
173 64 : print '(a)','| Past BLACS_Gridinfo.'
174 : end if
175 :
176 : ! All ELPA routines need MPI communicators for communicating within
177 : ! rows or columns of processes, these are set in elpa_get_communicators.
178 :
179 : mpierr = elpa_get_communicators(mpi_comm_world, my_prow, my_pcol, &
180 96 : mpi_comm_rows, mpi_comm_cols)
181 :
182 96 : if (myid==0) then
183 64 : print '(a)','| Past split communicator setup for rows and columns.'
184 : end if
185 :
186 : ! Determine the necessary size of the distributed matrices,
187 : ! we use the Scalapack tools routine NUMROC for that.
188 :
189 : call set_up_blacs_descriptor(na ,nblk, my_prow, my_pcol, np_rows, np_cols, &
190 96 : na_rows, na_cols, sc_desc, my_blacs_ctxt, info)
191 :
192 96 : if (myid==0) then
193 64 : print '(a)','| Past scalapack descriptor setup.'
194 : end if
195 : !-------------------------------------------------------------------------------
196 : ! Allocate matrices and set up a test matrix for the eigenvalue problem
197 :
198 96 : allocate(a (na_rows,na_cols))
199 96 : allocate(z (na_rows,na_cols))
200 96 : allocate(as(na_rows,na_cols))
201 :
202 96 : allocate(ev(na))
203 :
204 96 : call prepare_matrix_random(na, myid, sc_desc, a, z, as)
205 :
206 : ! set print flag in elpa1
207 96 : elpa_print_times = .true.
208 :
209 : !-------------------------------------------------------------------------------
210 : ! Calculate eigenvalues/eigenvectors
211 :
212 96 : if (myid==0) then
213 64 : print '(a)','| Entering one-stage ELPA solver ... '
214 64 : print *
215 : end if
216 :
217 : #ifdef WITH_MPI
218 64 : call mpi_barrier(mpi_comm_world, mpierr) ! for correct timings only
219 : #endif
220 : success = elpa_solve_evp_complex_single(na, nev, a, na_rows, ev, z, na_rows, nblk, &
221 : na_cols, &
222 96 : mpi_comm_rows, mpi_comm_cols, mpi_comm_world,method="1stage")
223 :
224 96 : if (.not.(success)) then
225 0 : write(error_unit,*) "elpa_solve_evp_complex produced an error! Aborting..."
226 : #ifdef WITH_MPI
227 0 : call MPI_ABORT(mpi_comm_world, 1, mpierr)
228 : #else
229 0 : call exit(1)
230 : #endif
231 : endif
232 :
233 96 : if (myid==0) then
234 64 : print '(a)','| One-step ELPA solver complete.'
235 64 : print *
236 : end if
237 :
238 96 : a = as
239 96 : z = as
240 :
241 96 : if (myid==0) then
242 64 : print '(a)','| Entering two-stage ELPA solver ... '
243 64 : print *
244 : end if
245 :
246 : #ifdef WITH_MPI
247 64 : call mpi_barrier(mpi_comm_world, mpierr) ! for correct timings only
248 : #endif
249 : success = elpa_solve_evp_complex_single(na, nev, a, na_rows, ev, z, na_rows, nblk, &
250 : na_cols, &
251 96 : mpi_comm_rows, mpi_comm_cols, mpi_comm_world,method="2stage")
252 :
253 96 : if (.not.(success)) then
254 0 : write(error_unit,*) "elpa_solve_evp_complex produced an error! Aborting..."
255 : #ifdef WITH_MPI
256 0 : call MPI_ABORT(mpi_comm_world, 1, mpierr)
257 : #else
258 0 : call exit(1)
259 : #endif
260 : endif
261 :
262 96 : if (myid==0) then
263 64 : print '(a)','| Two-step ELPA solver complete.'
264 64 : print *
265 : end if
266 :
267 96 : a = as
268 96 : z = as
269 :
270 96 : if (myid==0) then
271 64 : print '(a)','| Entering auto-chosen ELPA solver ... '
272 64 : print *
273 : end if
274 :
275 : #ifdef WITH_MPI
276 64 : call mpi_barrier(mpi_comm_world, mpierr) ! for correct timings only
277 : #endif
278 : success = elpa_solve_evp_complex_single(na, nev, a, na_rows, ev, z, na_rows, nblk, &
279 : na_cols, &
280 96 : mpi_comm_rows, mpi_comm_cols, mpi_comm_world,method="auto")
281 :
282 96 : if (.not.(success)) then
283 0 : write(error_unit,*) "elpa_solve_evp_complex produced an error! Aborting..."
284 : #ifdef WITH_MPI
285 0 : call MPI_ABORT(mpi_comm_world, 1, mpierr)
286 : #else
287 0 : call exit(1)
288 : #endif
289 : endif
290 :
291 96 : if (myid==0) then
292 64 : print '(a)','| Auto-chosen ELPA solver complete.'
293 64 : print *
294 : end if
295 :
296 :
297 96 : if(myid == 0) print *,'Time transform to tridi :',time_evp_fwd
298 96 : if(myid == 0) print *,'Time solve tridi :',time_evp_solve
299 96 : if(myid == 0) print *,'Time transform back EVs :',time_evp_back
300 96 : if(myid == 0) print *,'Total time (sum above) :',time_evp_back+time_evp_solve+time_evp_fwd
301 :
302 : ! if(write_to_file%eigenvectors) then
303 : ! write(unit = task_suffix, fmt = '(i8.8)') myid
304 : ! open(17,file="EVs_complex2_out_task_"//task_suffix(1:8)//".txt",form='formatted',status='new')
305 : ! write(17,*) "Part of eigenvectors: na_rows=",na_rows,"of na=",na," na_cols=",na_cols," of na=",na
306 : !
307 : ! do i=1,na_rows
308 : ! do j=1,na_cols
309 : ! write(17,*) "row=",i," col=",j," element of eigenvector=",z(i,j)
310 : ! enddo
311 : ! enddo
312 : ! close(17)
313 : ! endif
314 : ! if(write_to_file%eigenvalues) then
315 : ! if (myid == 0) then
316 : ! open(17,file="Eigenvalues_complex2_out.txt",form='formatted',status='new')
317 : ! do i=1,na
318 : ! write(17,*) i,ev(i)
319 : ! enddo
320 : ! close(17)
321 : ! endif
322 : ! endif
323 :
324 : !-------------------------------------------------------------------------------
325 : ! Test correctness of result (using plain scalapack routines)
326 96 : status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
327 :
328 96 : deallocate(a)
329 96 : deallocate(as)
330 :
331 96 : deallocate(z)
332 96 : deallocate(ev)
333 :
334 : #ifdef WITH_MPI
335 64 : call blacs_gridexit(my_blacs_ctxt)
336 64 : call mpi_finalize(mpierr)
337 : #endif
338 96 : call EXIT(STATUS)
339 : end
340 :
341 : !-------------------------------------------------------------------------------
|