test.F90 21.7 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 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
#ifdef WITH_SCALAPACK_TESTS
Pavel Kus's avatar
Pavel Kus committed
115
   use test_scalapack
116
#endif
117

118 119 120
#ifdef HAVE_REDIRECT
   use test_redirect
#endif
121 122 123
   implicit none

   ! matrix dimensions
124
   integer                     :: na, nev, nblk
125 126

   ! mpi
127 128 129 130 131
   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
132 133

   ! blacs
134
   integer                     :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol
135 136

   ! The Matrix
137
   MATRIX_TYPE, allocatable    :: a(:,:), as(:,:)
138
#if defined(TEST_HERMITIAN_MULTIPLY)
139
   MATRIX_TYPE, allocatable    :: b(:,:), c(:,:)
140
#endif
141
   ! eigenvectors
142
   MATRIX_TYPE, allocatable    :: z(:,:)
143
   ! eigenvalues
144
   EV_TYPE, allocatable        :: ev(:), ev_analytic(:)
145

146 147
   logical                     :: check_all_evals

Andreas Marek's avatar
Andreas Marek committed
148
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || TEST_QR_DECOMPOSITION == 1 || defined(TEST_HERMITIAN_MULTIPLY)
149 150
   EV_TYPE, allocatable        :: d(:), sd(:), ds(:), sds(:)
   EV_TYPE                     :: diagonalELement, subdiagonalElement
151
#endif
152
#if defined(TEST_CHOLESKY)
153 154
   MATRIX_TYPE, allocatable    :: d(:), sd(:), ds(:), sds(:)
   MATRIX_TYPE                 :: diagonalELement, subdiagonalElement
155
#endif
156

157
   integer                     :: error, status
158

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

176
   call read_input_parameters_traditional(na, nev, nblk, write_to_file)
177
   call setup_mpi(myid, nprocs)
178 179 180 181
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
     call redirect_stdout(myid)
182
#endif
183
#endif
184

185
   check_all_evals = .true.
186

Andreas Marek's avatar
Andreas Marek committed
187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203

   do_test_numeric_residual = .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.
#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
204 205 206 207 208
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

209 210 211 212 213
   if (myid == 0) then
     print '((a,i0))', 'Program ' // TEST_CASE
     print *, ""
   endif

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

   np_rows = nprocs/np_cols
229
   assert(nprocs == np_rows * np_cols)
230

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

Andreas Marek's avatar
Andreas Marek committed
243
#if TEST_QR_DECOMPOSITION == 1
244 245 246 247 248 249

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

266

267 268
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
269 270 271 272

   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
273 274
   allocate(a (na_rows,na_cols))
   allocate(as(na_rows,na_cols))
275 276 277
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

278 279 280 281 282
#ifdef TEST_HERMITIAN_MULTIPLY
   allocate(b (na_rows,na_cols))
   allocate(c (na_rows,na_cols))
#endif

Andreas Marek's avatar
Andreas Marek committed
283
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || TEST_QR_DECOMPOSITION == 1|| defined(TEST_CHOLESKY)
284 285 286 287 288
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
   allocate(ev_analytic(na))
#endif

289 290 291 292
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

Andreas Marek's avatar
Andreas Marek committed
293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312
#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
   ! RANDOM + TEST_CHOLESKY: no correctness check yet implemented
   ! 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)

    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
    do_test_frank_eigenvalues = .false.
    do_test_toeplitz_eigenvalues = .false.

#else /* TEST_EIGENVECTORS */

    if (nev .ge. 1) then
313
     call prepare_matrix_random(na, myid, sc_desc, a, z, as)
Andreas Marek's avatar
Andreas Marek committed
314 315 316 317 318 319 320 321

    do_test_analytic_eigenvalues = .false.
    do_test_analytic_eigenvalues_eigenvectors = .false.
    do_test_frank_eigenvalues = .false.
    do_test_toeplitz_eigenvalues = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
    do_test_numeric_residual = .true.
#endif
322
   else
Andreas Marek's avatar
Andreas Marek committed
323 324 325 326 327
     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)
328
#endif
Andreas Marek's avatar
Andreas Marek committed
329
     stop 77
330

Andreas Marek's avatar
Andreas Marek committed
331
   endif
332

Andreas Marek's avatar
Andreas Marek committed
333 334 335 336
#endif /* TEST_EIGENVECTORS */
#endif /* (TEST_MATRIX_RANDOM) */
#if defined(TEST_MATRIX_RANDOM) && (defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_CHOLESKY) || defined(TEST_EIGENVALUES))
#error "Random matrix is not allowed in this configuration"
337
#endif
338

Andreas Marek's avatar
Andreas Marek committed
339 340 341 342 343 344
#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
345

Andreas Marek's avatar
Andreas Marek committed
346 347
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
   as(:,:) = a
