test.F90 24.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
123
#ifdef HAVE_REDIRECT
   use test_redirect
#endif
124
125
126
   implicit none

   ! matrix dimensions
127
   integer                     :: na, nev, nblk
128
129

   ! mpi
130
131
132
133
134
   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
135
136

   ! blacs
137
   integer                     :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol
138
139

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

152
   logical                     :: check_all_evals, skip_check_correctness
153

154
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
155
156
   EV_TYPE, allocatable        :: d(:), sd(:), ds(:), sds(:)
   EV_TYPE                     :: diagonalELement, subdiagonalElement
157
#endif
158

159
   integer                     :: error, status
160

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

179
   call read_input_parameters_traditional(na, nev, nblk, write_to_file, skip_check_correctness)
180
   call setup_mpi(myid, nprocs)
181
182
183
184
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
     call redirect_stdout(myid)
185
#endif
186
#endif
187

188
   check_all_evals = .true.
189

Andreas Marek's avatar
Andreas Marek committed
190
191

   do_test_numeric_residual = .false.
192
   do_test_numeric_residual_generalized = .false.
Andreas Marek's avatar
Andreas Marek committed
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
   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
208
209
210
211
212
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

213
214
215
216
217
   if (myid == 0) then
     print '((a,i0))', 'Program ' // TEST_CASE
     print *, ""
   endif

218
#ifdef TEST_ALL_LAYOUTS
219
   do i_layout = 1, size(layouts)               ! layouts
220
     layout = layouts(i_layout)
221
     do np_cols = 1, nprocs                     ! factors
222
223
224
225
226
       if (mod(nprocs,np_cols) /= 0 ) then
         cycle
       endif
#else
   layout = 'C'
227
228
229
   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
230
#endif
231
232

   np_rows = nprocs/np_cols
233
   assert(nprocs == np_rows * np_cols)
234

235
236
237
238
   if (myid == 0) then
     print '((a,i0))', 'Matrix size: ', na
     print '((a,i0))', 'Num eigenvectors: ', nev
     print '((a,i0))', 'Blocksize: ', nblk
239
#ifdef WITH_MPI
240
241
     print '((a,i0))', 'Num MPI proc: ', nprocs
     print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
242
     print '(a)',      'Process layout: ' // layout
243
#endif
244
245
246
     print *,''
   endif

Andreas Marek's avatar
Andreas Marek committed
247
#if TEST_QR_DECOMPOSITION == 1
248
249
250
251
252
253

#if TEST_GPU == 1
#ifdef WITH_MPI
     call mpi_finalize(mpierr)
#endif
     stop 77
Andreas Marek's avatar
Andreas Marek committed
254
#endif /* TEST_GPU */
255
256
257
258
259
260
261
262
263
264
265
266
267
   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
268
269
#endif /* TEST_QR_DECOMPOSITION */

270

271
272
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
273
274
275
276

   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
277
278
   allocate(a (na_rows,na_cols))
   allocate(as(na_rows,na_cols))
279
280
281
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

282
283
284
285
286
#ifdef TEST_HERMITIAN_MULTIPLY
   allocate(b (na_rows,na_cols))
   allocate(c (na_rows,na_cols))
#endif

287
288
289
290
291
#ifdef TEST_GENERALIZED_EIGENPROBLEM
   allocate(b (na_rows,na_cols))
   allocate(bs (na_rows,na_cols))
#endif

292
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
293
294
295
296
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
#endif

297
298
299
300
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

Andreas Marek's avatar
Andreas Marek committed
301
302
303
304
305
#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
306
   ! RANDOM + TEST_CHOLESKY: wee need SPD matrix
Andreas Marek's avatar
Andreas Marek committed
307
308
309
310
311
312
313
   ! 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
314
     call prepare_matrix_random(na, myid, sc_desc, a, z, as)
Andreas Marek's avatar
Andreas Marek committed
315
316
317
#ifndef TEST_HERMITIAN_MULTIPLY
    do_test_numeric_residual = .true.
#endif
318
   else
