test.F90 26.1 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
#if !(defined(TEST_REAL) ^ defined(TEST_COMPLEX))
error: define exactly one of TEST_REAL or TEST_COMPLEX
53
54
55
56
57
58
#endif

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

59
60
#if !(defined(TEST_SOLVER_1STAGE) ^ defined(TEST_SOLVER_2STAGE) ^ defined(TEST_SCALAPACK_ALL) ^ defined(TEST_SCALAPACK_PART))
error: define exactly one of TEST_SOLVER_1STAGE or TEST_SOLVER_2STAGE or TEST_SCALAPACK_ALL or TEST_SCALAPACK_PART
61
62
#endif

63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
#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

78
79
80
#ifdef TEST_GENERALIZED_DECOMP_EIGENPROBLEM
#define TEST_GENERALIZED_EIGENPROBLEM
#endif
81

82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
#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

98
99
100
101
102
103
104
#ifdef TEST_REAL
#define KERNEL_KEY "real_kernel"
#endif
#ifdef TEST_COMPLEX
#define KERNEL_KEY "complex_kernel"
#endif

105
106
107
108
#include "assert.h"

program test
   use elpa
109
110
111
112
113
114
115

   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
116
   use test_analytic
117
#ifdef WITH_SCALAPACK_TESTS
Pavel Kus's avatar
Pavel Kus committed
118
   use test_scalapack
119
#endif
120

121
122
#ifdef HAVE_REDIRECT
   use test_redirect
123
124
125
#endif
#ifdef WITH_OPENMP
   use omp_lib
126
#endif
127
128
129
   implicit none

   ! matrix dimensions
130
   integer                     :: na, nev, nblk
131
132

   ! mpi
133
134
135
136
137
   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
138
139

   ! blacs
140
   integer                     :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol
141
142

   ! The Matrix
143
   MATRIX_TYPE, allocatable    :: a(:,:), as(:,:)
144
#if defined(TEST_HERMITIAN_MULTIPLY)
145
   MATRIX_TYPE, allocatable    :: b(:,:), c(:,:)
146
147
148
#endif
#if defined(TEST_GENERALIZED_EIGENPROBLEM)
   MATRIX_TYPE, allocatable    :: b(:,:), bs(:,:)
149
#endif
150
   ! eigenvectors
151
   MATRIX_TYPE, allocatable    :: z(:,:)
152
   ! eigenvalues
153
   EV_TYPE, allocatable        :: ev(:)
154

155
   logical                     :: check_all_evals, skip_check_correctness
156

157
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
158
159
   EV_TYPE, allocatable        :: d(:), sd(:), ds(:), sds(:)
   EV_TYPE                     :: diagonalELement, subdiagonalElement
160
#endif
161

162
   integer                     :: error, status
163

164
165
   type(output_t)              :: write_to_file
   class(elpa_t), pointer      :: e
166
#ifdef TEST_ALL_KERNELS
167
   integer                     :: i
168
169
170
#endif
#ifdef TEST_ALL_LAYOUTS
   character(len=1), parameter :: layouts(2) = [ 'C', 'R' ]
171
172
173
174
   integer                     :: i_layout
#endif
   integer                     :: kernel
   character(len=1)            :: layout
175
176
   logical                     :: do_test_numeric_residual, do_test_numeric_residual_generalized, &
                                  do_test_analytic_eigenvalues, &
Andreas Marek's avatar
Andreas Marek committed
177
178
179
180
181
                                  do_test_analytic_eigenvalues_eigenvectors,   &
                                  do_test_frank_eigenvalues,  &
                                  do_test_toeplitz_eigenvalues, do_test_cholesky,   &
                                  do_test_hermitian_multiply

182
#ifdef WITH_OPENMP
183
   integer                    :: max_threads, threads_caller
184
185
#endif

186
187
188
189
190
#ifdef SPLIT_COMM_MYSELF
   integer                    :: mpi_comm_rows, mpi_comm_cols, mpi_string_length, mpierr2
   character(len=MPI_MAX_ERROR_STRING) :: mpierr_string
#endif

191
   call read_input_parameters_traditional(na, nev, nblk, write_to_file, skip_check_correctness)
192
   call setup_mpi(myid, nprocs)
193
194
195
196
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
     call redirect_stdout(myid)
197
#endif
198
#endif
199

200
   check_all_evals = .true.
201

