test_check_correctness_template.F90 22.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
!    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.
!
! Author: A. Marek, MPCDF

44
    function check_correctness_evp_numeric_residuals_&
45
46
47
    &MATH_DATATYPE&
    &_&
    &PRECISION&
48
    & (na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol, bs) result(status)
49
      implicit none
50
#include "../../src/general/precision_kinds.F90"
51
      integer(kind=ik)                 :: status
Pavel Kus's avatar
Pavel Kus committed
52
      integer(kind=ik), intent(in)     :: na, nev, nblk, myid, np_rows, np_cols, my_prow, my_pcol
53
54
      MATH_DATATYPE(kind=rck), intent(in)           :: as(:,:), z(:,:)
      MATH_DATATYPE(kind=rck), intent(in), optional :: bs(:,:)
55
56
57
      real(kind=rk)                 :: ev(:)
      MATH_DATATYPE(kind=rck), dimension(size(as,dim=1),size(as,dim=2)) :: tmp1, tmp2
      MATH_DATATYPE(kind=rck)                :: xc
58

59
#ifndef WITH_MPI
60
61
#if REALCASE == 1
      real(kind=rck)                   :: dnrm2, snrm2
62
63
#endif
#if COMPLEXCASE == 1
64
      complex(kind=rck)                :: zdotc, cdotc
65
#endif /* COMPLEXCASE */
66
#endif
67
68
69

      integer(kind=ik)                 :: sc_desc(:)

Pavel Kus's avatar
Pavel Kus committed
70
      integer(kind=ik)                 :: i, rowLocal, colLocal
71
      real(kind=rck)                   :: err, errmax
72
73
74

      integer :: mpierr

75
76
      ! tolerance for the residual test for different math type/precision setups
      real(kind=rk), parameter       :: tol_res_real_double      = 5e-12_rk
77
      real(kind=rk), parameter       :: tol_res_real_single      = 3e-2_rk
78
      real(kind=rk), parameter       :: tol_res_complex_double   = 5e-12_rk
79
      real(kind=rk), parameter       :: tol_res_complex_single   = 3e-2_rk
80
      real(kind=rk)                  :: tol_res                  = tol_res_&
81
82
83
                                                                          &MATH_DATATYPE&
                                                                          &_&
                                                                          &PRECISION
84
85
86
      ! precision of generalized problem is lower
      real(kind=rk), parameter       :: generalized_penalty = 10.0_rk

87
      ! tolerance for the orthogonality test for different math type/precision setups
88
89
90
91
      real(kind=rk), parameter       :: tol_orth_real_double     = 5e-11_rk
      real(kind=rk), parameter       :: tol_orth_real_single     = 9e-2_rk
      real(kind=rk), parameter       :: tol_orth_complex_double  = 5e-11_rk
      real(kind=rk), parameter       :: tol_orth_complex_single  = 9e-3_rk
92
93
94
95
96
      real(kind=rk), parameter       :: tol_orth                 = tol_orth_&
                                                                          &MATH_DATATYPE&
                                                                          &_&
                                                                          &PRECISION

97
98
99
      if(present(bs)) then
          tol_res = generalized_penalty * tol_res
      endif
100
101
102
103
      status = 0

      ! 1. Residual (maximum of || A*Zi - Zi*EVi ||)

104
105
      ! tmp1 = Zi*EVi
      tmp1(:,:) = z(:,:)
106
      do i=1,nev
Pavel Kus's avatar
Pavel Kus committed
107
        xc = ev(i)
108
#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
109
110
        call p&
            &BLAS_CHAR&
111
            &scal(na, xc, tmp1, 1, i, sc_desc, 1)
112
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
113
        call BLAS_CHAR&
114
            &scal(na,xc,tmp1(:,i),1)
115
116
117
#endif /* WITH_MPI */
      enddo

118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
      ! for generalized EV problem, multiply by bs as well
      ! tmp2 = B * tmp1
      if(present(bs)) then
