test.F90 16.6 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
124
125
126
127
128
129
130
131
132
133
   implicit none

   ! matrix dimensions
   integer :: na, nev, nblk

   ! mpi
   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

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

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

#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY)
   EV_TYPE, allocatable :: d(:), sd(:), ds(:), sds(:)
Andreas Marek's avatar
Andreas Marek committed
148
   EV_TYPE              :: diagonalELement, subdiagonalElement
149
#endif
150
151
152
153
#if defined(TEST_CHOLESKY)
   MATRIX_TYPE, allocatable :: d(:), sd(:), ds(:), sds(:)
   MATRIX_TYPE              :: diagonalELement, subdiagonalElement
#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
164
#endif
#ifdef TEST_ALL_LAYOUTS
   character(len=1), parameter :: layouts(2) = [ 'C', 'R' ]
   integer :: i_layout
165
#endif
166
   integer :: kernel
167
   character(len=1) :: layout
Andreas Marek's avatar
Andreas Marek committed
168
169
170
171
172
173
174
175
176
177
!#ifdef TEST_COMPLEX
!   MATRIX_TYPE, allocatable :: tmp1(:,:), tmp2(:,:)
!   EV_TYPE :: norm, normmax
!#ifdef TEST_SINGLE
!   EV_TYPE :: pclange
!#else
!   EV_TYPE :: pzlange
!#endif
!   MATRIX_TYPE, parameter :: CONE = (1.0, 0.0) CZERO = (0.0, 0.0)
!#endif
178

179
   call read_input_parameters_traditional(na, nev, nblk, write_to_file)
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
189
190
191
192
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

193
194
195
196
197
   if (myid == 0) then
     print '((a,i0))', 'Program ' // TEST_CASE
     print *, ""
   endif

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

   np_rows = nprocs/np_cols
213
   assert(nprocs == np_rows * np_cols)
214

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

227
228
   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, my_prow, my_pcol)
229
230
231
232

   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
233
234
   allocate(a (na_rows,na_cols))
   allocate(as(na_rows,na_cols))
235
236
237
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

238
239
240
241
242
#ifdef TEST_HERMITIAN_MULTIPLY
   allocate(b (na_rows,na_cols))
   allocate(c (na_rows,na_cols))
#endif

243
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_CHOLESKY)
244
245
246
247
248
   allocate(d (na), ds(na))
   allocate(sd (na), sds(na))
   allocate(ev_analytic(na))
#endif

249
250
251
252
   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

253
#if defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY)
Pavel Kus's avatar
Pavel Kus committed
254
255
#ifdef TEST_MATRIX_ANALYTIC
   call prepare_matrix_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
256
   as(:,:) = a
257
258
259
260
261
262
263
264
#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
265
#else
266
267
268
269
270
271
272
     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
273
274
275
276
277
278
279
280
281
282

#ifdef TEST_HERMITIAN_MULTIPLY
#if REALCASE == 1

#ifdef DOUBLE_PRECISION_REAL
   b(:,:) = 2.0_rk8 * a(:,:)
   c(:,:) = 0.0_rk8
#else
   b(:,:) = 2.0_rk4 * a(:,:)
   c(:,:) = 0.0_rk4
283
#endif
284

Pavel Kus's avatar
Pavel Kus committed
285
#endif
286

287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
#if COMPLEXCASE == 1

#ifdef DOUBLE_PRECISION_COMPLEX
   b(:,:) = 2.0_ck8 * a(:,:)
   c(:,:) = 0.0_ck8
#else
   b(:,:) = 2.0_ck4 * a(:,:)
   c(:,:) = 0.0_ck4
#endif

#endif

#endif /* TEST_HERMITIAN_MULTIPLY */

#endif /* TEST_MATRIX_ANALYTIC */
#endif /* defined(TEST_EIGENVECTORS) || defined(TEST_HERMITIAN_MULTIPLY) */

304
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
305
306
307
308
309
310
311
312
313
314
315

#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)
316
#endif /* EIGENVALUES OR TRIDIAGONAL */
317

318
319
320
#if defined(TEST_CHOLESKY)

#ifdef TEST_SINGLE
321
322
   diagonalElement = (2.546_c_float, 0.0_c_float)
   subdiagonalElement =  (0.0_c_float, 0.0_c_float)