Andreas Marek's avatar
Andreas Marek committed
202
203

   do_test_numeric_residual = .false.
204
   do_test_numeric_residual_generalized = .false.
Andreas Marek's avatar
Andreas Marek committed
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
   do_test_analytic_eigenvalues = .false.
   do_test_analytic_eigenvalues_eigenvectors = .false.
   do_test_frank_eigenvalues = .false.
   do_test_toeplitz_eigenvalues = .false. 

   do_test_cholesky = .false.
#if defined(TEST_CHOLESKY)
   do_test_cholesky = .true.
#endif
   do_test_hermitian_multiply = .false.
#if defined(TEST_HERMITIAN_MULTIPLY)
   do_test_hermitian_multiply = .true.
#endif

   status = 0
220
221
222
223
224
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

225
226
227
228
229
   if (myid == 0) then
     print '((a,i0))', 'Program ' // TEST_CASE
     print *, ""
   endif

230
#ifdef TEST_ALL_LAYOUTS
231
   do i_layout = 1, size(layouts)               ! layouts
232
     layout = layouts(i_layout)
233
     do np_cols = 1, nprocs                     ! factors
234
235
236
237
238
       if (mod(nprocs,np_cols) /= 0 ) then
         cycle
       endif
#else
   layout = 'C'
239
240
241
   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
242
#endif
243
244

   np_rows = nprocs/np_cols
245
   assert(nprocs == np_rows * np_cols)
246

247
248
249
250
   if (myid == 0) then
     print '((a,i0))', 'Matrix size: ', na
     print '((a,i0))', 'Num eigenvectors: ', nev
     print '((a,i0))', 'Blocksize: ', nblk
251
#ifdef WITH_MPI
252
253
     print '((a,i0))', 'Num MPI proc: ', nprocs
     print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
254
     print '(a)',      'Process layout: ' // layout
255
#endif
256
257
258
     print *,''
   endif

Andreas Marek's avatar
Andreas Marek committed
259
#if TEST_QR_DECOMPOSITION == 1
260
261
262
263
264
265

#if TEST_GPU == 1
#ifdef WITH_MPI
     call mpi_finalize(mpierr)
#endif
     stop 77
Andreas Marek's avatar
Andreas Marek committed
266
#endif /* TEST_GPU */
267
268
269
270
271
272
273
274
275
276
277
278
279
   if (nblk .lt. 64) then
     if (myid .eq. 0) then
       print *,"At the moment QR decomposition need blocksize of at least 64"
     endif
     if ((na .lt. 64) .and. (myid .eq. 0)) then
       print *,"This is why the matrix size must also be at least 64 or only 1 MPI task can be used"
     endif

#ifdef WITH_MPI
     call mpi_finalize(mpierr)
#endif
     stop 77
   endif
Andreas Marek's avatar
Andreas Marek committed
280
281
#endif /* TEST_QR_DECOMPOSITION */

282

283
284
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
285
286
287
288

   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
289
290
   allocate(a (na_rows,na_cols))
   allocate(as(na_rows,na_cols))
291
292
293
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

294
295
296
297
298
#ifdef TEST_HERMITIAN_MULTIPLY
   allocate(b (na_rows,na_cols))
   allocate(c (na_rows,na_cols))
#endif

299
300
301
302
303
#ifdef TEST_GENERALIZED_EIGENPROBLEM
   allocate(b (na_rows,na_cols))
   allocate(bs (na_rows,na_cols))
#endif

304
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
305
306
307
308
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
#endif

309
310
311
312
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

Andreas Marek's avatar
Andreas Marek committed
313
314
315
316
317
#if defined(TEST_MATRIX_RANDOM) && !defined(TEST_SOLVE_TRIDIAGONAL) && !defined(TEST_CHOLESKY) && !defined(TEST_EIGENVALUES)
   ! the random matrix can be used in allmost all tests; but for some no
   ! correctness checks have been implemented; do not allow these
   ! combinations
   ! RANDOM + TEST_SOLVE_TRIDIAGONAL: we need a TOEPLITZ MATRIX
318
   ! RANDOM + TEST_CHOLESKY: wee need SPD matrix
Andreas Marek's avatar
Andreas Marek committed
319
320
321
322
323
324
325
   ! RANDOM + TEST_EIGENVALUES: no correctness check known

   ! We also have to take care of special case in TEST_EIGENVECTORS