Andreas Marek's avatar
Andreas Marek committed
319
320
321
322
323
     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)
324
#endif
Andreas Marek's avatar
Andreas Marek committed
325
326
327
     stop 77
   endif
#endif /* TEST_EIGENVECTORS */
328
329
330
331
    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
332
#endif /* (TEST_MATRIX_RANDOM) */
333
334
335

#if defined(TEST_MATRIX_RANDOM) && defined(TEST_CHOLESKY)
     call prepare_matrix_random_spd(na, myid, sc_desc, a, z, as, &
336
                 nblk, np_rows, np_cols, my_prow, my_pcol)
337
338
339
340
341
342
    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 */

343
344
345
346
347
348
349
350
351
352
353
354
#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 */

355
#if defined(TEST_MATRIX_RANDOM) && (defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVALUES))
Andreas Marek's avatar
Andreas Marek committed
356
#error "Random matrix is not allowed in this configuration"
357
#endif
358

Andreas Marek's avatar
Andreas Marek committed
359
360
361
362
363
364
#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
365

Andreas Marek's avatar
Andreas Marek committed
366
367
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
   as(:,:) = a
368

Andreas Marek's avatar
Andreas Marek committed
369
370
371
372
   do_test_numeric_residual = .false.
   do_test_analytic_eigenvalues_eigenvectors = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
   do_test_analytic_eigenvalues = .true.
373
#endif
Andreas Marek's avatar
Andreas Marek committed
374
375
376
#if defined(TEST_EIGENVECTORS)
   if (nev .ge. 1) then
     do_test_analytic_eigenvalues_eigenvectors = .true.
Pavel Kus's avatar
Pavel Kus committed
377
     do_test_analytic_eigenvalues = .false.
Andreas Marek's avatar
Andreas Marek committed
378
379
380
   else
     do_test_analytic_eigenvalues_eigenvectors = .false.
   endif
381
#endif
Andreas Marek's avatar
Andreas Marek committed
382
383
   do_test_frank_eigenvalues = .false.
   do_test_toeplitz_eigenvalues = .false.
384
#endif /* TEST_MATRIX_ANALYTIC */
Andreas Marek's avatar
Andreas Marek committed
385
386
387
#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
388

Andreas Marek's avatar
Andreas Marek committed
389
390
#if defined(TEST_MATRIX_TOEPLITZ)
   ! The Toeplitz matrix works in each test
Andreas Marek's avatar
Andreas Marek committed
391
392
393
394
395
396
397
#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
398

399
! actually we test cholesky for diagonal matrix only
400
401
#if defined(TEST_CHOLESKY)
#ifdef TEST_SINGLE
402
403
  diagonalElement = (2.546_c_float, 0.0_c_float)
  subdiagonalElement =  (0.0_c_float, 0.0_c_float)
404
#else
405
406
  diagonalElement = (2.546_c_double, 0.0_c_double)
  subdiagonalElement =  (0.0_c_double, 0.0_c_double)
407
#endif
Andreas Marek's avatar
Andreas Marek committed
408
409
#endif /* TEST_CHOLESKY */

410
   call prepare_matrix_toeplitz(na, diagonalElement, subdiagonalElement, &
411
412
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
413
414


Andreas Marek's avatar
Andreas Marek committed
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
   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
432

Andreas Marek's avatar
Andreas Marek committed
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
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
#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 */
510

511
512
513
514
515
516
517
518
519
520
521
522
! 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

523
524
   e => elpa_allocate()

525
526
527
528
529
530
531
532
533
534
   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)
535
536

#ifdef WITH_MPI
537
538
539
540
541
542
   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)
543
544
545
546
#endif
#ifdef TEST_GENERALIZED_EIGENPROBLEM
   call e%set("blacs_context", my_blacs_ctxt, error)
   assert_elpa_ok(error)
547
#endif
Andreas Marek's avatar
Andreas Marek committed
548
   call e%set("timings",1,error)
549
550
551
   assert_elpa_ok(e%setup())

