test.F90 13.2 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
118
119
120
121
122
123
124
125
126
127

   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
128
   integer :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol 
129
130
131
132
133
134
135

   ! The Matrix
   MATRIX_TYPE, allocatable :: a(:,:), as(:,:)
   ! eigenvectors
   MATRIX_TYPE, allocatable :: z(:,:)
   ! eigenvalues
   EV_TYPE, allocatable :: ev(:)
136
137
#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)
   EV_TYPE, allocatable :: d(:), sd(:), ev_analytic(:), ds(:), sds(:)
Andreas Marek's avatar
Andreas Marek committed
138
   EV_TYPE              :: diagonalELement, subdiagonalElement
139
140
#endif

141

142
   integer :: error, status
143
144
145

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

156
#if defined(TEST_COMPLEX) && defined(__SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
157
158
159
#ifdef WITH_MPI
   call MPI_finalize(mpierr)
#endif
160
161
   stop 77
#endif
162
   call read_input_parameters_traditional(na, nev, nblk, write_to_file)
Pavel Kus's avatar
Pavel Kus committed
163

164
165
   call setup_mpi(myid, nprocs)

166
167
168
169
170
171
172
173
174
175
176
177
178
179
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

#ifdef TEST_ALL_LAYOUTS
   do i_layout = 1, size(layouts)               ! layout loop
     layout = layouts(i_layout)
     do np_cols = 1, nprocs                     ! factor loop
       if (mod(nprocs,np_cols) /= 0 ) then
         cycle
       endif
#else
   layout = 'C'
180
181
182
   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
183
#endif
184
185

   np_rows = nprocs/np_cols
186
   assert(nprocs == np_rows * np_cols)
187

188
189
190
191
192
193
   if (myid == 0) then
     print '((a,i0))', 'Matrix size: ', na
     print '((a,i0))', 'Num eigenvectors: ', nev
     print '((a,i0))', 'Blocksize: ', nblk
     print '((a,i0))', 'Num MPI proc: ', nprocs
     print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
194
     print '(a)',      'Process layout: ' // layout
195
196
197
     print *,''
   endif

198
199
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
200
201
202
203

   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
204
205
206
207
   allocate(a (na_rows,na_cols))
#ifdef TEST_MATRIX_RANDOM
   allocate(as(na_rows,na_cols))
#endif
208
209
210
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

211
212
213
214
215
216
#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
   allocate(ev_analytic(na))
#endif

217
218
219
220
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

221
#ifdef __EIGENVECTORS
Pavel Kus's avatar
Pavel Kus committed
222
223
224
#ifdef TEST_MATRIX_ANALYTIC
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
#else
225
   call prepare_matrix(na, myid, sc_desc, a, z, as)
226
#endif
Pavel Kus's avatar
Pavel Kus committed
227
#endif
228
229

#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
230
231
232
233
234
235
236
237
238
239
240

#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)
241
#endif
242
243
244

   e => elpa_allocate()

245
246
247
248
249
250
251
252
253
254
   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)
255
256

#ifdef WITH_MPI
257
258
259
260
261
262
   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)
263
#endif
264

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

267
268
269
270
271
272
273
   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
274
   assert_elpa_ok(error)
275

276
277
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
278

279
280
281
#ifdef TEST_ALL_KERNELS
   do i = 0, elpa_option_cardinality(KERNEL_KEY)
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
282
#endif
283
#ifdef TEST_KERNEL
284
     kernel = TEST_KERNEL
285
#endif
286
287

#ifdef TEST_SOLVER_2STAGE
288
     call e%set(KERNEL_KEY, kernel, error)
289
290
291
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
292
293
294
     if (error /= ELPA_OK) then
       cycle
     endif
295
#endif
296
     if (myid == 0) print *, elpa_int_value_to_string(KERNEL_KEY, kernel), " kernel"
297
298

     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
Andreas Marek's avatar
Andreas Marek committed
299
#else /* ALL_KERNELS */
300
301

#ifdef __EIGENVECTORS
302
     call e%timer_start("e%eigenvectors()")
303
304
305
306
307
308
309
#endif
#ifdef __EIGENVALUES
     call e%timer_start("e%eigenvalues()")
#endif
#ifdef __SOLVE_TRIDIAGONAL
     call e%timer_start("e%solve_tridiagonal()")
#endif
310
#endif
311

312
     ! The actual solve step
313
#ifdef __EIGENVECTORS
314
     call e%eigenvectors(a, ev, z, error)
315
316
317
318
319
320
321
322
323
#endif
#ifdef __EIGENVALUES
     call e%eigenvalues(a, ev, error)
#endif
#if defined(__SOLVE_TRIDIAGONAL) && !defined(TEST_COMPLEX)
     call e%solve_tridiagonal(d, sd, z, error)
     ev(:) = d(:)
#endif

324
325
     assert_elpa_ok(error)

326
327
328
#ifdef TEST_SOLVER_2STAGE
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#else
329
#ifdef __EIGENVECTORS
330
     call e%timer_stop("e%eigenvectors()")