#if !defined(TEST_EIGENVECTORS)
    call prepare_matrix_random(na, myid, sc_desc, a, z, as)
#else /* TEST_EIGENVECTORS */
    if (nev .ge. 1) then
326
     call prepare_matrix_random(na, myid, sc_desc, a, z, as)
Andreas Marek's avatar
Andreas Marek committed
327
328
329
#ifndef TEST_HERMITIAN_MULTIPLY
    do_test_numeric_residual = .true.
#endif
330
   else
Andreas Marek's avatar
Andreas Marek committed
331
332
333
334
335
     if (myid .eq. 0) then
       print *,"At the moment with the random matrix you need nev >=1"
     endif
#ifdef WITH_MPI
     call mpi_finalize(mpierr)
336
#endif
Andreas Marek's avatar
Andreas Marek committed
337
338
339
     stop 77
   endif
#endif /* TEST_EIGENVECTORS */
340
341
342
343
    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
    do_test_frank_eigenvalues = .false.
    do_test_toeplitz_eigenvalues = .false.
Andreas Marek's avatar
Andreas Marek committed
344
#endif /* (TEST_MATRIX_RANDOM) */
345
346
347

#if defined(TEST_MATRIX_RANDOM) && defined(TEST_CHOLESKY)
     call prepare_matrix_random_spd(na, myid, sc_desc, a, z, as, &
348
                 nblk, np_rows, np_cols, my_prow, my_pcol)
349
350
351
352
353
354
    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
    do_test_frank_eigenvalues = .false.
    do_test_toeplitz_eigenvalues = .false.
#endif /* TEST_MATRIX_RANDOM and TEST_CHOLESKY */

355
356
357
358
359
360
361
362
363
364
365
366
#if defined(TEST_MATRIX_RANDOM) && defined(TEST_GENERALIZED_EIGENPROBLEM)
   ! call prepare_matrix_random(na, myid, sc_desc, a, z, as)
    call prepare_matrix_random_spd(na, myid, sc_desc, b, z, bs, &
                 nblk, np_rows, np_cols, my_prow, my_pcol)
    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
    do_test_frank_eigenvalues = .false.
    do_test_toeplitz_eigenvalues = .false.
    do_test_numeric_residual = .false.
    do_test_numeric_residual_generalized = .true.
#endif /* TEST_MATRIX_RANDOM and TEST_GENERALIZED_EIGENPROBLEM */

367
#if defined(TEST_MATRIX_RANDOM) && (defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVALUES))
Andreas Marek's avatar
Andreas Marek committed
368
#error "Random matrix is not allowed in this configuration"
369
#endif
370

Andreas Marek's avatar
Andreas Marek committed
371
372
373
374
375
376
#if defined(TEST_MATRIX_ANALYTIC)  && !defined(TEST_SOLVE_TRIDIAGONAL) && !defined(TEST_CHOLESKY)
   ! the analytic matrix can be used in allmost all tests; but for some no
   ! correctness checks have been implemented; do not allow these
   ! combinations
   ! ANALYTIC + TEST_SOLVE_TRIDIAGONAL: we need a TOEPLITZ MATRIX
   ! ANALTIC  + TEST_CHOLESKY: no correctness check yet implemented
377

Andreas Marek's avatar
Andreas Marek committed
378
379
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
   as(:,:) = a
380

Andreas Marek's avatar
Andreas Marek committed
381
382
383
384
   do_test_numeric_residual = .false.
   do_test_analytic_eigenvalues_eigenvectors = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
   do_test_analytic_eigenvalues = .true.
385
#endif
Andreas Marek's avatar
Andreas Marek committed
386
387
388
#if defined(TEST_EIGENVECTORS)
   if (nev .ge. 1) then
     do_test_analytic_eigenvalues_eigenvectors = .true.
Pavel Kus's avatar
Pavel Kus committed
389
     do_test_analytic_eigenvalues = .false.
Andreas Marek's avatar
Andreas Marek committed
390
391
392
   else
     do_test_analytic_eigenvalues_eigenvectors = .false.
   endif
393
#endif
Andreas Marek's avatar
Andreas Marek committed
394
395
   do_test_frank_eigenvalues = .false.
   do_test_toeplitz_eigenvalues = .false.