#ifdef WITH_MPI
      call scal_PRECISION_GEMM('N', 'N', na, nev, na, ONE, bs, 1, 1, sc_desc, &
                  tmp1, 1, 1, sc_desc, ZERO, tmp2, 1, 1, sc_desc)
#else /* WITH_MPI */
      call PRECISION_GEMM('N','N',na,nev,na,ONE,bs,na,tmp1,na,ZERO,tmp2,na)
#endif /* WITH_MPI */
      else
        ! normal eigenvalue problem .. no need to multiply
        tmp2(:,:) = tmp1(:,:)
      end if

      ! tmp1 =  A * Z
      ! as is original stored matrix, Z are the EVs
#ifdef WITH_MPI
      call scal_PRECISION_GEMM('N', 'N', na, nev, na, ONE, as, 1, 1, sc_desc, &
                  z, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
#else /* WITH_MPI */
      call PRECISION_GEMM('N','N',na,nev,na,ONE,as,na,z,na,ZERO,tmp1,na)
#endif /* WITH_MPI */

141
142
143
144
      !  tmp1 = A*Zi - Zi*EVi
      tmp1(:,:) =  tmp1(:,:) - tmp2(:,:)

      ! Get maximum norm of columns of tmp1
145
      errmax = 0.0_rk
146
147
148

      do i=1,nev
#if REALCASE == 1
149
        err = 0.0_rk
150
#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
151
        call scal_PRECISION_NRM2(na, err, tmp1, 1, i, sc_desc, 1)
152
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
153
        err = PRECISION_NRM2(na,tmp1(1,i),1)
154
155
156
157
158
159
160
#endif /* WITH_MPI */
        errmax = max(errmax, err)
#endif /* REALCASE */

#if COMPLEXCASE == 1
        xc = 0
#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
161
        call scal_PRECISION_DOTC(na, xc, tmp1, 1, i, sc_desc, 1, tmp1, 1, i, sc_desc, 1)
162
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
163
        xc = PRECISION_DOTC(na,tmp1,1,tmp1,1)
164
#endif /* WITH_MPI */
165
        errmax = max(errmax, sqrt(real(xc,kind=REAL_DATATYPE)))
166
167
168
169
170
171
#endif /* COMPLEXCASE */
      enddo

      ! Get maximum error norm over all processors
      err = errmax
#ifdef WITH_MPI
172
      call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
173
174
175
#else /* WITH_MPI */
      errmax = err
#endif /* WITH_MPI */
176
      if (myid==0) print *,'Results of numerical residual checks:'
177
      if (myid==0) print *,'Error Residual     :',errmax
178
      if (nev .ge. 2) then
179
        if (errmax .gt. tol_res .or. errmax .eq. 0.0_rk) then
180
181
182
          status = 1
        endif
      else
183
        if (errmax .gt. tol_res) then
184
185
          status = 1
        endif
186
187
188
      endif

      ! 2. Eigenvector orthogonality
189
190
191
      !TODO for the generalized eigenvector problem, the orthogonality test has to be altered
      !TODO at the moment, it is skipped
      if(.not. present(bs)) then
192
193
194
      ! tmp1 = Z**T * Z
      tmp1 = 0
#ifdef WITH_MPI
195
196
      call scal_PRECISION_GEMM(BLAS_TRANS_OR_CONJ, 'N', nev, nev, na, ONE, z, 1, 1, &
                        sc_desc, z, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
197
#else /* WITH_MPI */
198
      call PRECISION_GEMM(BLAS_TRANS_OR_CONJ,'N',nev,nev,na,ONE,z,na,z,na,ZERO,tmp1,na)
199
#endif /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
200
201
      !TODO for the C interface, not all information is passed (zeros instead)
      !TODO than this part of the test cannot be done
202
      !TODO either we will not have this part of test at all, or it will be improved
Pavel Kus's avatar
Pavel Kus committed
203
204
      if(nblk > 0) then
        ! First check, whether the elements on diagonal are 1 .. "normality" of the vectors
205
        err = 0.0_rk
Pavel Kus's avatar
Pavel Kus committed
206
207
        do i=1, nev
          if (map_global_array_index_to_local_index(i, i, rowLocal, colLocal, nblk, np_rows, np_cols, my_prow, my_pcol)) then
208
             err = max(err, abs(tmp1(rowLocal,colLocal) - 1.0_rk))
Pavel Kus's avatar
Pavel Kus committed
209
210
211
212
213
214
215
216
217
           endif
        end do
#ifdef WITH_MPI
        call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
#else /* WITH_MPI */
        errmax = err
#endif /* WITH_MPI */
        if (myid==0) print *,'Maximal error in eigenvector lengths:',errmax
      end if
218

Pavel Kus's avatar
Pavel Kus committed
219
      ! Second, find the maximal error in the whole Z**T * Z matrix (its diference from identity matrix)
220
221
222
      ! Initialize tmp2 to unit matrix
      tmp2 = 0
#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
223
      call scal_PRECISION_LASET('A', nev, nev, ZERO, ONE, tmp2, 1, 1, sc_desc)
224
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
225
      call PRECISION_LASET('A',nev,nev,ZERO,ONE,tmp2,na)
226
227
228
229
230
231
232
233
#endif /* WITH_MPI */

      !      ! tmp1 = Z**T * Z - Unit Matrix
      tmp1(:,:) =  tmp1(:,:) - tmp2(:,:)

      ! Get maximum error (max abs value in tmp1)
      err = maxval(abs(tmp1))
#ifdef WITH_MPI
234
      call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
235
236
237
238
#else /* WITH_MPI */
      errmax = err
#endif /* WITH_MPI */
      if (myid==0) print *,'Error Orthogonality:',errmax
239
240
241
242
243
244
245
246
        if (nev .ge. 2) then
          if (errmax .gt. tol_orth .or. errmax .eq. 0.0_rk) then
            status = 1
          endif
        else
          if (errmax .gt. tol_orth) then
            status = 1
          endif
247
        endif
248
      endif  ! skiping test of orthogonality for generalized eigenproblem
249
250
251
252
253
254
    end function


#if REALCASE == 1

#ifdef DOUBLE_PRECISION_REAL
255
    !c> int check_correctness_evp_numeric_residuals_real_double_f(int na, int nev, int na_rows, int na_cols,
256
257
258
    !c>                                         double *as, double *z, double *ev,
    !c>                                         int sc_desc[9], int myid);
#else
259
    !c> int check_correctness_evp_numeric_residuals_real_single_f(int na, int nev, int na_rows, int na_cols,
260
261
262
263
264
265
266
267
    !c>                                         float *as, float *z, float *ev,
    !c>                                         int sc_desc[9], int myid);
#endif

#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
268
    !c> int check_correctness_evp_numeric_residuals_complex_double_f(int na, int nev, int na_rows, int na_cols,
269
270
271
    !c>                                         complex double *as, complex double *z, double *ev,
    !c>                                         int sc_desc[9], int myid);
#else
272
    !c> int check_correctness_evp_numeric_residuals_complex_single_f(int na, int nev, int na_rows, int na_cols,
Andreas Marek's avatar
Andreas Marek committed
273
    !c>                                         complex float *as, complex float *z, float *ev,
274
275
276
277
    !c>                                         int sc_desc[9], int myid);
#endif
#endif /* COMPLEXCASE */

278
function check_correctness_evp_numeric_residuals_&
279
280
281
282
&MATH_DATATYPE&
&_&
&PRECISION&
&_f (na, nev, na_rows, na_cols, as, z, ev, sc_desc, myid) result(status) &
283
      bind(C,name="check_correctness_evp_numeric_residuals_&
284
285
286
287
288
289
290
291
      &MATH_DATATYPE&
      &_&
      &PRECISION&
      &_f")

      use iso_c_binding

      implicit none
292
#include "../../src/general/precision_kinds.F90"
293
294
295

      integer(kind=c_int)            :: status
      integer(kind=c_int), value     :: na, nev, myid, na_rows, na_cols
296
297
      MATH_DATATYPE(kind=rck)     :: as(1:na_rows,1:na_cols), z(1:na_rows,1:na_cols)
      real(kind=rck)    :: ev(1:na)
298
299
      integer(kind=c_int)            :: sc_desc(1:9)

Pavel Kus's avatar
Pavel Kus committed
300
301
302
      ! TODO: I did not want to add all the variables to the C interface as well
      ! TODO: I think that we should find a better way to pass this information
      ! TODO: to all the functions anyway (get it from sc_desc, pass elpa_t, etc..)
303
      status = check_correctness_evp_numeric_residuals_&
304
305
306
      &MATH_DATATYPE&
      &_&
      &PRECISION&
Pavel Kus's avatar
Pavel Kus committed
307
      & (na, nev, as, z, ev, sc_desc, 0, myid, 0, 0, 0, 0)
308
309
310

    end function

Andreas Marek's avatar
Andreas Marek committed
311
312
313
314
315
316
317
    function check_correctness_eigenvalues_toeplitz_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, diagonalElement, subdiagonalElement, ev, z, myid) result(status)
      use iso_c_binding
      implicit none
318
#include "../../src/general/precision_kinds.F90"
Andreas Marek's avatar
Andreas Marek committed
319
320
321

      integer               :: status, ii, j, myid
      integer, intent(in)   :: na
322
323
324
      real(kind=rck) :: diagonalElement, subdiagonalElement
      real(kind=rck) :: ev_analytic(na), ev(na)
      MATH_DATATYPE(kind=rck) :: z(:,:)
Andreas Marek's avatar
Andreas Marek committed
325
326

#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
327
      real(kind=rck), parameter   :: pi = 3.141592653589793238462643383279_c_double
Andreas Marek's avatar
Andreas Marek committed
328
#else
329
      real(kind=rck), parameter   :: pi = 3.1415926535897932_c_float
Andreas Marek's avatar
Andreas Marek committed
330
#endif
331
      real(kind=rck)              :: tmp, maxerr
332
      integer                     :: loctmp
Andreas Marek's avatar
Andreas Marek committed
333
334
335
336
      status = 0

     ! analytic solution
     do ii=1, na
337
338
339
       ev_analytic(ii) = diagonalElement + 2.0_rk * &
                         subdiagonalElement *cos( pi*real(ii,kind=rk)/ &
                         real(na+1,kind=rk) )
Andreas Marek's avatar
Andreas Marek committed
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
     enddo

     ! sort analytic solution:

     ! this hack is neither elegant, nor optimized: for huge matrixes it might be expensive
     ! a proper sorting algorithmus might be implemented here

     tmp    = minval(ev_analytic)
     loctmp = minloc(ev_analytic, 1)

     ev_analytic(loctmp) = ev_analytic(1)
     ev_analytic(1) = tmp
     do ii=2, na
       tmp = ev_analytic(ii)
       do j= ii, na
         if (ev_analytic(j) .lt. tmp) then
           tmp    = ev_analytic(j)
           loctmp = j
         endif
       enddo
       ev_analytic(loctmp) = ev_analytic(ii)
       ev_analytic(ii) = tmp
     enddo

     ! compute a simple error max of eigenvalues
     maxerr = 0.0
     maxerr = maxval( (ev(:) - ev_analytic(:))/ev_analytic(:) , 1)

#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
369
     if (maxerr .gt. 8.e-13_c_double .or. maxerr .eq. 0.0_c_double) then
Andreas Marek's avatar
Andreas Marek committed
370
#else
371
     if (maxerr .gt. 8.e-4_c_float .or. maxerr .eq. 0.0_c_float) then
Andreas Marek's avatar
Andreas Marek committed
372
373
374
#endif
       status = 1
       if (myid .eq. 0) then
Andreas Marek's avatar
Andreas Marek committed
375
         print *,"Result of Toeplitz matrix test: "
Andreas Marek's avatar
Andreas Marek committed
376
377
378
         print *,"Eigenvalues differ from analytic solution: maxerr = ",maxerr
       endif
     endif
Andreas Marek's avatar
Andreas Marek committed
379
380
381
382
383
384
385

    if (status .eq. 0) then
       if (myid .eq. 0) then
         print *,"Result of Toeplitz matrix test: test passed"
         print *,"Eigenvalues differ from analytic solution: maxerr = ",maxerr
       endif
    endif
Andreas Marek's avatar
Andreas Marek committed
386
    end function
387

388
389
390
391
392
393
394
395
396
397
    function check_correctness_cholesky_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, a, as, na_rows, sc_desc, myid) result(status)
      implicit none
#include "../../src/general/precision_kinds.F90"
      integer(kind=ik)                 :: status
      integer(kind=ik), intent(in)     :: na, myid, na_rows

398
399
400
      MATH_DATATYPE(kind=rck), intent(in)       :: a(:,:), as(:,:)
      MATH_DATATYPE(kind=rck), dimension(size(as,dim=1),size(as,dim=2)) :: tmp1, tmp2
      real(kind=rk)                   :: norm, normmax
401
402

#ifdef WITH_MPI
403
404
405
      real(kind=rck)                   :: p&
                                           &BLAS_CHAR&
                                           &lange
406
#else /* WITH_MPI */
407
408
      real(kind=rck)                   :: BLAS_CHAR&
                                          &lange
409
410
411
412
413
414
415
#endif /* WITH_MPI */

      integer(kind=ik)                 :: sc_desc(:)
      real(kind=rck)                   :: err, errmax
      integer :: mpierr

      status = 0
416
      tmp1(:,:) = 0.0_rck
417
418
419
420
421


#if REALCASE == 1
      ! tmp1 = a**T
#ifdef WITH_MPI
422
423
424
      call p&
          &BLAS_CHAR&
          &tran(na, na, 1.0_rck, a, 1, 1, sc_desc, 0.0_rck, tmp1, 1, 1, sc_desc)
425
426
427
428
429
430
431
432
#else /* WITH_MPI */
      tmp1 = transpose(a)
#endif /* WITH_MPI */
#endif /* REALCASE == 1 */

#if COMPLEXCASE == 1
      ! tmp1 = a**H
#ifdef WITH_MPI
433
434
435
      call p&
            &BLAS_CHAR&
            &tranc(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
436
437
438
#else /* WITH_MPI */
      tmp1 = transpose(conjg(a))
#endif /* WITH_MPI */
439
#endif /* COMPLEXCASE == 1 */
440

Pavel Kus's avatar
Pavel Kus committed
441
      ! tmp2 = a**T * a
442
#ifdef WITH_MPI
443
444
      call p&
            &BLAS_CHAR&
Pavel Kus's avatar
Pavel Kus committed
445
            &gemm("N","N", na, na, na, ONE, tmp1, 1, 1, sc_desc, a, 1, 1, &
446
               sc_desc, ZERO, tmp2, 1, 1, sc_desc)
447
#else /* WITH_MPI */
448
      call BLAS_CHAR&
Pavel Kus's avatar
Pavel Kus committed
449
                    &gemm("N","N", na, na, na, ONE, tmp1, na, a, na, ZERO, tmp2, na)
450
451
452
453
#endif /* WITH_MPI */

      ! compare tmp2 with original matrix
      tmp2(:,:) = tmp2(:,:) - as(:,:)
454

455
#ifdef WITH_MPI
456
457
458
      norm = p&
              &BLAS_CHAR&
              &lange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
459
#else /* WITH_MPI */
460
461
      norm = BLAS_CHAR&
             &lange("M", na, na, tmp2, na_rows, tmp1)
462
463
464
465
#endif /* WITH_MPI */


#ifdef WITH_MPI
466
      call mpi_allreduce(norm,normmax,1,MPI_REAL_PRECISION,MPI_MAX,MPI_COMM_WORLD,mpierr)
467
468
469
470
471
472
473
474
#else /* WITH_MPI */
      normmax = norm
#endif /* WITH_MPI */

      if (myid .eq. 0) then
        print *," Maximum error of result: ", normmax
      endif

475
#if REALCASE == 1
476
#ifdef DOUBLE_PRECISION_REAL
477
!      if (normmax .gt. 5e-12_rk8 .or. normmax .eq. 0.0_rk8) then
478
479
480
481
      if (normmax .gt. 5e-12_rk8) then
        status = 1
      endif
#else
482
483
!      if (normmax .gt. 5e-4_rk4 .or. normmax .eq. 0.0_rk4) then
      if (normmax .gt. 5e-4_rk4 ) then
484
485
486
        status = 1
      endif
#endif
487
#endif
488

489
490
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
491
492
!      if (normmax .gt. 5e-11_rk8 .or. normmax .eq. 0.0_rk8) then
      if (normmax .gt. 5e-11_rk8 ) then
493
494
495
        status = 1
      endif
#else
496
!      if (normmax .gt. 5e-3_rk4 .or. normmax .eq. 0.0_rk4) then
497
498
499
500
501
      if (normmax .gt. 5e-3_rk4) then
        status = 1
      endif
#endif
#endif
502
503
    end function

504
505
506
507
508
509
510
511
512
    function check_correctness_hermitian_multiply_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, a, b, c, na_rows, sc_desc, myid) result(status)
      implicit none
#include "../../src/general/precision_kinds.F90"
      integer(kind=ik)                 :: status
      integer(kind=ik), intent(in)     :: na, myid, na_rows
513
514
      MATH_DATATYPE(kind=rck), intent(in)       :: a(:,:), b(:,:), c(:,:)
      MATH_DATATYPE(kind=rck), dimension(size(a,dim=1),size(a,dim=2)) :: tmp1, tmp2
515
516
517
      real(kind=rck)                   :: norm, normmax

#ifdef WITH_MPI
518
519
520
      real(kind=rck)                   :: p&
                                           &BLAS_CHAR&
                                           &lange
521
#else /* WITH_MPI */
522
523
      real(kind=rck)                   :: BLAS_CHAR&
                                          &lange
524
525
526
527
528
529
530
#endif /* WITH_MPI */

      integer(kind=ik)                 :: sc_desc(:)
      real(kind=rck)                   :: err, errmax
      integer :: mpierr

      status = 0
531
      tmp1(:,:) = ZERO
532
533
534
535

#if REALCASE == 1
      ! tmp1 = a**T
#ifdef WITH_MPI
536
537
538
      call p&
            &BLAS_CHAR&
            &tran(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
539
540
541
542
543
544
545
546
547
#else /* WITH_MPI */
      tmp1 = transpose(a)
#endif /* WITH_MPI */

#endif /* REALCASE == 1 */

#if COMPLEXCASE == 1
      ! tmp1 = a**H
#ifdef WITH_MPI
548
549
550
      call p&
            &BLAS_CHAR&
            &tranc(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
551
552
553
#else /* WITH_MPI */
      tmp1 = transpose(conjg(a))
#endif /* WITH_MPI */
554
#endif /* COMPLEXCASE == 1 */
555
556
557

   ! tmp2 = tmp1 * b
#ifdef WITH_MPI
558
559
560
561
   call p&
         &BLAS_CHAR&
         &gemm("N","N", na, na, na, ONE, tmp1, 1, 1, sc_desc, b, 1, 1, &
               sc_desc, ZERO, tmp2, 1, 1, sc_desc)
562
#else
563
564
   call BLAS_CHAR&
        &gemm("N","N", na, na, na, ONE, tmp1, na, b, na, ZERO, tmp2, na)
565
566
567
568
569
570
#endif

      ! compare tmp2 with c
      tmp2(:,:) = tmp2(:,:) - c(:,:)

#ifdef WITH_MPI
571
572
573
      norm = p&
              &BLAS_CHAR&
              &lange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
574
#else /* WITH_MPI */
575
576
      norm = BLAS_CHAR&
             &lange("M", na, na, tmp2, na_rows, tmp1)
577
578
579
#endif /* WITH_MPI */

#ifdef WITH_MPI
580
      call mpi_allreduce(norm,normmax,1,MPI_REAL_PRECISION,MPI_MAX,MPI_COMM_WORLD,mpierr)
581
582
583
584
585
586
587
588
589
#else /* WITH_MPI */
      normmax = norm
#endif /* WITH_MPI */

      if (myid .eq. 0) then
        print *," Maximum error of result: ", normmax
      endif

#ifdef DOUBLE_PRECISION_REAL
590
      if (normmax .gt. 5e-11_rk8 ) then
591
592
593
        status = 1
      endif
#else
594
      if (normmax .gt. 5e-3_rk4 ) then
595
596
597
598
        status = 1
      endif
#endif

599
600
601
602
603
604
605
606
607
#ifdef DOUBLE_PRECISION_COMPLEX
      if (normmax .gt. 5e-11_rk8 ) then
        status = 1
      endif
#else
      if (normmax .gt. 5e-3_rk4 ) then
        status = 1
      endif
#endif
608
609
    end function

Andreas Marek's avatar
Andreas Marek committed
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
    function check_correctness_eigenvalues_frank_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, ev, z, myid) result(status)
      use iso_c_binding
      implicit none
#include "../../src/general/precision_kinds.F90"

      integer                   :: status, i, j, myid
      integer, intent(in)       :: na
      real(kind=rck)            :: ev_analytic(na), ev(na)
      MATH_DATATYPE(kind=rck)   :: z(:,:)

#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
      real(kind=rck), parameter :: pi = 3.141592653589793238462643383279_c_double
#else
      real(kind=rck), parameter :: pi = 3.1415926535897932_c_float
#endif
      real(kind=rck)            :: tmp, maxerr
      integer                   :: loctmp
      status = 0

     ! analytic solution
     do i = 1, na
       j = na - i
#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
       ev_analytic(i) = pi * (2.0_c_double * real(j,kind=c_double) + 1.0_c_double) / &
           (2.0_c_double * real(na,kind=c_double) + 1.0_c_double)
       ev_analytic(i) = 0.5_c_double / (1.0_c_double - cos(ev_analytic(i)))
#else
       ev_analytic(i) = pi * (2.0_c_float * real(j,kind=c_float) + 1.0_c_float) / &
           (2.0_c_float * real(na,kind=c_float) + 1.0_c_float)
       ev_analytic(i) = 0.5_c_float / (1.0_c_float - cos(ev_analytic(i)))
#endif
     enddo

     ! sort analytic solution:

     ! this hack is neither elegant, nor optimized: for huge matrixes it might be expensive
     ! a proper sorting algorithmus might be implemented here
651

Andreas Marek's avatar
Andreas Marek committed
652
653
     tmp    = minval(ev_analytic)
     loctmp = minloc(ev_analytic, 1)
654

Andreas Marek's avatar
Andreas Marek committed
655
656
657
658
659
660
661
662
663
664
665
666
667
     ev_analytic(loctmp) = ev_analytic(1)
     ev_analytic(1) = tmp
     do i=2, na
       tmp = ev_analytic(i)
       do j= i, na
         if (ev_analytic(j) .lt. tmp) then
           tmp    = ev_analytic(j)
           loctmp = j
         endif
       enddo
       ev_analytic(loctmp) = ev_analytic(i)
       ev_analytic(i) = tmp
     enddo
668

Andreas Marek's avatar
Andreas Marek committed
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
     ! compute a simple error max of eigenvalues
     maxerr = 0.0
     maxerr = maxval( (ev(:) - ev_analytic(:))/ev_analytic(:) , 1)

#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
     if (maxerr .gt. 8.e-13_c_double) then
#else
     if (maxerr .gt. 8.e-4_c_float) then
#endif
       status = 1
       if (myid .eq. 0) then
         print *,"Result of Frank matrix test: "
         print *,"Eigenvalues differ from analytic solution: maxerr = ",maxerr
       endif
     endif
    end function
685

686
! vim: syntax=fortran