331
332
333
334
335
336
337
#endif
#ifdef __EIGENVALUES
     call e%timer_stop("e%eigenvalues()")
#endif
#ifdef __SOLVE_TRIDIAGONAL
     call e%timer_stop("e%solve_tridiagonal()")
#endif
338
339
#endif

340
     if (myid .eq. 0) then
341
342
343
#ifdef TEST_SOLVER_2STAGE
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
#else
344
#ifdef __EIGENVECTORS
345
       call e%print_times("e%eigenvectors()")
346
347
348
349
350
351
352
#endif
#ifdef __EIGENVALUES
       call e%print_times("e%eigenvalues()")
#endif
#ifdef __SOLVE_TRIDIAGONAL
     call e%print_times("e%solve_tridiagonal()")
#endif
353
#endif
354
     endif
355
356

#ifdef __EIGENVECTORS
Pavel Kus's avatar
Pavel Kus committed
357
358
359
#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
360
     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
361
#endif
362
     if (status /= 0) then
Andreas Marek's avatar
Andreas Marek committed
363
       if (myid == 0) print *, "Result incorrect!"
364
365
       call exit(status)
     endif
Andreas Marek's avatar
Andreas Marek committed
366
     if (myid == 0) print *, ""
367
368
369
#endif
#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)

Andreas Marek's avatar
Andreas Marek committed
370
371
     status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
372
373
374
375
376
377

     if (status /= 0) then
       call exit(status)
     endif
#ifdef __SOLVE_TRIDIAGONAL
     ! check eigenvectors
Pavel Kus's avatar
Pavel Kus committed
378
     status = check_correctness(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
379
380
381
382
383
384
385
386
     if (status /= 0) then
       if (myid == 0) print *, "Result incorrect!"
       call exit(status)
     endif
     if (myid == 0) print *, ""
#endif

#endif
387
388

#ifdef TEST_ALL_KERNELS
Pavel Kus's avatar
Pavel Kus committed
389
390
391
#ifdef TEST_MATRIX_ANALYTIC
     call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
#else
392
     a(:,:) = as(:,:)
Pavel Kus's avatar
Pavel Kus committed
393
#endif
394
395
396
397
#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)
     d = ds
     sd = sds
#endif
398
399
   end do
#endif
Andreas Marek's avatar
Andreas Marek committed
400

401
402
403
   call elpa_deallocate(e)

   deallocate(a)
Pavel Kus's avatar
Pavel Kus committed
404
#ifdef TEST_MATRIX_RANDOM
405
   deallocate(as)
Pavel Kus's avatar
Pavel Kus committed
406
#endif
407
408
409
   deallocate(z)
   deallocate(ev)

410
#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)
411
412
413
414
415
   deallocate(d, ds)
   deallocate(sd, sds)
   deallocate(ev_analytic)
#endif

416
417
418
419
420
421
422
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif

   call elpa_uninit()

423
424
425
426
427
428
429
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif

   call exit(status)

430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
!#if defined(__EIGENVALUES) || defined(__SOLVE_TRIDIAGONAL)
!   contains
!
!     !Processor col for global col number
!     pure function pcol(global_col, nblk, np_cols) result(local_col)
!       implicit none
!       integer(kind=c_int), intent(in) :: global_col, nblk, np_cols
!       integer(kind=c_int)             :: local_col
!       local_col = MOD((global_col-1)/nblk,np_cols)
!     end function
!
!     !Processor row for global row number
!     pure function prow(global_row, nblk, np_rows) result(local_row)
!       implicit none
!       integer(kind=c_int), intent(in) :: global_row, nblk, np_rows
!       integer(kind=c_int)             :: local_row
!       local_row = MOD((global_row-1)/nblk,np_rows)
!     end function
!
!     function map_global_array_index_to_local_index(iGLobal, jGlobal, iLocal, jLocal , nblk, np_rows, np_cols, my_prow, my_pcol) &
!       result(possible)
!       implicit none
!
!       integer(kind=c_int)              :: pi, pj, li, lj, xi, xj
!       integer(kind=c_int), intent(in)  :: iGlobal, jGlobal, nblk, np_rows, np_cols, my_prow, my_pcol
!       integer(kind=c_int), intent(out) :: iLocal, jLocal
!       logical                       :: possible
!
!       possible = .true.
!       iLocal = 0
!       jLocal = 0
!
!       pi = prow(iGlobal, nblk, np_rows)
!
!       if (my_prow .ne. pi) then
!         possible = .false.
!         return
!       endif
!
!       pj = pcol(jGlobal, nblk, np_cols)
!
!       if (my_pcol .ne. pj) then
!         possible = .false.
!         return
!       endif
!       li = (iGlobal-1)/(np_rows*nblk) ! block number for rows
!       lj = (jGlobal-1)/(np_cols*nblk) ! block number for columns
!
!       xi = mod( (iGlobal-1),nblk)+1   ! offset in block li
!       xj = mod( (jGlobal-1),nblk)+1   ! offset in block lj
!
!       iLocal = li * nblk + xi
!       jLocal = lj * nblk + xj
!
!     end function
!#endif
486
end program