test.F90 19.9 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

148
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || 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
170
171
172
173
174
175
176
177
178
179
180
181
182
   integer                     :: i_layout
#endif
   integer                     :: kernel
   character(len=1)            :: layout
#ifdef TEST_COMPLEX
   EV_TYPE                     :: norm, normmax
   MATRIX_TYPE, allocatable    :: tmp1(:,:), tmp2(:,:)
#ifdef TEST_DOUBLE
   MATRIX_TYPE, parameter      :: CONE = (1.0_c_double, 0.0_c_double), &
                                  CZERO = (0.0_c_double, 0.0_c_double)
   EV_TYPE :: pzlange, zlange
#else
   MATRIX_TYPE, parameter      :: CONE = (1.0_c_float, 0.0_c_float), &
                                  CZERO = (0.0_c_float, 0.0_c_float)
   EV_TYPE :: pclange, clange
#endif
#endif
183
   call read_input_parameters_traditional(na, nev, nblk, write_to_file)
184
   call setup_mpi(myid, nprocs)
185
186
187
188
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
     call redirect_stdout(myid)
189
#endif
190
#endif
191

192
   check_all_evals = .true.
193

194
195
196
197
198
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

199
200
201
202
203
   if (myid == 0) then
     print '((a,i0))', 'Program ' // TEST_CASE
     print *, ""
   endif

204
#ifdef TEST_ALL_LAYOUTS
205
   do i_layout = 1, size(layouts)               ! layouts
206
     layout = layouts(i_layout)
207
     do np_cols = 1, nprocs                     ! factors
208
209
210
211
212
       if (mod(nprocs,np_cols) /= 0 ) then
         cycle
       endif
#else
   layout = 'C'
213
214
215
   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
216
#endif
217
218

   np_rows = nprocs/np_cols
219
   assert(nprocs == np_rows * np_cols)
220

221
222
223
224
   if (myid == 0) then
     print '((a,i0))', 'Matrix size: ', na
     print '((a,i0))', 'Num eigenvectors: ', nev
     print '((a,i0))', 'Blocksize: ', nblk
225
#ifdef WITH_MPI
226
227
     print '((a,i0))', 'Num MPI proc: ', nprocs
     print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
228
     print '(a)',      'Process layout: ' // layout
229
#endif
230
231
232
     print *,''
   endif

233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
#ifdef TEST_QR_DECOMPOSITION

#if TEST_GPU == 1
#ifdef WITH_MPI
     call mpi_finalize(mpierr)
#endif
     stop 77
#endif
   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
#endif

256
257
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
258
259
260
261

   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
262
263
   allocate(a (na_rows,na_cols))
   allocate(as(na_rows,na_cols))
264
265
266
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

267
268
269
270
271
#ifdef TEST_HERMITIAN_MULTIPLY
   allocate(b (na_rows,na_cols))
   allocate(c (na_rows,na_cols))
#endif

272
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || defined(TEST_CHOLESKY)
273
274
275
276
277
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
   allocate(ev_analytic(na))
#endif

278
279
280
281
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

282
#if defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY) || defined(TEST_QR_DECOMPOSITION)
Pavel Kus's avatar
Pavel Kus committed
283
284
#ifdef TEST_MATRIX_ANALYTIC
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
285
   as(:,:) = a
286
287
#else
   if (nev .ge. 1) then
288
     call prepare_matrix_random(na, myid, sc_desc, a, z, as)
289
290
291
292
293
   else
     ! zero eigenvectors and not analytic test => toeplitz matrix
#ifdef TEST_SINGLE
     diagonalElement = 0.45_c_float
     subdiagonalElement =  0.78_c_float
Pavel Kus's avatar
Pavel Kus committed
294
#else
295
296
297
     diagonalElement = 0.45_c_double
     subdiagonalElement =  0.78_c_double
#endif
298
     call prepare_matrix_toeplitz(na, diagonalElement, subdiagonalElement, &
299
300
301
                                  d, sd, ds, sds, a, as, nblk, np_rows, &
                                  np_cols, my_prow, my_pcol)
   endif