323
#else
324
325
   diagonalElement = (2.546_c_double, 0.0_c_double)
   subdiagonalElement =  (0.0_c_double, 0.0_c_double)
326
327
328
329
#endif
   call prepare_toeplitz_matrix(na, diagonalElement, subdiagonalElement, &
                                d, sd, ds, sds, a, as, nblk, np_rows, &
                                np_cols, my_prow, my_pcol)
330
331


332
333
#endif /* TEST_CHOLESKY */

334
335
   e => elpa_allocate()

336
337
338
339
340
341
342
343
344
345
   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)
346
347

#ifdef WITH_MPI
348
349
350
351
352
353
   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)
354
#endif
355

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

358
359
360
361
362
363
364
   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
365
   assert_elpa_ok(error)
366

367
368
   call e%set("gpu", TEST_GPU, error)
   assert_elpa_ok(error)
369

370
371
   if (myid == 0) print *, ""

372
#ifdef TEST_ALL_KERNELS
373
   do i = 0, elpa_option_cardinality(KERNEL_KEY)  ! kernels
374
     kernel = elpa_option_enumerate(KERNEL_KEY, i)
375
#endif
376
#ifdef TEST_KERNEL
377
     kernel = TEST_KERNEL
378
#endif
379
380

#ifdef TEST_SOLVER_2STAGE
381
     call e%set(KERNEL_KEY, kernel, error)
382
383
384
#ifdef TEST_KERNEL
     assert_elpa_ok(error)
#else
385
386
387
     if (error /= ELPA_OK) then
       cycle
     endif
388
389
390
391
392
393
     ! 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
394
395
#endif

396
#ifdef TEST_ALL_KERNELS
397
     call e%timer_start(elpa_int_value_to_string(KERNEL_KEY, kernel))
398
#endif
399

400
     ! The actual solve step
401
402
#ifdef TEST_EIGENVECTORS
     call e%timer_start("e%eigenvectors()")
Pavel Kus's avatar
Pavel Kus committed
403
404
405
#ifdef TEST_SCALAPACK_ALL
     call solve_scalapack_all(na, a, sc_desc, ev, z)
#else
406
     call e%eigenvectors(a, ev, z, error)
Pavel Kus's avatar
Pavel Kus committed
407
#endif
408
     call e%timer_stop("e%eigenvectors()")
409
#endif /* TEST_EIGENVECTORS */
410
411
412

#ifdef TEST_EIGENVALUES
     call e%timer_start("e%eigenvalues()")
413
     call e%eigenvalues(a, ev, error)
414
     call e%timer_stop("e%eigenvalues()")
415
#endif
416
417
418

#if defined(TEST_SOLVE_TRIDIAGONAL)
     call e%timer_start("e%solve_tridiagonal()")
419
     call e%solve_tridiagonal(d, sd, z, error)
420
     call e%timer_stop("e%solve_tridiagonal()")
421
422
423
     ev(:) = d(:)
#endif

424
425
426
427
428
429
#if defined(TEST_CHOLESKY)
     call e%timer_start("e%cholesky()")
     call e%cholesky(a, error)
     call e%timer_stop("e%cholesky()")
#endif

430
431
432
433
434
#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
435

436
437
     assert_elpa_ok(error)

438
#ifdef TEST_ALL_KERNELS
439
440
441
     call e%timer_stop(elpa_int_value_to_string(KERNEL_KEY, kernel))
#endif

442
     if (myid .eq. 0) then
443
#ifdef TEST_ALL_KERNELS
444
       call e%print_times(elpa_int_value_to_string(KERNEL_KEY, kernel))
445
446
#else /* TEST_ALL_KERNELS */

447
#ifdef TEST_EIGENVECTORS
448
       call e%print_times("e%eigenvectors()")
449
#endif
450
#ifdef TEST_EIGENVALUES
451
452
       call e%print_times("e%eigenvalues()")
#endif
453
454
#ifdef TEST_SOLVE_TRIDIAGONAL
       call e%print_times("e%solve_tridiagonal()")
455
#endif
456
457
#ifdef TEST_CHOLESKY
       call e%print_times("e%cholesky()")
458
#endif
459
460
461
#ifdef TEST_HERMITIAN_MULTIPLY
       call e%print_times("e%hermitian_multiply()")