#ifdef TEST_SOLVER_1STAGE
Andreas Marek's avatar
Andreas Marek committed
552
   call e%set("solver", ELPA_SOLVER_1STAGE,error)
553
#else
Andreas Marek's avatar
Andreas Marek committed
554
   call e%set("solver", ELPA_SOLVER_2STAGE,error)
555
#endif
556
   assert_elpa_ok(error)
557

558
559
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
560

Andreas Marek's avatar
Andreas Marek committed
561
#if TEST_QR_DECOMPOSITION == 1
562
563
564
565
   call e%set("qr", 1, error)
   assert_elpa_ok(error)
#endif

566
567
   if (myid == 0) print *, ""

568
#ifdef TEST_ALL_KERNELS
569
   do i = 0, elpa_option_cardinality(KERNEL_KEY)  ! kernels
570
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
571
#endif
572
#ifdef TEST_KERNEL
573
     kernel = TEST_KERNEL
574
#endif
575
576

#ifdef TEST_SOLVER_2STAGE
577
     call e%set(KERNEL_KEY, kernel, error)
578
579
580
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
581
582
583
     if (error /= ELPA_OK) then
       cycle
     endif
584
     ! actually used kernel might be different if forced via environment variables
Andreas Marek's avatar
Andreas Marek committed
585
     call e%get(KERNEL_KEY, kernel, error)
586
587
588
589
#endif
     if (myid == 0) then
       print *, elpa_int_value_to_string(KERNEL_KEY, kernel) // " kernel"
     endif
590
591
#endif

592
#ifdef TEST_ALL_KERNELS
593
     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
594
#endif
595

596
     ! The actual solve step
Andreas Marek's avatar
Andreas Marek committed
597
598
599
600
#if defined(TEST_EIGENVECTORS)
#if TEST_QR_DECOMPOSITION == 1
     call e%timer_start("e%eigenvectors_qr()")
#else
601
     call e%timer_start("e%eigenvectors()")
Andreas Marek's avatar
Andreas Marek committed
602
#endif
Pavel Kus's avatar
Pavel Kus committed
603
604
#ifdef TEST_SCALAPACK_ALL
     call solve_scalapack_all(na, a, sc_desc, ev, z)
605
606
#elif TEST_SCALAPACK_PART
     call solve_scalapack_part(na, a, sc_desc, nev, ev, z)
607
     check_all_evals = .false. ! scalapack does not compute all eigenvectors
Pavel Kus's avatar
Pavel Kus committed
608
#else
609
     call e%eigenvectors(a, ev, z, error)
Pavel Kus's avatar
Pavel Kus committed
610
#endif
Andreas Marek's avatar
Andreas Marek committed
611
612
613
#if TEST_QR_DECOMPOSITION == 1
     call e%timer_stop("e%eigenvectors_qr()")
#else
614
     call e%timer_stop("e%eigenvectors()")
Andreas Marek's avatar
Andreas Marek committed
615
616
#endif
#endif /* TEST_EIGENVECTORS  */
617
618
619

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
620
     call e%eigenvalues(a, ev, error)
621
     call e%timer_stop("e%eigenvalues()")
622
#endif
623
624
625

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
626
     call e%solve_tridiagonal(d, sd, z, error)
627
     call e%timer_stop("e%solve_tridiagonal()")
628
629
630
     ev(:) = d(:)
#endif

631
632
633
#if defined(TEST_CHOLESKY)
     call e%timer_start("e%cholesky()")
     call e%cholesky(a, error)
634
     assert_elpa_ok(error)
635
636
637
     call e%timer_stop("e%cholesky()")
#endif

638
639
640
641
642
#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
643

644
645
#if defined(TEST_GENERALIZED_EIGENPROBLEM)
     call e%timer_start("e%generalized_eigenvectors()")
646
647
648
#if defined(TEST_GENERALIZED_DECOMP_EIGENPROBLEM)
     call e%timer_start("is_already_decomposed=.false.")
#endif
649
     call e%generalized_eigenvectors(a, b, ev, z, .false., error)
