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

Pavel Kus's avatar
Pavel Kus committed
59
#if !(defined(TEST_SOLVER_1STAGE) ^ defined(TEST_SOLVER_2STAGE) ^ defined(TEST_SCALAPACK_ALL))
60
error: define exactly one of TEST_SOLVER_1STAGE or TEST_SOLVER_2STAGE or TEST_SCALAPACK_ALL
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
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || defined(TEST_HERMITIAN_MULTIPLY)
147
148
   EV_TYPE, allocatable        :: d(:), sd(:), ds(:), sds(:)
   EV_TYPE                     :: diagonalELement, subdiagonalElement
149
#endif
150
#if defined(TEST_CHOLESKY)
151
152
   MATRIX_TYPE, allocatable    :: d(:), sd(:), ds(:), sds(:)
   MATRIX_TYPE                 :: diagonalELement, subdiagonalElement
153
#endif
154

155
   integer                     :: error, status
156

157
158
   type(output_t)              :: write_to_file
   class(elpa_t), pointer      :: e
159
#ifdef TEST_ALL_KERNELS
160
   integer                     :: i
161
162
163
#endif
#ifdef TEST_ALL_LAYOUTS
   character(len=1), parameter :: layouts(2) = [ 'C', 'R' ]
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
   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
181
   call read_input_parameters_traditional(na, nev, nblk, write_to_file)
182
   call setup_mpi(myid, nprocs)
183
184
185
186
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
     call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
     call redirect_stdout(myid)
187
#endif
188
#endif
189

190
191


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

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

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

   np_rows = nprocs/np_cols
217
   assert(nprocs == np_rows * np_cols)
218

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

231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
#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

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

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

265
266
267
268
269
#ifdef TEST_HERMITIAN_MULTIPLY
   allocate(b (na_rows,na_cols))
   allocate(c (na_rows,na_cols))
#endif

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

276
277
278
279
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

280
#if defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY) || defined(TEST_QR_DECOMPOSITION)
Pavel Kus's avatar
Pavel Kus committed
281
282
#ifdef TEST_MATRIX_ANALYTIC
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
283
   as(:,:) = a
284
285
286
287
288
289
290
291
#else
   if (nev .ge. 1) then
     call prepare_matrix(na, myid, sc_desc, a, z, as)
   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
292
#else
293
294
295
296
297
298
299
     diagonalElement = 0.45_c_double
     subdiagonalElement =  0.78_c_double
#endif
     call prepare_toeplitz_matrix(na, diagonalElement, subdiagonalElement, &
                                  d, sd, ds, sds, a, as, nblk, np_rows, &
                                  np_cols, my_prow, my_pcol)
   endif
300
301

#ifdef TEST_HERMITIAN_MULTIPLY
302
#ifdef TEST_REAL
303

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

Pavel Kus's avatar
Pavel Kus committed
312
#endif
313

314
#ifdef TEST_COMPLEX
315

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

#endif

#endif /* TEST_HERMITIAN_MULTIPLY */

#endif /* TEST_MATRIX_ANALYTIC */
329
#endif /* defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY) || defined(TEST_QR_DECOMPOSITION) */
330

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

#ifdef TEST_SINGLE
   diagonalElement = 0.45_c_float
   subdiagonalElement =  0.78_c_float
#else
   diagonalElement = 0.45_c_double
   subdiagonalElement =  0.78_c_double