396
#endif /* TEST_MATRIX_ANALYTIC */
Andreas Marek's avatar
Andreas Marek committed
397
398
399
#if defined(TEST_MATRIX_ANALYTIC) && (defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_CHOLESKY))
#error "Analytic matrix is not allowd in this configuration"
#endif
Andreas Marek's avatar
Andreas Marek committed
400

Andreas Marek's avatar
Andreas Marek committed
401
402
#if defined(TEST_MATRIX_TOEPLITZ)
   ! The Toeplitz matrix works in each test
Andreas Marek's avatar
Andreas Marek committed
403
404
405
406
407
408
409
#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
410

411
! actually we test cholesky for diagonal matrix only
412
413
#if defined(TEST_CHOLESKY)
#ifdef TEST_SINGLE
414
415
  diagonalElement = (2.546_c_float, 0.0_c_float)
  subdiagonalElement =  (0.0_c_float, 0.0_c_float)
416
#else
417
418
  diagonalElement = (2.546_c_double, 0.0_c_double)
  subdiagonalElement =  (0.0_c_double, 0.0_c_double)
419
#endif
Andreas Marek's avatar
Andreas Marek committed
420
421
#endif /* TEST_CHOLESKY */

422
   call prepare_matrix_toeplitz(na, diagonalElement, subdiagonalElement, &
423
424
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
425
426


Andreas Marek's avatar
Andreas Marek committed
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
   do_test_numeric_residual = .false.
#if defined(TEST_EIGENVECTORS)
   if (nev .ge. 1) then
     do_test_numeric_residual = .true.
   else
     do_test_numeric_residual = .false.
   endif
#endif

   do_test_analytic_eigenvalues = .false.
   do_test_analytic_eigenvalues_eigenvectors = .false.
   do_test_frank_eigenvalues = .false.
#if defined(TEST_CHOLESKY)
   do_test_toeplitz_eigenvalues = .false.
#else
   do_test_toeplitz_eigenvalues = .true.
#endif
444

Andreas Marek's avatar
Andreas Marek committed
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
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
#endif /* TEST_MATRIX_TOEPLITZ */


#if defined(TEST_MATRIX_FRANK) && !defined(TEST_SOLVE_TRIDIAGONAL) && !defined(TEST_CHOLESKY)
   ! the random matrix can be used in allmost all tests; but for some no
   ! correctness checks have been implemented; do not allow these
   ! combinations
   ! FRANK + TEST_SOLVE_TRIDIAGONAL: we need a TOEPLITZ MATRIX
   ! FRANK + TEST_CHOLESKY: no correctness check yet implemented

   ! We also have to take care of special case in TEST_EIGENVECTORS
#if !defined(TEST_EIGENVECTORS)
    call prepare_matrix_frank(na, a, z, as, nblk, np_rows, np_cols, my_prow, my_pcol)

    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
    do_test_frank_eigenvalues = .true.
#endif
    do_test_toeplitz_eigenvalues = .false.

#else /* TEST_EIGENVECTORS */

    if (nev .ge. 1) then
     call prepare_matrix_frank(na, a, z, as, nblk, np_rows, np_cols, my_prow, my_pcol)

    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
    do_test_frank_eigenvalues = .true.
#endif
    do_test_toeplitz_eigenvalues = .false.
    do_test_numeric_residual = .false.
   else
    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
    do_test_frank_eigenvalues = .true.
#endif
    do_test_toeplitz_eigenvalues = .false.
    do_test_numeric_residual = .false.

   endif

#endif /* TEST_EIGENVECTORS */
#endif /* (TEST_MATRIX_FRANK) */
#if defined(TEST_MATRIX_FRANK) && (defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_CHOLESKY))
#error "FRANK matrix is not allowed in this configuration"
#endif


#ifdef TEST_HERMITIAN_MULTIPLY
#ifdef TEST_REAL

#ifdef TEST_DOUBLE
   b(:,:) = 2.0_c_double * a(:,:)
   c(:,:) = 0.0_c_double
#else
   b(:,:) = 2.0_c_float * a(:,:)
   c(:,:) = 0.0_c_float
#endif

#endif /* TEST_REAL */

#ifdef TEST_COMPLEX

#ifdef TEST_DOUBLE
   b(:,:) = 2.0_c_double * a(:,:)
   c(:,:) = (0.0_c_double, 0.0_c_double)
#else
   b(:,:) = 2.0_c_float * a(:,:)
   c(:,:) = (0.0_c_float, 0.0_c_float)
#endif

#endif /* TEST_COMPLEX */

#endif /* TEST_HERMITIAN_MULTIPLY */
522

523
524
525
526
527
528
529
530
531
532
533
534
! if the test is used for (repeated) performacne tests, one might want to skip the checking
! of the results, which might be time-consuming and not necessary.
   if(skip_check_correctness) then
     do_test_numeric_residual = .false.
     do_test_numeric_residual_generalized = .false.
     do_test_analytic_eigenvalues = .false.
     do_test_analytic_eigenvalues_eigenvectors = .false.
     do_test_frank_eigenvalues = .false.
     do_test_toeplitz_eigenvalues = .false. 
     do_test_cholesky = .false.
   endif

535
536
537
538
539
540
541
542

#ifdef WITH_OPENMP
   threads_caller = omp_get_max_threads()
   if (myid == 0) then
     print *,"The calling program uses ",threads_caller," threads"
   endif
#endif

543
544
   e => elpa_allocate(error)
   assert_elpa_ok(error)
545

546
547
548
549
550
551
552
553
554
555
   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)
