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 complex 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 complex ELPA 2 kernel is set as the default kernel.
64 : !> However, this can be overriden by setting
65 : !> the environment variable "COMPLEX_ELPA_KERNEL" to an
66 : !> appropiate value.
67 : !>
68 96 : program test_complex2_single_banded
69 :
70 : !-------------------------------------------------------------------------------
71 : ! Standard eigenvalue problem - COMPLEX 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 96 : use elpa
82 :
83 : use test_util
84 : use test_read_input_parameters
85 : use test_check_correctness
86 : use test_setup_mpi
87 : use test_blacs_infrastructure
88 : use test_prepare_matrix
89 : #ifdef HAVE_REDIRECT
90 : use test_redirect
91 : #endif
92 :
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 : #ifdef WITH_MPI
111 : integer(kind=ik), external :: numroc
112 : #endif
113 : complex(kind=ck4), parameter :: CZERO = (0.0_rk4,0.0_rk4), CONE = (1.0_rk4,0.0_rk4)
114 96 : real(kind=rk4), allocatable :: ev(:)
115 :
116 192 : complex(kind=ck4), allocatable :: a(:,:), z(:,:), as(:,:)
117 :
118 : integer(kind=ik) :: STATUS
119 : #ifdef WITH_OPENMP
120 : integer(kind=ik) :: omp_get_max_threads, required_mpi_thread_level, provided_mpi_thread_level
121 : #endif
122 : type(output_t) :: write_to_file
123 : integer(kind=ik) :: success
124 : character(len=8) :: task_suffix
125 : integer(kind=ik) :: j
126 :
127 :
128 : integer(kind=ik) :: global_row, global_col, local_row, local_col
129 : integer(kind=ik) :: bandwidth
130 : class(elpa_t), pointer :: e
131 :
132 : #define COMPLEXCASE
133 : #define DOUBLE_PRECISION_COMPLEX 1
134 :
135 96 : call read_input_parameters_traditional(na, nev, nblk, write_to_file)
136 : !-------------------------------------------------------------------------------
137 : ! MPI Initialization
138 96 : call setup_mpi(myid, nprocs)
139 :
140 96 : STATUS = 0
141 :
142 : !-------------------------------------------------------------------------------
143 : ! Selection of number of processor rows/columns
144 : ! We try to set up the grid square-like, i.e. start the search for possible
145 : ! divisors of nprocs with a number next to the square root of nprocs
146 : ! and decrement it until a divisor is found.
147 :
148 96 : do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
149 0 : if(mod(nprocs,np_cols) == 0 ) exit
150 : enddo
151 : ! at the end of the above loop, nprocs is always divisible by np_cols
152 :
153 96 : np_rows = nprocs/np_cols
154 :
155 96 : if(myid==0) then
156 64 : print *
157 64 : print '(a)','Standard eigenvalue problem - COMPLEX version'
158 64 : print *
159 64 : print '(3(a,i0))','Matrix size=',na,', Number of eigenvectors=',nev,', Block size=',nblk
160 64 : print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
161 64 : print *
162 : endif
163 :
164 : !-------------------------------------------------------------------------------
165 : ! Set up BLACS context and MPI communicators
166 : !
167 : ! The BLACS context is only necessary for using Scalapack.
168 : !
169 : ! For ELPA, the MPI communicators along rows/cols are sufficient,
170 : ! and the grid setup may be done in an arbitrary way as long as it is
171 : ! consistent (i.e. 0<=my_prow<np_rows, 0<=my_pcol<np_cols and every
172 : ! process has a unique (my_prow,my_pcol) pair).
173 :
174 : call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, 'C', &
175 96 : my_blacs_ctxt, my_prow, my_pcol)
176 :
177 96 : if (myid==0) then
178 64 : print '(a)','| Past BLACS_Gridinfo.'
179 : end if
180 :
181 : ! Determine the necessary size of the distributed matrices,
182 : ! we use the Scalapack tools routine NUMROC for that.
183 :
184 : call set_up_blacs_descriptor(na ,nblk, my_prow, my_pcol, np_rows, np_cols, &
185 96 : na_rows, na_cols, sc_desc, my_blacs_ctxt, info)
186 :
187 96 : if (myid==0) then
188 64 : print '(a)','| Past scalapack descriptor setup.'
189 : end if
190 : !-------------------------------------------------------------------------------
191 : ! Allocate matrices and set up a test matrix for the eigenvalue problem
192 :
193 96 : allocate(a (na_rows,na_cols))
194 96 : allocate(z (na_rows,na_cols))
195 96 : allocate(as(na_rows,na_cols))
196 :
197 96 : allocate(ev(na))
198 :
199 96 : call prepare_matrix_random(na, myid, sc_desc, a, z, as)
200 :
201 : ! set values outside of the bandwidth to zero
202 96 : bandwidth = nblk
203 :
204 9696 : do local_row = 1, na_rows
205 9600 : global_row = index_l2g( local_row, nblk, my_prow, np_rows )
206 1449600 : do local_col = 1, na_cols
207 1440000 : global_col = index_l2g( local_col, nblk, my_pcol, np_cols )
208 :
209 1440000 : if (ABS(global_row-global_col) > bandwidth) then
210 1140608 : a(local_row, local_col) = 0
211 1140608 : as(local_row, local_col) = 0
212 : end if
213 : end do
214 : end do
215 :
216 96 : if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
217 0 : print *, "ELPA API version not supported"
218 0 : stop 1
219 : endif
220 :
221 96 : e => elpa_allocate()
222 :
223 96 : call e%set("na", na, success)
224 96 : assert_elpa_ok(success)
225 96 : call e%set("nev", nev, success)
226 96 : assert_elpa_ok(success)
227 96 : call e%set("local_nrows", na_rows, success)
228 96 : assert_elpa_ok(success)
229 96 : call e%set("local_ncols", na_cols, success)
230 96 : assert_elpa_ok(success)
231 96 : call e%set("nblk", nblk, success)
232 96 : assert_elpa_ok(success)
233 96 : call e%set("mpi_comm_parent", MPI_COMM_WORLD, success)
234 96 : assert_elpa_ok(success)
235 96 : call e%set("process_row", my_prow, success)
236 96 : assert_elpa_ok(success)
237 96 : call e%set("process_col", my_pcol, success)
238 96 : assert_elpa_ok(success)
239 :
240 96 : call e%set("bandwidth", bandwidth, success)
241 96 : assert_elpa_ok(success)
242 :
243 96 : assert(e%setup() .eq. ELPA_OK)
244 :
245 96 : call e%set("solver", ELPA_SOLVER_2STAGE, success)
246 96 : assert_elpa_ok(success)
247 96 : call e%eigenvectors(a, ev, z, success)
248 96 : assert_elpa_ok(success)
249 96 : call elpa_deallocate(e)
250 :
251 96 : call elpa_uninit()
252 :
253 : !-------------------------------------------------------------------------------
254 : ! Test correctness of result (using plain scalapack routines)
255 96 : status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
256 :
257 96 : deallocate(a)
258 96 : deallocate(as)
259 :
260 96 : deallocate(z)
261 96 : deallocate(ev)
262 :
263 : #ifdef WITH_MPI
264 64 : call blacs_gridexit(my_blacs_ctxt)
265 64 : call mpi_finalize(mpierr)
266 : #endif
267 96 : call EXIT(STATUS)
268 : end
269 :
270 : !-------------------------------------------------------------------------------
|