#endif
   call prepare_toeplitz_matrix(na, diagonalElement, subdiagonalElement, &
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
343
#endif /* EIGENVALUES OR TRIDIAGONAL */
344

345
346
347
#if defined(TEST_CHOLESKY)

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


359
360
#endif /* TEST_CHOLESKY */

361
362
   e => elpa_allocate()

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

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

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

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

394
395
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
396

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

402
403
   if (myid == 0) print *, ""

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

#ifdef TEST_SOLVER_2STAGE
413
     call e%set(KERNEL_KEY, kernel, error)
414
415
416
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
417
418
419
     if (error /= ELPA_OK) then
       cycle
     endif
420
421
422
423
424
425
     ! 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
426
427
#endif

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

432
     ! The actual solve step
433
#if defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION)
434
     call e%timer_start("e%eigenvectors()")
Pavel Kus's avatar
Pavel Kus committed
435
436
437
#ifdef TEST_SCALAPACK_ALL
     call solve_scalapack_all(na, a, sc_desc, ev, z)
#else
438
     call e%eigenvectors(a, ev, z, error)
Pavel Kus's avatar
Pavel Kus committed
439
#endif
440
     call e%timer_stop("e%eigenvectors()")
441
#endif /* TEST_EIGENVECTORS || defined(TEST_QR_DECOMPOSITION) */
442
443
444

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
445
     call e%eigenvalues(a, ev, error)
446
     call e%timer_stop("e%eigenvalues()")
447
#endif
448
449
450

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
451
     call e%solve_tridiagonal(d, sd, z, error)
452
     call e%timer_stop("e%solve_tridiagonal()")
453
454
455
     ev(:) = d(:)
#endif

456
457
458
459
460
461
#if defined(TEST_CHOLESKY)
     call e%timer_start("e%cholesky()")
     call e%cholesky(a, error)
     call e%timer_stop("e%cholesky()")
#endif

462
463
464
465
466
#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
467

468
469
     assert_elpa_ok(error)

470
#ifdef TEST_ALL_KERNELS
471
472
473
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

474
     if (myid .eq. 0) then
475
#ifdef TEST_ALL_KERNELS
476
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
477
478
#else /* TEST_ALL_KERNELS */

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

497
#if defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION)
Pavel Kus's avatar
Pavel Kus committed
498
499
500
#ifdef TEST_MATRIX_ANALYTIC
     status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
#else
501
502
503
504
505
506
507
     if (nev .ge. 1) then
       status = check_correctness(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows,np_cols, my_prow, my_pcol)
     else
       ! zero eigenvectors and no analytic test => toeplitz
       status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
     endif
Andreas Marek's avatar
Andreas Marek committed
508
#endif
509
     call check_status(status, myid)
510
511
#endif

512
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
513
514
     status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
515
     call check_status(status, myid)
516

517
#ifdef TEST_SOLVE_TRIDIAGONAL
518
     ! check eigenvectors
Pavel Kus's avatar
Pavel Kus committed
519
     status = check_correctness(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
520
     call check_status(status, myid)
521
522
#endif
#endif
523

524
525
526
527
528
#if defined(TEST_CHOLESKY)
     status = check_correctness_cholesky(na, a, as, na_rows, sc_desc, myid )
     call check_status(status, myid)
#endif

529
#if defined(TEST_HERMITIAN_MULTIPLY)
530
531

#ifdef TEST_REAL
532
533
534
     status = check_correctness_hermitian_multiply(na, a, b, c, na_rows, sc_desc, myid )
     call check_status(status, myid)
#endif
535
536
537
538
539
540
541
542
#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
543
   tmp1(:,:) = (0.0_c_double, 0.0_c_double)
544
#else
545
   tmp1(:,:) = (0.0_c_float, 0.0_c_float)
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
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
#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
603
   if (normmax .gt. 5e-11_rk8) then
604
#else
Andreas Marek's avatar
Andreas Marek committed
605
   if (normmax .gt. 5e-3_rk4) then
606
607
608
609
610
611
612
613
614
615
#endif
        print *,"norm= ",normmax
        status = 1
   endif

   deallocate(tmp1)
   deallocate(tmp2)

#endif
#endif /* TEST_HERMITIAN_MULTIPLY */
616
617


618
619
620
621
     if (myid == 0) then
       print *, ""
     endif

622
623
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
624
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || defined(TEST_CHOLESKY)
625
626
627
     d = ds
     sd = sds
#endif
628
   end do ! kernels
629
#endif
Andreas Marek's avatar
Andreas Marek committed
630

631
632
633
634
635
636
637
   call elpa_deallocate(e)

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

638
639
640
641
642
#ifdef TEST_HERMITIAN_MULTIPLY
   deallocate(b)
   deallocate(c)
#endif

643
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_QR_DECOMPOSITION) || defined(TEST_CHOLESKY)
644
645
646
647
648
   deallocate(d, ds)
   deallocate(sd, sds)
   deallocate(ev_analytic)
#endif

649
650
651
652
653
654
655
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif

   call elpa_uninit()

656
657
658
659
660
661
662
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif

   call exit(status)

663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
   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

678
end program