556
557

#ifdef WITH_MPI
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
#ifdef SPLIT_COMM_MYSELF
   call mpi_comm_split(MPI_COMM_WORLD,my_pcol,my_prow,mpi_comm_rows,mpierr)
   if (mpierr .ne. MPI_SUCCESS) then
     call MPI_ERROR_STRING(mpierr,mpierr_string, mpi_string_length, mpierr2)
     write(error_unit,*) "MPI ERROR occured during mpi_comm_split for row communicator: ", trim(mpierr_string)
     stop 1
   endif

   call mpi_comm_split(MPI_COMM_WORLD,my_prow,my_pcol,mpi_comm_cols, mpierr)
   if (mpierr .ne. MPI_SUCCESS) then
     call MPI_ERROR_STRING(mpierr,mpierr_string, mpi_string_length, mpierr2)
     write(error_unit,*) "MPI ERROR occured during mpi_comm_split for col communicator: ", trim(mpierr_string)
     stop 1
   endif

573
574
   call e%set("mpi_comm_parent", MPI_COMM_WORLD, error)
   assert_elpa_ok(error)
575
576
577
578
579
   call e%set("mpi_comm_rows", mpi_comm_rows, error)
   assert_elpa_ok(error)
   call e%set("mpi_comm_cols", mpi_comm_cols, error)
   assert_elpa_ok(error)
#else
580
581
582
583
584
585
   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)
586
#endif
587
#endif
588
589
590
#ifdef TEST_GENERALIZED_EIGENPROBLEM
   call e%set("blacs_context", my_blacs_ctxt, error)
   assert_elpa_ok(error)
591
#endif
Andreas Marek's avatar
Andreas Marek committed
592
   call e%set("timings",1,error)
593
594
595
   assert_elpa_ok(e%setup())

#ifdef TEST_SOLVER_1STAGE
Andreas Marek's avatar
Andreas Marek committed
596
   call e%set("solver", ELPA_SOLVER_1STAGE,error)
597
#else
Andreas Marek's avatar
Andreas Marek committed
598
   call e%set("solver", ELPA_SOLVER_2STAGE,error)
599
#endif
600
   assert_elpa_ok(error)
601

602
603
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
604

Andreas Marek's avatar
Andreas Marek committed
605
#if TEST_QR_DECOMPOSITION == 1
606
607
608
609
   call e%set("qr", 1, error)
   assert_elpa_ok(error)
#endif

610
611
612
613
614
615
#ifdef WITH_OPENMP
   max_threads=omp_get_max_threads()
   call e%set("omp_threads", max_threads, error)
   assert_elpa_ok(error)
#endif

616
617
   if (myid == 0) print *, ""

618
#ifdef TEST_ALL_KERNELS
619
   do i = 0, elpa_option_cardinality(KERNEL_KEY)  ! kernels
620
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
621
622
     if (kernel .eq. ELPA_2STAGE_REAL_GPU) continue
     if (kernel .eq. ELPA_2STAGE_COMPLEX_GPU) continue
623
#endif
624
#ifdef TEST_KERNEL
625
     kernel = TEST_KERNEL
626
#endif
627
628

#ifdef TEST_SOLVER_2STAGE
629
     call e%set(KERNEL_KEY, kernel, error)
630
631
632
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
633
634
635
     if (error /= ELPA_OK) then
       cycle
     endif
