test.F90 11.6 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - Max Planck Computing and Data Facility (MPCDF), formerly known as
!      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
!      Schwerpunkt Wissenschaftliches Rechnen ,
!    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
!    http://elpa.mpcdf.mpg.de/
!
!    ELPA is free software: you can redistribute it and/or modify
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
!    Software Foundation.
!
!    ELPA is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!    GNU Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public License
!    along with ELPA.  If not, see <http://www.gnu.org/licenses/>
!
!    ELPA reflects a substantial effort on the part of the original
!    ELPA consortium, and we ask you to respect the spirit of the
!    license that we chose: i.e., please contribute any changes you
!    may have back to the original ELPA library distribution, and keep
!    any derivatives of ELPA under the same license that we chose for
!    the original distribution, the GNU Lesser General Public License.
!
!
#include "config-f90.h"

! Define one of TEST_REAL or TEST_COMPLEX
! Define one of TEST_SINGLE or TEST_DOUBLE
! Define one of TEST_SOLVER_1STAGE or TEST_SOLVER_2STAGE
! Define TEST_GPU \in [0, 1]
49
! Define either TEST_ALL_KERNELS or a TEST_KERNEL \in [any valid kernel]
50
51
52
53
54
55
56
57
58
59
60
61
62

#if !(defined(TEST_REAL) ^ defined(TEST_COMPLEX))
error: define exactly one of TEST_REAL or TEST_COMPLEX
#endif

#if !(defined(TEST_SINGLE) ^ defined(TEST_DOUBLE))
error: define exactly one of TEST_SINGLE or TEST_DOUBLE
#endif

#if !(defined(TEST_SOLVER_1STAGE) ^ defined(TEST_SOLVER_2STAGE))
error: define exactly one of TEST_SOLVER_1STAGE or TEST_SOLVER_2STAGE
#endif

63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
#ifdef TEST_SOLVER_1STAGE
#ifdef TEST_ALL_KERNELS
error: TEST_ALL_KERNELS cannot be defined for TEST_SOLVER_1STAGE
#endif
#ifdef TEST_KERNEL
error: TEST_KERNEL cannot be defined for TEST_SOLVER_1STAGE
#endif
#endif

#ifdef TEST_SOLVER_2STAGE
#if !(defined(TEST_KERNEL) ^ defined(TEST_ALL_KERNELS))
error: define either TEST_ALL_KERNELS or a valid TEST_KERNEL
#endif
#endif


79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
#ifdef TEST_SINGLE
#  define EV_TYPE real(kind=C_FLOAT)
#  ifdef TEST_REAL
#    define MATRIX_TYPE real(kind=C_FLOAT)
#  else
#    define MATRIX_TYPE complex(kind=C_FLOAT_COMPLEX)
#  endif
#else
#  define EV_TYPE real(kind=C_DOUBLE)
#  ifdef TEST_REAL
#    define MATRIX_TYPE real(kind=C_DOUBLE)
#  else
#    define MATRIX_TYPE complex(kind=C_DOUBLE_COMPLEX)
#  endif
#endif

95
96
97
98
99
100
101
#ifdef TEST_REAL
#define KERNEL_KEY "real_kernel"
#endif
#ifdef TEST_COMPLEX
#define KERNEL_KEY "complex_kernel"
#endif

102
103
104
105
#include "assert.h"

program test
   use elpa
106
107
108
109
110
111
112

   use test_util
   use test_setup_mpi
   use test_prepare_matrix
   use test_read_input_parameters
   use test_blacs_infrastructure
   use test_check_correctness
Pavel Kus's avatar
Pavel Kus committed
113
   use test_analytic
114

115
116
117
#ifdef HAVE_REDIRECT
   use test_redirect
#endif
118
119
120
121
122
123
124
125
126
127
128
129
130
   implicit none

   ! matrix dimensions
   integer :: na, nev, nblk

   ! mpi
   integer :: myid, nprocs
   integer :: na_cols, na_rows  ! local matrix size
   integer :: np_cols, np_rows  ! number of MPI processes per column/row
   integer :: my_prow, my_pcol  ! local MPI task position (my_prow, my_pcol) in the grid (0..np_cols -1, 0..np_rows -1)
   integer :: mpierr

   ! blacs
131
   integer :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol 
132
133
134
135
136
137
138

   ! The Matrix
   MATRIX_TYPE, allocatable :: a(:,:), as(:,:)
   ! eigenvectors
   MATRIX_TYPE, allocatable :: z(:,:)
   ! eigenvalues
   EV_TYPE, allocatable :: ev(:)