348

Andreas Marek's avatar
Andreas Marek committed
349 350 351 352
   do_test_numeric_residual = .false.
   do_test_analytic_eigenvalues_eigenvectors = .false.
#ifndef TEST_HERMITIAN_MULTIPLY
   do_test_analytic_eigenvalues = .true.
353
#endif
Andreas Marek's avatar
Andreas Marek committed
354 355 356 357 358 359 360 361
#if defined(TEST_EIGENVECTORS)
   if (nev .ge. 1) then
     do_test_analytic_eigenvalues_eigenvectors = .true.
     do_test_numeric_residual = .true.
   else
     do_test_analytic_eigenvalues_eigenvectors = .false.
     do_test_numeric_residual = .false.
   endif
362
#endif
Andreas Marek's avatar
Andreas Marek committed
363 364
   do_test_frank_eigenvalues = .false.
   do_test_toeplitz_eigenvalues = .false.
365
#endif /* TEST_MATRIX_ANALYTIC */
Andreas Marek's avatar
Andreas Marek committed
366 367 368
#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
369

Andreas Marek's avatar
Andreas Marek committed
370 371
#if defined(TEST_MATRIX_TOEPLITZ)
   ! The Toeplitz matrix works in each test
Andreas Marek's avatar
Andreas Marek committed
372 373 374 375 376 377 378
#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
379

380 381
#if defined(TEST_CHOLESKY)
#ifdef TEST_SINGLE
382 383
   diagonalElement = (2.546_c_float, 0.0_c_float)
   subdiagonalElement =  (0.0_c_float, 0.0_c_float)
384
#else
385 386
   diagonalElement = (2.546_c_double, 0.0_c_double)
   subdiagonalElement =  (0.0_c_double, 0.0_c_double)
387
#endif
Andreas Marek's avatar
Andreas Marek committed
388 389
#endif /* TEST_CHOLESKY */

390
   call prepare_matrix_toeplitz(na, diagonalElement, subdiagonalElement, &
391 392
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
393 394


Andreas Marek's avatar
Andreas Marek committed
395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 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 486 487 488
   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
#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 */
489

490 491
   e => elpa_allocate()

492 493 494 495 496 497 498 499 500 501
   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)
502 503

#ifdef WITH_MPI
504 505 506 507 508 509
   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)
510
#endif
Andreas Marek's avatar
Andreas Marek committed
511
   call e%set("timings",1,error)
512 513 514
   assert_elpa_ok(e%setup())

#ifdef TEST_SOLVER_1STAGE
Andreas Marek's avatar
Andreas Marek committed
515
   call e%set("solver", ELPA_SOLVER_1STAGE,error)
516
#else
Andreas Marek's avatar
Andreas Marek committed
517
   call e%set("solver", ELPA_SOLVER_2STAGE,error)
518
#endif
519
   assert_elpa_ok(error)
520

521 522
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
523

Andreas Marek's avatar
Andreas Marek committed
524
#if TEST_QR_DECOMPOSITION == 1
525 526 527 528
   call e%set("qr", 1, error)
   assert_elpa_ok(error)
#endif

529 530
   if (myid == 0) print *, ""

531
#ifdef TEST_ALL_KERNELS
532
   do i = 0, elpa_option_cardinality(KERNEL_KEY)  ! kernels
533
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
534
#endif
535
#ifdef TEST_KERNEL
536
     kernel = TEST_KERNEL
537
#endif
538 539

#ifdef TEST_SOLVER_2STAGE
540
     call e%set(KERNEL_KEY, kernel, error)
541 542 543
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
544 545 546
     if (error /= ELPA_OK) then
       cycle
     endif
547
     ! actually used kernel might be different if forced via environment variables
Andreas Marek's avatar
Andreas Marek committed
548
     call e%get(KERNEL_KEY, kernel, error)
549 550 551 552
#endif
     if (myid == 0) then
       print *, elpa_int_value_to_string(KERNEL_KEY, kernel) // " kernel"
     endif
553 554
#endif

555
#ifdef TEST_ALL_KERNELS
556
     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
557
#endif
558

559
     ! The actual solve step
Andreas Marek's avatar
Andreas Marek committed
560 561 562 563
#if defined(TEST_EIGENVECTORS)
#if TEST_QR_DECOMPOSITION == 1
     call e%timer_start("e%eigenvectors_qr()")
#else
564
     call e%timer_start("e%eigenvectors()")
Andreas Marek's avatar
Andreas Marek committed
565
#endif
Pavel Kus's avatar
Pavel Kus committed
566 567
#ifdef TEST_SCALAPACK_ALL
     call solve_scalapack_all(na, a, sc_desc, ev, z)
568 569
#elif TEST_SCALAPACK_PART
     call solve_scalapack_part(na, a, sc_desc, nev, ev, z)