302
303

#ifdef TEST_HERMITIAN_MULTIPLY
304
#ifdef TEST_REAL
305

306
307
308
#ifdef TEST_DOUBLE
   b(:,:) = 2.0_c_double * a(:,:)
   c(:,:) = 0.0_c_double
309
#else
310
311
   b(:,:) = 2.0_c_float * a(:,:)
   c(:,:) = 0.0_c_float
312
#endif
313

Pavel Kus's avatar
Pavel Kus committed
314
#endif
315

316
#ifdef TEST_COMPLEX
317

318
319
320
#ifdef TEST_DOUBLE
   b(:,:) = 2.0_c_double * a(:,:)
   c(:,:) = (0.0_c_double, 0.0_c_double)
321
#else
322
323
   b(:,:) = 2.0_c_float * a(:,:)
   c(:,:) = (0.0_c_float, 0.0_c_float)
324
325
326
327
328
329
330
#endif

#endif

#endif /* TEST_HERMITIAN_MULTIPLY */

#endif /* TEST_MATRIX_ANALYTIC */
331
#endif /* defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY) || defined(TEST_QR_DECOMPOSITION) */
332

333
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
334
335
336
337
338
339
340
341

#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
342
   call prepare_matrix_toeplitz(na, diagonalElement, subdiagonalElement, &
Andreas Marek's avatar
Andreas Marek committed
343
344
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
345
#endif /* EIGENVALUES OR TRIDIAGONAL */
346

347
348
349
#if defined(TEST_CHOLESKY)

#ifdef TEST_SINGLE
350
351
   diagonalElement = (2.546_c_float, 0.0_c_float)
   subdiagonalElement =  (0.0_c_float, 0.0_c_float)
352
#else
353
354
   diagonalElement = (2.546_c_double, 0.0_c_double)
   subdiagonalElement =  (0.0_c_double, 0.0_c_double)
355
#endif
356
   call prepare_matrix_toeplitz(na, diagonalElement, subdiagonalElement, &
357
358
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
359
360


361
362
#endif /* TEST_CHOLESKY */

363
364
   e => elpa_allocate()

365
366
367
368
369
370
371
372
373
374
   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)
375
376

#ifdef WITH_MPI
377
378
379
380
381
382
   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)
383
#endif
384

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

387
388
389
390
391
392
393
   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
394
   assert_elpa_ok(error)
395

396
397
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
398

399
400
401
402
403
#ifdef TEST_QR_DECOMPOSITION
   call e%set("qr", 1, error)
   assert_elpa_ok(error)
#endif

404
405
   if (myid == 0) print *, ""

406
#ifdef TEST_ALL_KERNELS
407
   do i = 0, elpa_option_cardinality(KERNEL_KEY)  ! kernels
408
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
409
#endif
410
#ifdef TEST_KERNEL
411
     kernel = TEST_KERNEL
412
#endif
413
414

#ifdef TEST_SOLVER_2STAGE
415
     call e%set(KERNEL_KEY, kernel, error)
416
417
418
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
419
420
421
     if (error /= ELPA_OK) then
       cycle
     endif
422
423
424
425
426
427
     ! actually used kernel might be different if forced via environment variables
     call e%get(KERNEL_KEY, kernel)
#endif
     if (myid == 0) then
       print *, elpa_int_value_to_string(KERNEL_KEY, kernel) // " kernel"
     endif
428
429
#endif

430
#ifdef TEST_ALL_KERNELS
431
     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
432
#endif
433

434
     ! The actual solve step
435
#if defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION)
436
     call e%timer_start("e%eigenvectors()")
Pavel Kus's avatar
Pavel Kus committed
437
438
#ifdef TEST_SCALAPACK_ALL
     call solve_scalapack_all(na, a, sc_desc, ev, z)
439
440
#elif TEST_SCALAPACK_PART
     call solve_scalapack_part(na, a, sc_desc, nev, ev, z)