636
     ! actually used kernel might be different if forced via environment variables
Andreas Marek's avatar
Andreas Marek committed
637
     call e%get(KERNEL_KEY, kernel, error)
638
     assert_elpa_ok(error)
639
640
641
642
#endif
     if (myid == 0) then
       print *, elpa_int_value_to_string(KERNEL_KEY, kernel) // " kernel"
     endif
643
644
#endif

Pavel Kus's avatar
Pavel Kus committed
645
646

! print all parameters
647
648
     call e%print_settings(error)
     assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
649

650
#ifdef TEST_ALL_KERNELS
651
     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
652
#endif
653

654
     ! The actual solve step
Andreas Marek's avatar
Andreas Marek committed
655
656
657
658
#if defined(TEST_EIGENVECTORS)
#if TEST_QR_DECOMPOSITION == 1
     call e%timer_start("e%eigenvectors_qr()")
#else
659
     call e%timer_start("e%eigenvectors()")
Andreas Marek's avatar
Andreas Marek committed
660
#endif
Pavel Kus's avatar
Pavel Kus committed
661
662
#ifdef TEST_SCALAPACK_ALL
     call solve_scalapack_all(na, a, sc_desc, ev, z)
663
664
#elif TEST_SCALAPACK_PART
     call solve_scalapack_part(na, a, sc_desc, nev, ev, z)
665
     check_all_evals = .false. ! scalapack does not compute all eigenvectors
Pavel Kus's avatar
Pavel Kus committed
666
#else
667
     call e%eigenvectors(a, ev, z, error)
Pavel Kus's avatar
Pavel Kus committed
668
#endif
Andreas Marek's avatar
Andreas Marek committed
669
670
671
#if TEST_QR_DECOMPOSITION == 1
     call e%timer_stop("e%eigenvectors_qr()")
#else
672
     call e%timer_stop("e%eigenvectors()")
Andreas Marek's avatar
Andreas Marek committed
673
674
#endif
#endif /* TEST_EIGENVECTORS  */
675
676
677

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
678
     call e%eigenvalues(a, ev, error)
679
     call e%timer_stop("e%eigenvalues()")
680
#endif
681
682
683

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
684
     call e%solve_tridiagonal(d, sd, z, error)
685
     call e%timer_stop("e%solve_tridiagonal()")
686
687
688
     ev(:) = d(:)
#endif

689
690
691
#if defined(TEST_CHOLESKY)
     call e%timer_start("e%cholesky()")
     call e%cholesky(a, error)
692
     assert_elpa_ok(error)
693
694
695
     call e%timer_stop("e%cholesky()")
#endif

696
697
698
699
700
#if defined(TEST_HERMITIAN_MULTIPLY)
     call e%timer_start("e%hermitian_multiply()")
     call e%hermitian_multiply('F','F', na, a, b, na_rows, na_cols, c, na_rows, na_cols, error)
     call e%timer_stop("e%hermitian_multiply()")
#endif
Pavel Kus's avatar
Pavel Kus committed
701

702
703
#if defined(TEST_GENERALIZED_EIGENPROBLEM)
     call e%timer_start("e%generalized_eigenvectors()")
704
705
706
#if defined(TEST_GENERALIZED_DECOMP_EIGENPROBLEM)
     call e%timer_start("is_already_decomposed=.false.")
#endif
707
     call e%generalized_eigenvectors(a, b, ev, z, .false., error)
708
709
710
711
#if defined(TEST_GENERALIZED_DECOMP_EIGENPROBLEM)
     call e%timer_stop("is_already_decomposed=.false.")
     a = as
     call e%timer_start("is_already_decomposed=.true.")
712
     call e%generalized_eigenvectors(a, b, ev, z, .true., error)
713
714
     call e%timer_stop("is_already_decomposed=.true.")
#endif
715
716
717
     call e%timer_stop("e%generalized_eigenvectors()")
#endif

718
719
     assert_elpa_ok(error)

720
#ifdef TEST_ALL_KERNELS
721
722
723
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

724
     if (myid .eq. 0) then
725
#ifdef TEST_ALL_KERNELS
726
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
727
728
#else /* TEST_ALL_KERNELS */

Andreas Marek's avatar
Andreas Marek committed
729
730
731
732
#if defined(TEST_EIGENVECTORS)
#if TEST_QR_DECOMPOSITION == 1
       call e%print_times("e%eigenvectors_qr()")