570
     check_all_evals = .false. ! scalapack does not compute all eigenvectors
Pavel Kus's avatar
Pavel Kus committed
571
#else
572
     call e%eigenvectors(a, ev, z, error)
Pavel Kus's avatar
Pavel Kus committed
573
#endif
Andreas Marek's avatar
Andreas Marek committed
574 575 576
#if TEST_QR_DECOMPOSITION == 1
     call e%timer_stop("e%eigenvectors_qr()")
#else
577
     call e%timer_stop("e%eigenvectors()")
Andreas Marek's avatar
Andreas Marek committed
578 579
#endif
#endif /* TEST_EIGENVECTORS  */
580 581 582

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
583
     call e%eigenvalues(a, ev, error)
584
     call e%timer_stop("e%eigenvalues()")
585
#endif
586 587 588

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
589
     call e%solve_tridiagonal(d, sd, z, error)
590
     call e%timer_stop("e%solve_tridiagonal()")
591 592 593
     ev(:) = d(:)
#endif

594 595 596 597 598 599
#if defined(TEST_CHOLESKY)
     call e%timer_start("e%cholesky()")
     call e%cholesky(a, error)
     call e%timer_stop("e%cholesky()")
#endif

600 601 602 603 604
#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
605

606 607
     assert_elpa_ok(error)

608
#ifdef TEST_ALL_KERNELS
609 610 611
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

612
     if (myid .eq. 0) then
613
#ifdef TEST_ALL_KERNELS
614
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
615 616
#else /* TEST_ALL_KERNELS */

Andreas Marek's avatar
Andreas Marek committed
617 618 619 620
#if defined(TEST_EIGENVECTORS)
#if TEST_QR_DECOMPOSITION == 1
       call e%print_times("e%eigenvectors_qr()")
#else
621
       call e%print_times("e%eigenvectors()")
622
#endif
Andreas Marek's avatar
Andreas Marek committed
623
#endif
624
#ifdef TEST_EIGENVALUES
625 626
       call e%print_times("e%eigenvalues()")
#endif
627 628
#ifdef TEST_SOLVE_TRIDIAGONAL
       call e%print_times("e%solve_tridiagonal()")
629
#endif
630 631
#ifdef TEST_CHOLESKY
       call e%print_times("e%cholesky()")
632
#endif
633 634 635
#ifdef TEST_HERMITIAN_MULTIPLY
       call e%print_times("e%hermitian_multiply()")
#endif
636
#endif /* TEST_ALL_KERNELS */
637
     endif
638

Andreas Marek's avatar
Andreas Marek committed
639 640 641 642 643 644 645 646 647 648
     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
     if(do_test_numeric_residual) then
649
       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
650 651 652 653 654 655
       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)
656
     endif
657

Andreas Marek's avatar
Andreas Marek committed
658
     if (do_test_toeplitz_eigenvalues) then
659
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
660
       status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
Andreas Marek's avatar
Andreas Marek committed
661
         subdiagonalElement, ev, z, myid)
Andreas Marek's avatar
Andreas Marek committed
662
       call check_status(status, myid)
663
#endif
Andreas Marek's avatar
Andreas Marek committed
664
     endif
665

Andreas Marek's avatar
Andreas Marek committed
666 667 668 669
     if (do_test_cholesky) then
       status = check_correctness_cholesky(na, a, as, na_rows, sc_desc, myid )
       call check_status(status, myid)
     endif
670

Andreas Marek's avatar
Andreas Marek committed
671 672 673 674 675
#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
676
#endif
677

678 679 680 681
     if (myid == 0) then
       print *, ""
     endif

682 683
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
Andreas Marek's avatar
Andreas Marek committed
684
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || TEST_QR_DECOMPOSITION == 1 || defined(TEST_CHOLESKY)
685 686 687
     d = ds
     sd = sds
#endif
688
   end do ! kernels
689
#endif
Andreas Marek's avatar
Andreas Marek committed
690

691 692 693 694 695 696
   call elpa_deallocate(e)

   deallocate(a)
   deallocate(as)
   deallocate(z)
   deallocate(ev)
697 698 699 700
#ifdef TEST_HERMITIAN_MULTIPLY
   deallocate(b)
   deallocate(c)
#endif
Andreas Marek's avatar
Andreas Marek committed
701
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || TEST_QR_DECOMPOSITION == 1 || defined(TEST_CHOLESKY)
702 703 704 705 706
   deallocate(d, ds)
   deallocate(sd, sds)
   deallocate(ev_analytic)
#endif

707 708 709 710 711 712
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif
   call elpa_uninit()

713 714 715 716 717 718
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif
   call exit(status)

719 720 721 722 723 724 725 726 727 728 729 730 731 732 733
   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

734
end program