#endif
462
#endif /* TEST_ALL_KERNELS */
463
     endif
464

465
#ifdef TEST_EIGENVECTORS
Pavel Kus's avatar
Pavel Kus committed
466
467
468
#ifdef TEST_MATRIX_ANALYTIC
     status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
#else
469
470
471
472
473
474
475
     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
476
#endif
477
     call check_status(status, myid)
478
479
#endif

480
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL)
Andreas Marek's avatar
Andreas Marek committed
481
482
     status = check_correctness_eigenvalues_toeplitz(na, diagonalElement, &
         subdiagonalElement, ev, z, myid)
483
     call check_status(status, myid)
484

485
#ifdef TEST_SOLVE_TRIDIAGONAL
486
     ! check eigenvectors
Pavel Kus's avatar
Pavel Kus committed
487
     status = check_correctness(na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
488
     call check_status(status, myid)
489
490
#endif
#endif
491

492
#if defined(TEST_CHOLESKY)
Andreas Marek's avatar
Andreas Marek committed
493
494

!#ifdef TEST_REAL
495
496
     status = check_correctness_cholesky(na, a, as, na_rows, sc_desc, myid )
     call check_status(status, myid)
Andreas Marek's avatar
Andreas Marek committed
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
!#endif
        !-------------------------------------------------------------------------------
!#ifdef TEST_COMPLEX
!   ! Test correctness of result (using plain scalapack routines)
!   allocate(tmp1(na_rows,na_cols))
!   allocate(tmp2(na_rows,na_cols))
!
!   tmp1(:,:) = 0.0_ck8
!
!   ! tmp1 = a**H
!#ifdef WITH_MPI
!   call pztranc(na, na, CONE, a, 1, 1, sc_desc, CZERO, tmp1, 1, 1, sc_desc)
!#else
!   tmp1 = transpose(conjg(a))
!#endif
!   ! tmp2 = a * a**H
!#ifdef WITH_MPI
!   call pzgemm("N","N", na, na, na, CONE, a, 1, 1, sc_desc, tmp1, 1, 1, &
!               sc_desc, CZERO, tmp2, 1, 1, sc_desc)
!#else
!   call zgemm("N","N", na, na, na, CONE, a, na, tmp1, na, CZERO, tmp2, na)
!#endif
!
!   ! compare tmp2 with c
!   tmp2(:,:) = tmp2(:,:) - as(:,:)
!
!#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
!#ifdef WITH_MPI
!   call mpi_allreduce(norm,normmax,1,MPI_REAL8,MPI_MAX,MPI_COMM_WORLD,mpierr)
!#else
!   normmax = norm
!#endif
!   if (myid .eq. 0) then
!     print *," Maximum error of result: ", normmax
!   endif
!
!   if (normmax .gt. 5e-11_rk8) then
!        status = 1
!   endif
!
!   deallocate(tmp1, tmp2)
!#endif

544
545
#endif

546
547
548
549
550
551
#if defined(TEST_HERMITIAN_MULTIPLY)
     status = check_correctness_hermitian_multiply(na, a, b, c, na_rows, sc_desc, myid )
     call check_status(status, myid)
#endif


552
553
554
555
     if (myid == 0) then
       print *, ""
     endif

556
557
#ifdef TEST_ALL_KERNELS
     a(:,:) = as(:,:)
558
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_CHOLESKY)
559
560
561
     d = ds
     sd = sds
#endif
562
   end do ! kernels
563
#endif
Andreas Marek's avatar
Andreas Marek committed
564

565
566
567
568
569
570
571
   call elpa_deallocate(e)

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

572
573
574
575
576
#ifdef TEST_HERMITIAN_MULTIPLY
   deallocate(b)
   deallocate(c)
#endif

577
#if defined(TEST_EIGENVALUES) || defined(TEST_SOLVE_TRIDIAGONAL) || defined(TEST_EIGENVECTORS) || defined(TEST_CHOLESKY)
578
579
580
581
582
   deallocate(d, ds)
   deallocate(sd, sds)
   deallocate(ev_analytic)
#endif

583
584
585
586
587
588
589
#ifdef TEST_ALL_LAYOUTS
   end do ! factors
   end do ! layouts
#endif

   call elpa_uninit()

590
591
592
593
594
595
596
#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif

   call exit(status)

597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
   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

612
end program