#else
733
       call e%print_times("e%eigenvectors()")
734
#endif
Andreas Marek's avatar
Andreas Marek committed
735
#endif
736
#ifdef TEST_EIGENVALUES
737
738
       call e%print_times("e%eigenvalues()")
#endif
739
740
#ifdef TEST_SOLVE_TRIDIAGONAL
       call e%print_times("e%solve_tridiagonal()")
741
#endif
742
743
#ifdef TEST_CHOLESKY
       call e%print_times("e%cholesky()")
744
#endif
745
746
747
#ifdef TEST_HERMITIAN_MULTIPLY
       call e%print_times("e%hermitian_multiply()")
#endif
748
749
750
#ifdef TEST_GENERALIZED_EIGENPROBLEM
      call e%print_times("e%generalized_eigenvectors()")
#endif
751
#endif /* TEST_ALL_KERNELS */
752
     endif
753

Andreas Marek's avatar
Andreas Marek committed
754
755
756
757
758
759
760
761
762
     if (do_test_analytic_eigenvalues) then
       status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol, check_all_evals, .false.)
       call check_status(status, myid)
     endif

     if (do_test_analytic_eigenvalues_eigenvectors) then
       status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol, check_all_evals, .true.)
       call check_status(status, myid)
     endif
763

Andreas Marek's avatar
Andreas Marek committed
764
     if(do_test_numeric_residual) then
765
       status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, my_pcol)
Andreas Marek's avatar
Andreas Marek committed
766
767
768
769
770
771
       call check_status(status, myid)
     endif

     if (do_test_frank_eigenvalues) then
       status = check_correctness_eigenvalues_frank(na, ev, z, myid)
       call check_status(status, myid)
772
     endif
773

Andreas Marek's avatar
Andreas Marek committed
774
     if (do_test_toeplitz_eigenvalues) then
775
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
776
       status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
777
                                                       subdiagonalElement, ev, z, myid)
Andreas Marek's avatar
Andreas Marek committed
778
       call check_status(status, myid)
779
#endif
Andreas Marek's avatar
Andreas Marek committed
780
     endif
781

Andreas Marek's avatar
Andreas Marek committed
782
783
784
785
     if (do_test_cholesky) then
       status = check_correctness_cholesky(na, a, as, na_rows, sc_desc, myid )
       call check_status(status, myid)
     endif
786

Andreas Marek's avatar
Andreas Marek committed
787
788
789
790
791
#ifdef TEST_HERMITIAN_MULTIPLY
     if (do_test_hermitian_multiply) then
       status = check_correctness_hermitian_multiply(na, a, b, c, na_rows, sc_desc, myid )
       call check_status(status, myid)
     endif
792
#endif
793

794
795
796
797
798
799
800
801
#ifdef TEST_GENERALIZED_EIGENPROBLEM
     if(do_test_numeric_residual_generalized) then
       status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, &
       my_pcol, bs)
       call check_status(status, myid)
     endif
#endif

802
803
804
805
806
807
808
809
810

#ifdef WITH_OPENMP
     if (threads_caller .ne. omp_get_max_threads()) then
       if (myid .eq. 0) then
         print *, " ERROR! the number of OpenMP threads has not been restored correctly"
       endif
       status = 1
     endif
#endif
811
812
813
814
     if (myid == 0) then
       print *, ""
     endif

815
816
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
817
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
818
819
820
     d = ds
     sd = sds
#endif
821
   end do ! kernels
822
#endif
Andreas Marek's avatar
Andreas Marek committed
823

824
825
   call elpa_deallocate(e, error)
   assert_elpa_ok(error)
826
827
828
829
830

   deallocate(a)
   deallocate(as)
   deallocate(z)
   deallocate(ev)
831
832
833
834
#ifdef TEST_HERMITIAN_MULTIPLY
   deallocate(b)
   deallocate(c)
#endif
835
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
836
837
838
   deallocate(d, ds)
   deallocate(sd, sds)
#endif
839
840
841
#if defined(TEST_GENERALIZED_EIGENPROBLEM)
  deallocate(b, bs)
#endif
842

843
844
845
846
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif
847
848
   call elpa_uninit(error)
   assert_elpa_ok(error)
849

850
851
852
853
854
855
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif
   call exit(status)

856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
   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

871
end program