test_check_correctness_template.F90 25.5 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 */
200
201
202
203
204
205
206
      ! First check, whether the elements on diagonal are 1 .. "normality" of the vectors
      err = 0.0_rk
      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
           err = max(err, abs(tmp1(rowLocal,colLocal) - 1.0_rk))
         endif
      end do
Pavel Kus's avatar
Pavel Kus committed
207
#ifdef WITH_MPI
208
      call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
Pavel Kus's avatar
Pavel Kus committed
209
#else /* WITH_MPI */
210
      errmax = err
Pavel Kus's avatar
Pavel Kus committed
211
#endif /* WITH_MPI */
212
      if (myid==0) print *,'Maximal error in eigenvector lengths:',errmax
213

Pavel Kus's avatar
Pavel Kus committed
214
      ! Second, find the maximal error in the whole Z**T * Z matrix (its diference from identity matrix)
215
216
217
      ! Initialize tmp2 to unit matrix
      tmp2 = 0
#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
218
      call scal_PRECISION_LASET('A', nev, nev, ZERO, ONE, tmp2, 1, 1, sc_desc)
219
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
220
      call PRECISION_LASET('A',nev,nev,ZERO,ONE,tmp2,na)
221
222
223
224
225
226
227
228
#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
229
      call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
230
231
232
233
#else /* WITH_MPI */
      errmax = err
#endif /* WITH_MPI */
      if (myid==0) print *,'Error Orthogonality:',errmax
234
235
236
237
238
239
240
241
        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
242
        endif
243
      endif  ! skiping test of orthogonality for generalized eigenproblem
244
245
246
247
248
    end function


#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
249
    !c> int check_correctness_evp_numeric_residuals_real_double_f(int na, int nev, int na_rows, int na_cols,
250
251
    !c>                                         double *as, double *z, double *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol);
252
#else
253
    !c> int check_correctness_evp_numeric_residuals_real_single_f(int na, int nev, int na_rows, int na_cols,
254
255
    !c>                                         float *as, float *z, float *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol);
256
257
258
259
260
#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
261
    !c> int check_correctness_evp_numeric_residuals_complex_double_f(int na, int nev, int na_rows, int na_cols,
262
263
    !c>                                         complex double *as, complex double *z, double *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol);
264
#else
265
    !c> int check_correctness_evp_numeric_residuals_complex_single_f(int na, int nev, int na_rows, int na_cols,
266
267
    !c>                                         complex float *as, complex float *z, float *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol);
268
269
270
#endif
#endif /* COMPLEXCASE */

271
function check_correctness_evp_numeric_residuals_&
272
273
274
&MATH_DATATYPE&
&_&
&PRECISION&
275
&_f (na, nev, na_rows, na_cols, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol) result(status) &
276
      bind(C,name="check_correctness_evp_numeric_residuals_&
277
278
279
280
281
282
283
284
      &MATH_DATATYPE&
      &_&
      &PRECISION&
      &_f")

      use iso_c_binding

      implicit none
285
#include "../../src/general/precision_kinds.F90"
286
287

      integer(kind=c_int)            :: status
288
      integer(kind=c_int), value     :: na, nev, myid, na_rows, na_cols, nblk, np_rows, np_cols, my_prow, my_pcol
289
290
      MATH_DATATYPE(kind=rck)     :: as(1:na_rows,1:na_cols), z(1:na_rows,1:na_cols)
      real(kind=rck)    :: ev(1:na)
291
292
      integer(kind=c_int)            :: sc_desc(1:9)

293
      status = check_correctness_evp_numeric_residuals_&
294
295
296
      &MATH_DATATYPE&
      &_&
      &PRECISION&
297
      & (na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol)
298
299
300

    end function

301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
!---- variant for the generalized eigenproblem
!---- unlike in Fortran, we cannot use optional parameter
!---- we thus define a different function
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
    !c> int check_correctness_evp_gen_numeric_residuals_real_double_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         double *as, double *z, double *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol,
    !c>                                         double *bs);
#else
    !c> int check_correctness_evp_gen_numeric_residuals_real_single_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         float *as, float *z, float *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol, 
    !c>                                         float *bs);
#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
    !c> int check_correctness_evp_gen_numeric_residuals_complex_double_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         complex double *as, complex double *z, double *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol,
    !c>                                         complex double *bs);
#else
    !c> int check_correctness_evp_gen_numeric_residuals_complex_single_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         complex float *as, complex float *z, float *ev, int sc_desc[9],
    !c>                                         int nblk, int myid, int np_rows, int np_cols, int my_prow, int my_pcol,
    !c>                                         complex float *bs);
#endif
#endif /* COMPLEXCASE */