139
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
140
   EV_TYPE, allocatable :: d(:), sd(:), ev_analytic(:), ds(:), sds(:)
Andreas Marek's avatar
Andreas Marek committed
141
   EV_TYPE              :: diagonalELement, subdiagonalElement
142
143
#endif

144

145
   integer :: error, status
146
147
148

   type(output_t) :: write_to_file
   class(elpa_t), pointer :: e
149
#ifdef TEST_ALL_KERNELS
150
   integer :: i
151
152
153
154
#endif
#ifdef TEST_ALL_LAYOUTS
   character(len=1), parameter :: layouts(2) = [ 'C', 'R' ]
   integer :: i_layout
155
#endif
156
   integer :: kernel
157
   character(len=1) :: layout
158

159
   call read_input_parameters_traditional(na, nev, nblk, write_to_file)
160
   call setup_mpi(myid, nprocs)
161
162
163
164
165
166
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
#endif
     call redirect_stdout(myid)
#endif
167

168
169
170
171
172
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

173
174
175
176
177
   if (myid == 0) then
     print '((a,i0))', 'Program ' // TEST_CASE
     print *, ""
   endif

178
#ifdef TEST_ALL_LAYOUTS
179
   do i_layout = 1, size(layouts)               ! layouts
180
     layout = layouts(i_layout)
181
     do np_cols = 1, nprocs                     ! factors
182
183
184
185
186
       if (mod(nprocs,np_cols) /= 0 ) then
         cycle
       endif
#else
   layout = 'C'
187
188
189
   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
190
#endif
191
192

   np_rows = nprocs/np_cols
193
   assert(nprocs == np_rows * np_cols)
194

195
196
197
198
   if (myid == 0) then
     print '((a,i0))', 'Matrix size: ', na
     print '((a,i0))', 'Num eigenvectors: ', nev
     print '((a,i0))', 'Blocksize: ', nblk
199
#ifdef WITH_MPI
200
201
     print '((a,i0))', 'Num MPI proc: ', nprocs
     print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
202
     print '(a)',      'Process layout: ' // layout
203
#endif
204
205
206
     print *,''
   endif

207
208
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
209
210
211
212

   call set_up_blacs_descriptor(na, nblk, my_prow, my_pcol, np_rows, np_cols, &
                                na_rows, na_cols, sc_desc, my_blacs_ctxt, info)

Pavel Kus's avatar
Pavel Kus committed
213
214
   allocate(a (na_rows,na_cols))
   allocate(as(na_rows,na_cols))
215
216
217
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

218
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
219
220
221
222
223
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
   allocate(ev_analytic(na))
#endif

224
225
226
227
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

228
#ifdef TEST_EIGENVECTORS
Pavel Kus's avatar
Pavel Kus committed
229
230
#ifdef TEST_MATRIX_ANALYTIC
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
231
   as(:,:) = a
Pavel Kus's avatar
Pavel Kus committed
232
#else
233
   call prepare_matrix(na, myid, sc_desc, a, z, as)
234
#endif
Pavel Kus's avatar
Pavel Kus committed
235
#endif
236

237
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
238
239
240
241
242
243
244
245
246
247
248

#ifdef TEST_SINGLE
   diagonalElement = 0.45_c_float
   subdiagonalElement =  0.78_c_float
#else
   diagonalElement = 0.45_c_double
   subdiagonalElement =  0.78_c_double