441
     check_all_evals = .false. ! scalapack does not compute all eigenvectors
Pavel Kus's avatar
Pavel Kus committed
442
#else
443
     call e%eigenvectors(a, ev, z, error)
Pavel Kus's avatar
Pavel Kus committed
444
#endif
445
     call e%timer_stop("e%eigenvectors()")
446
#endif /* TEST_EIGENVECTORS || defined(TEST_QR_DECOMPOSITION) */
447
448
449

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
450
     call e%eigenvalues(a, ev, error)
451
     call e%timer_stop("e%eigenvalues()")
452
#endif
453
454
455

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
456
     call e%solve_tridiagonal(d, sd, z, error)
457
     call e%timer_stop("e%solve_tridiagonal()")
458
459
460
     ev(:) = d(:)
#endif

461
462
463
464
465
466
#if defined(TEST_CHOLESKY)
     call e%timer_start("e%cholesky()")
     call e%cholesky(a, error)
     call e%timer_stop("e%cholesky()")
#endif

467
468
469
470
471
#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
472

473
474
     assert_elpa_ok(error)

475
#ifdef TEST_ALL_KERNELS
476
477
478
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

479
     if (myid .eq. 0) then
480
#ifdef TEST_ALL_KERNELS
481
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
482
483
#else /* TEST_ALL_KERNELS */

484
#if defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION)
485
       call e%print_times("e%eigenvectors()")
486
#endif
487
#ifdef TEST_EIGENVALUES
488
489
       call e%print_times("e%eigenvalues()")
#endif
490
491
#ifdef TEST_SOLVE_TRIDIAGONAL
       call e%print_times("e%solve_tridiagonal()")
492
#endif
493
494
#ifdef TEST_CHOLESKY
       call e%print_times("e%cholesky()")
495
#endif
496
497
498
#ifdef TEST_HERMITIAN_MULTIPLY
       call e%print_times("e%hermitian_multiply()")
#endif
499
#endif /* TEST_ALL_KERNELS */
500
     endif
501

502
#if defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION)
Pavel Kus's avatar
Pavel Kus committed
503
#ifdef TEST_MATRIX_ANALYTIC
504
505
!
!#if defined(TEST_MATRIX_ANALYTIC)
506
     status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol, check_all_evals)
507
508
509
510
511
512
     call check_status(status, myid)
     if (.true.) then
       ! also check residuals
       status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, my_pcol)
       call check_status(status, myid)
     endif
Pavel Kus's avatar
Pavel Kus committed
513
#else
514
515
516
!#elif defined(TEST_MATRIX_FRANK)
!     status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, my_pcol)
!#elif defined(TEST_MATRIX_RANDOM)
517
     if (nev .ge. 1) then
518
       status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, my_pcol)
519
520
521
522
523
     else
       ! zero eigenvectors and no analytic test => toeplitz
       status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
     endif
524
     call check_status(status, myid)
525
526
527
!#else
!#error "MATRIX TYPE"
!#endif
528
#endif
529
#endif /* defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) */
530

531
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
532
533
     status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
534
     call check_status(status, myid)
535

536
#ifdef TEST_SOLVE_TRIDIAGONAL
537
     ! check eigenvectors