function check_correctness_evp_gen_numeric_residuals_&
&MATH_DATATYPE&
&_&
&PRECISION&
&_f (na, nev, na_rows, na_cols, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol, bs) result(status) &
      bind(C,name="check_correctness_evp_gen_numeric_residuals_&
      &MATH_DATATYPE&
      &_&
      &PRECISION&
      &_f")

      use iso_c_binding

      implicit none
#include "../../src/general/precision_kinds.F90"

      integer(kind=c_int)            :: status
      integer(kind=c_int), value     :: na, nev, myid, na_rows, na_cols, nblk, np_rows, np_cols, my_prow, my_pcol
      MATH_DATATYPE(kind=rck)     :: as(1:na_rows,1:na_cols), z(1:na_rows,1:na_cols), bs(1:na_rows,1:na_cols)
      real(kind=rck)    :: ev(1:na)
      integer(kind=c_int)            :: sc_desc(1:9)

      status = check_correctness_evp_numeric_residuals_&
      &MATH_DATATYPE&
      &_&
      &PRECISION&
      & (na, nev, as, z, ev, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol, bs)

    end function

    !-----------------------------------------------------------------------------------------------------------

Andreas Marek's avatar
Andreas Marek committed
364
365
366
367
368
369
370
    function check_correctness_eigenvalues_toeplitz_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, diagonalElement, subdiagonalElement, ev, z, myid) result(status)
      use iso_c_binding
      implicit none
371
#include "../../src/general/precision_kinds.F90"
Andreas Marek's avatar
Andreas Marek committed
372
373
374

      integer               :: status, ii, j, myid
      integer, intent(in)   :: na
375
376
377
      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
378
379

#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
380
      real(kind=rck), parameter   :: pi = 3.141592653589793238462643383279_c_double
Andreas Marek's avatar
Andreas Marek committed
381
#else
382
      real(kind=rck), parameter   :: pi = 3.1415926535897932_c_float
Andreas Marek's avatar
Andreas Marek committed
383
#endif
384
      real(kind=rck)              :: tmp, maxerr
385
      integer                     :: loctmp
Andreas Marek's avatar
Andreas Marek committed
386
387
388
389
      status = 0

     ! analytic solution
     do ii=1, na
390
391
392
       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
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
     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)
422
     if (maxerr .gt. 8.e-13_c_double .or. maxerr .eq. 0.0_c_double) then
Andreas Marek's avatar
Andreas Marek committed
423
#else
424
     if (maxerr .gt. 8.e-4_c_float .or. maxerr .eq. 0.0_c_float) then
Andreas Marek's avatar
Andreas Marek committed
425
426
427
#endif
       status = 1
       if (myid .eq. 0) then
Andreas Marek's avatar
Andreas Marek committed
428
         print *,"Result of Toeplitz matrix test: "
Andreas Marek's avatar
Andreas Marek committed
429
430
431
         print *,"Eigenvalues differ from analytic solution: maxerr = ",maxerr
       endif
     endif
Andreas Marek's avatar
Andreas Marek committed
432
433
434
435
436
437
438

    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
439
    end function
440

441
442
443
444
445
446
447
448
449
450
    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

451
452
453
      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
454
455

#ifdef WITH_MPI
456
457
458
      real(kind=rck)                   :: p&
                                           &BLAS_CHAR&
                                           &lange
459
#else /* WITH_MPI */
460
461
      real(kind=rck)                   :: BLAS_CHAR&
                                          &lange
462
463
464
465
466
467
468
#endif /* WITH_MPI */

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

      status = 0
469
      tmp1(:,:) = 0.0_rck
470
471
472
473
474


#if REALCASE == 1
      ! tmp1 = a**T
#ifdef WITH_MPI
475
476
477
      call p&
          &BLAS_CHAR&
          &tran(na, na, 1.0_rck, a, 1, 1, sc_desc, 0.0_rck, tmp1, 1, 1, sc_desc)
478
479
480
481
482
483
484
485
#else /* WITH_MPI */
      tmp1 = transpose(a)
#endif /* WITH_MPI */
#endif /* REALCASE == 1 */

#if COMPLEXCASE == 1
      ! tmp1 = a**H