#endif
   call prepare_toeplitz_matrix(na, diagonalElement, subdiagonalElement, &
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
249
#endif
250
251
252

   e => elpa_allocate()

253
254
255
256
257
258
259
260
261
262
   call e%set("na", na, error)
   assert_elpa_ok(error)
   call e%set("nev", nev, error)
   assert_elpa_ok(error)
   call e%set("local_nrows", na_rows, error)
   assert_elpa_ok(error)
   call e%set("local_ncols", na_cols, error)
   assert_elpa_ok(error)
   call e%set("nblk", nblk, error)
   assert_elpa_ok(error)
263
264

#ifdef WITH_MPI
265
266
267
268
269
270
   call e%set("mpi_comm_parent", MPI_COMM_WORLD, error)
   assert_elpa_ok(error)
   call e%set("process_row", my_prow, error)
   assert_elpa_ok(error)
   call e%set("process_col", my_pcol, error)
   assert_elpa_ok(error)
271
#endif
272

Andreas Marek's avatar
Andreas Marek committed
273
274
   call e%set("timings",1)

275
276
277
278
279
280
281
   assert_elpa_ok(e%setup())

#ifdef TEST_SOLVER_1STAGE
   call e%set("solver", ELPA_SOLVER_1STAGE)
#else
   call e%set("solver", ELPA_SOLVER_2STAGE)
#endif
282
   assert_elpa_ok(error)
283

284
285
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
286

287
288
   if (myid == 0) print *, ""

289
#ifdef TEST_ALL_KERNELS
290
   do i = 0, elpa_option_cardinality(KERNEL_KEY)  ! kernels
291
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
292
#endif
293
#ifdef TEST_KERNEL
294
     kernel = TEST_KERNEL
295
#endif
296
297

#ifdef TEST_SOLVER_2STAGE
298
     call e%set(KERNEL_KEY, kernel, error)
299
300
301
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
302
303
304
     if (error /= ELPA_OK) then
       cycle
     endif
305
306
307
308
309
310
     ! actually used kernel might be different if forced via environment variables
     call e%get(KERNEL_KEY, kernel)
#endif
     if (myid == 0) then
       print *, elpa_int_value_to_string(KERNEL_KEY, kernel) // " kernel"
     endif
311
312
#endif

313
#ifdef TEST_ALL_KERNELS
314
     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
315
#endif
316

317
     ! The actual solve step
318
319
#ifdef TEST_EIGENVECTORS
     call e%timer_start("e%eigenvectors()")
320
     call e%eigenvectors(a, ev, z, error)
321
     call e%timer_stop("e%eigenvectors()")
322
#endif
323
324
325

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
326
     call e%eigenvalues(a, ev, error)
327
     call e%timer_stop("e%eigenvalues()")
328
#endif
329
330
331

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
332
     call e%solve_tridiagonal(d, sd, z, error)
333
     call e%timer_stop("e%solve_tridiagonal()")
334
335
336
     ev(:) = d(:)
#endif

337
338
     assert_elpa_ok(error)

339
#ifdef TEST_ALL_KERNELS
340
341
342
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

343
     if (myid .eq. 0) then
344
#ifdef TEST_ALL_KERNELS
345
346
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
#else
347
#ifdef TEST_EIGENVECTORS
348
       call e%print_times("e%eigenvectors()")
349
#endif
350
#ifdef TEST_EIGENVALUES
351
352
       call e%print_times("e%eigenvalues()")
#endif
353
354
#ifdef TEST_SOLVE_TRIDIAGONAL
       call e%print_times("e%solve_tridiagonal()")
355
#endif
356
#endif
357
     endif
358

359
#ifdef TEST_EIGENVECTORS
Pavel Kus's avatar
Pavel Kus committed
360
361
362
#ifdef TEST_MATRIX_ANALYTIC
     status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
#else
Pavel Kus's avatar
Pavel Kus committed
363
     status = check_correctness(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, my_pcol)
Andreas Marek's avatar
Andreas Marek committed
364
#endif
365
     call check_status(status, myid)
366
367
#endif

368
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
369
370
     status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
371
     call check_status(status, myid)
372

373
#ifdef TEST_SOLVE_TRIDIAGONAL
374
     ! check eigenvectors
Pavel Kus's avatar
Pavel Kus committed
375
     status = check_correctness(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
376
     call check_status(status, myid)
377
378
#endif
#endif
379

380
381
382
383
     if (myid == 0) then
       print *, ""
     endif

384
385
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
386
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
387
388
389
     d = ds
     sd = sds
#endif
390
   end do ! kernels
391
#endif
Andreas Marek's avatar
Andreas Marek committed
392

393
394
395
396
397
398
399
   call elpa_deallocate(e)

   deallocate(a)
   deallocate(as)
   deallocate(z)
   deallocate(ev)

400
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
401
402
403
404
405
   deallocate(d, ds)
   deallocate(sd, sds)
   deallocate(ev_analytic)
#endif

406
407
408
409
410
411
412
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif

   call elpa_uninit()

413
414
415
416
417
418
419
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif

   call exit(status)

420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
   contains

     subroutine check_status(status, myid)
       implicit none
       integer, intent(in) :: status, myid
       integer :: mpierr
       if (status /= 0) then
         if (myid == 0) print *, "Result incorrect!"
#ifdef WITH_MPI
         call mpi_finalize(mpierr)
#endif
         call exit(status)
       endif
     end subroutine

435
end program