650
651
652
653
#if defined(TEST_GENERALIZED_DECOMP_EIGENPROBLEM)
     call e%timer_stop("is_already_decomposed=.false.")
     a = as
     call e%timer_start("is_already_decomposed=.true.")
654
     call e%generalized_eigenvectors(a, b, ev, z, .true., error)
655
656
     call e%timer_stop("is_already_decomposed=.true.")
#endif
657
658
659
     call e%timer_stop("e%generalized_eigenvectors()")
#endif

660
661
     assert_elpa_ok(error)

662
#ifdef TEST_ALL_KERNELS
663
664
665
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

666
     if (myid .eq. 0) then
667
#ifdef TEST_ALL_KERNELS
668
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
669
670
#else /* TEST_ALL_KERNELS */

Andreas Marek's avatar
Andreas Marek committed
671
672
673
674
#if defined(TEST_EIGENVECTORS)
#if TEST_QR_DECOMPOSITION == 1
       call e%print_times("e%eigenvectors_qr()")
#else
675
       call e%print_times("e%eigenvectors()")
676
#endif
Andreas Marek's avatar
Andreas Marek committed
677
#endif
678
#ifdef TEST_EIGENVALUES
679
680
       call e%print_times("e%eigenvalues()")
#endif
681
682
#ifdef TEST_SOLVE_TRIDIAGONAL
       call e%print_times("e%solve_tridiagonal()")
683
#endif
684
685
#ifdef TEST_CHOLESKY
       call e%print_times("e%cholesky()")
686
#endif
687
688
689
#ifdef TEST_HERMITIAN_MULTIPLY
       call e%print_times("e%hermitian_multiply()")
#endif
690
691
692
#ifdef TEST_GENERALIZED_EIGENPROBLEM
      call e%print_times("e%generalized_eigenvectors()")
#endif
693
#endif /* TEST_ALL_KERNELS */
694
     endif
695

Andreas Marek's avatar
Andreas Marek committed
696
697
698
699
700
701
702
703
704
     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
705

Andreas Marek's avatar
Andreas Marek committed
706
     if(do_test_numeric_residual) then
707
       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
708
709
710
711
712
713
       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)
714
     endif
715

Andreas Marek's avatar
Andreas Marek committed
716
     if (do_test_toeplitz_eigenvalues) then
717
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
718
       status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
Andreas Marek's avatar
Andreas Marek committed
719
         subdiagonalElement, ev, z, myid)
Andreas Marek's avatar
Andreas Marek committed
720
       call check_status(status, myid)
721
#endif
Andreas Marek's avatar
Andreas Marek committed
722
     endif
723

Andreas Marek's avatar
Andreas Marek committed
724
725
726
727
     if (do_test_cholesky) then
       status = check_correctness_cholesky(na, a, as, na_rows, sc_desc, myid )
       call check_status(status, myid)
     endif
728

Andreas Marek's avatar
Andreas Marek committed
729
730
731
732
733
#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
734
#endif
735

736
737
738
739
740
741
742
743
#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

744
745
746
747
     if (myid == 0) then
       print *, ""
     endif

748
749
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
750
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
751
752
753
     d = ds
     sd = sds
#endif
754
   end do ! kernels
755
#endif
Andreas Marek's avatar
Andreas Marek committed
756

757
758
759
760
761
762
   call elpa_deallocate(e)

   deallocate(a)
   deallocate(as)
   deallocate(z)
   deallocate(ev)
763
764
765
766
#ifdef TEST_HERMITIAN_MULTIPLY
   deallocate(b)
   deallocate(c)
#endif
767
#if defined(TEST_MATRIX_TOEPLITZ) || defined(TEST_MATRIX_FRANK)
768
769
770
   deallocate(d, ds)
   deallocate(sd, sds)
#endif
771
772
773
#if defined(TEST_GENERALIZED_EIGENPROBLEM)
  deallocate(b, bs)
#endif
774

775
776
777
778
779
780
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif
   call elpa_uninit()

781
782
783
784
785
786
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif
   call exit(status)

787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
   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

802
end program