#ifdef WITH_MPI
486
487
488
      call p&
            &BLAS_CHAR&
            &tranc(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
489
490
491
#else /* WITH_MPI */
      tmp1 = transpose(conjg(a))
#endif /* WITH_MPI */
492
#endif /* COMPLEXCASE == 1 */
493

Pavel Kus's avatar
Pavel Kus committed
494
      ! tmp2 = a**T * a
495
#ifdef WITH_MPI
496
497
      call p&
            &BLAS_CHAR&
Pavel Kus's avatar
Pavel Kus committed
498
            &gemm("N","N", na, na, na, ONE, tmp1, 1, 1, sc_desc, a, 1, 1, &
499
               sc_desc, ZERO, tmp2, 1, 1, sc_desc)
500
#else /* WITH_MPI */
501
      call BLAS_CHAR&
Pavel Kus's avatar
Pavel Kus committed
502
                    &gemm("N","N", na, na, na, ONE, tmp1, na, a, na, ZERO, tmp2, na)
503
504
505
506
#endif /* WITH_MPI */

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

508
#ifdef WITH_MPI
509
510
511
      norm = p&
              &BLAS_CHAR&
              &lange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
512
#else /* WITH_MPI */
513
514
      norm = BLAS_CHAR&
             &lange("M", na, na, tmp2, na_rows, tmp1)
515
516
517
518
#endif /* WITH_MPI */


#ifdef WITH_MPI
519
      call mpi_allreduce(norm,normmax,1,MPI_REAL_PRECISION,MPI_MAX,MPI_COMM_WORLD,mpierr)
520
521
522
523
524
525
526
527
#else /* WITH_MPI */
      normmax = norm
#endif /* WITH_MPI */

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

528
#if REALCASE == 1
529
#ifdef DOUBLE_PRECISION_REAL
530
!      if (normmax .gt. 5e-12_rk8 .or. normmax .eq. 0.0_rk8) then
531
532
533
534
      if (normmax .gt. 5e-12_rk8) then
        status = 1
      endif
#else
535
536
!      if (normmax .gt. 5e-4_rk4 .or. normmax .eq. 0.0_rk4) then
      if (normmax .gt. 5e-4_rk4 ) then
537
538
539
        status = 1
      endif
#endif
540
#endif
541

542
543
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
544
545
!      if (normmax .gt. 5e-11_rk8 .or. normmax .eq. 0.0_rk8) then
      if (normmax .gt. 5e-11_rk8 ) then
546
547
548
        status = 1
      endif
#else
549
!      if (normmax .gt. 5e-3_rk4 .or. normmax .eq. 0.0_rk4) then
550
551
552
553
554
      if (normmax .gt. 5e-3_rk4) then
        status = 1
      endif
#endif
#endif
555
556
    end function

557
558
559
560
561
562
563
564
565
    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
566
567
      MATH_DATATYPE(kind=rck), intent(in)       :: a(:,:), b(:,:), c(:,:)
      MATH_DATATYPE(kind=rck), dimension(size(a,dim=1),size(a,dim=2)) :: tmp1, tmp2
568
569
570
      real(kind=rck)                   :: norm, normmax

#ifdef WITH_MPI
571
572
573
      real(kind=rck)                   :: p&
                                           &BLAS_CHAR&
                                           &lange
574
#else /* WITH_MPI */
575
576
      real(kind=rck)                   :: BLAS_CHAR&
                                          &lange
577
578
579
580
581
582
583
#endif /* WITH_MPI */

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

      status = 0
584
      tmp1(:,:) = ZERO
585
586
587
588

#if REALCASE == 1
      ! tmp1 = a**T
#ifdef WITH_MPI
589
590
591
      call p&
            &BLAS_CHAR&
            &tran(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
592
593
594
595
596
597
598
599
600
#else /* WITH_MPI */
      tmp1 = transpose(a)
#endif /* WITH_MPI */

#endif /* REALCASE == 1 */

#if COMPLEXCASE == 1
      ! tmp1 = a**H
#ifdef WITH_MPI
601
602
603
      call p&
            &BLAS_CHAR&
            &tranc(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
604
605
606
#else /* WITH_MPI */
      tmp1 = transpose(conjg(a))
#endif /* WITH_MPI */
607
#endif /* COMPLEXCASE == 1 */
608
609
610

   ! tmp2 = tmp1 * b
#ifdef WITH_MPI
611
612
613
614
   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)
615
#else
616
617
   call BLAS_CHAR&
        &gemm("N","N", na, na, na, ONE, tmp1, na, b, na, ZERO, tmp2, na)
618
619
620
621
622
623
#endif

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

#ifdef WITH_MPI
624
625
626
      norm = p&
              &BLAS_CHAR&
              &lange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
627
#else /* WITH_MPI */
628
629
      norm = BLAS_CHAR&
             &lange("M", na, na, tmp2, na_rows, tmp1)
630
631
632
#endif /* WITH_MPI */

#ifdef WITH_MPI
633
      call mpi_allreduce(norm,normmax,1,MPI_REAL_PRECISION,MPI_MAX,MPI_COMM_WORLD,mpierr)
634
635
636
637
638
639
640
641
642
#else /* WITH_MPI */
      normmax = norm
#endif /* WITH_MPI */

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

#ifdef DOUBLE_PRECISION_REAL
643
      if (normmax .gt. 5e-11_rk8 ) then
644
645
646
        status = 1
      endif
#else
647
      if (normmax .gt. 5e-3_rk4 ) then
648
649
650
651
        status = 1
      endif
#endif

652
653
654
655
656
657
658
659
660
#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
661
662
    end function

Andreas Marek's avatar
Andreas Marek committed
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
    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
704

Andreas Marek's avatar
Andreas Marek committed
705
706
     tmp    = minval(ev_analytic)
     loctmp = minloc(ev_analytic, 1)
707

Andreas Marek's avatar
Andreas Marek committed
708
709
710
711
712
713
714
715
716
717
718
719
720
     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
721

Andreas Marek's avatar
Andreas Marek committed
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
     ! 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
738

739
! vim: syntax=fortran