538
     status = check_correctness_evp_numeric_residuals(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
539
     call check_status(status, myid)
540
541
#endif
#endif
542

543
544
545
546
547
#if defined(TEST_CHOLESKY)
     status = check_correctness_cholesky(na, a, as, na_rows, sc_desc, myid )
     call check_status(status, myid)
#endif

548
#if defined(TEST_HERMITIAN_MULTIPLY)
549
550

#ifdef TEST_REAL
551
552
553
     status = check_correctness_hermitian_multiply(na, a, b, c, na_rows, sc_desc, myid )
     call check_status(status, myid)
#endif
554
555
556
557
558
559
560
561
#ifdef TEST_COMPLEX
   status = 0

   !-------------------------------------------------------------------------------
   ! Test correctness of result (using plain scalapack routines)
   allocate(tmp1(na_rows,na_cols))
   allocate(tmp2(na_rows,na_cols))
#ifdef TEST_DOUBLE
562
   tmp1(:,:) = (0.0_c_double, 0.0_c_double)
563
#else
564
   tmp1(:,:) = (0.0_c_float, 0.0_c_float)
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
#endif
   ! tmp1 = a**T
#ifdef WITH_MPI
#ifdef TEST_DOUBLE
   call pztranc(na, na, CONE, a, 1, 1, sc_desc, CZERO, tmp1, 1, 1, sc_desc)
#else
   call pctranc(na, na, CONE, a, 1, 1, sc_desc, CZERO, tmp1, 1, 1, sc_desc)
#endif
#else
   tmp1 = transpose(conjg(a))
#endif
   ! tmp2 = tmp1 * b
#ifdef TEST_DOUBLE
#ifdef WITH_MPI
   call pzgemm("N","N", na, na, na, CONE, tmp1, 1, 1, sc_desc, b, 1, 1, &
               sc_desc, CZERO, tmp2, 1, 1, sc_desc)
#else
   call zgemm("N","N", na, na, na, CONE, tmp1, na, b, na, CZERO, tmp2, na)
#endif
#else
#ifdef WITH_MPI
   call pcgemm("N","N", na, na, na, CONE, tmp1, 1, 1, sc_desc, b, 1, 1, &
               sc_desc, CZERO, tmp2, 1, 1, sc_desc)
#else
   call cgemm("N","N", na, na, na, CONE, tmp1, na, b, na, CZERO, tmp2, na)
#endif
#endif

   ! compare tmp2 with c
   tmp2(:,:) = tmp2(:,:) - c(:,:)
#ifdef TEST_DOUBLE
#ifdef WITH_MPI
   norm = pzlange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
#else
   norm = zlange("M",na, na, tmp2, na_rows, tmp1)
#endif
#else
#ifdef WITH_MPI
   norm = pclange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
#else
   norm = clange("M",na, na, tmp2, na_rows, tmp1)
#endif
#endif
#ifdef WITH_MPI
#ifdef TEST_DOUBLE
   call mpi_allreduce(norm,normmax,1,MPI_REAL8,MPI_MAX,MPI_COMM_WORLD,mpierr)
#else
   call mpi_allreduce(norm,normmax,1,MPI_REAL4,MPI_MAX,MPI_COMM_WORLD,mpierr)
#endif
#else
   normmax = norm
#endif
   if (myid .eq. 0) then
     print *," Maximum error of result: ", normmax
   endif

#ifdef TEST_DOUBLE
Andreas Marek's avatar
Andreas Marek committed
622
   if (normmax .gt. 5e-11_rk8) then
623
#else
Andreas Marek's avatar
Andreas Marek committed
624
   if (normmax .gt. 5e-3_rk4) then
625
626
627
628
629
630
631
632
633
634
#endif
        print *,"norm= ",normmax
        status = 1
   endif

   deallocate(tmp1)
   deallocate(tmp2)

#endif
#endif /* TEST_HERMITIAN_MULTIPLY */
635
636


637
638
639
640
     if (myid == 0) then
       print *, ""
     endif

641
642
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
643
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || defined(TEST_CHOLESKY)
644
645
646
     d = ds
     sd = sds
#endif
647
   end do ! kernels
648
#endif
Andreas Marek's avatar
Andreas Marek committed
649

650
651
652
653
654
655
656
   call elpa_deallocate(e)

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

657
658
659
660
661
#ifdef TEST_HERMITIAN_MULTIPLY
   deallocate(b)
   deallocate(c)
#endif

662
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || defined(TEST_CHOLESKY)
663
664
665
666
667
   deallocate(d, ds)
   deallocate(sd, sds)
   deallocate(ev_analytic)
#endif

668
669
670
671
672
673
674
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif

   call elpa_uninit()

675
676
677
678
679
680
681
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif

   call exit(status)

682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
   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

697
end program