test_check_correctness_template.F90 19.9 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
!    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
77
78
79
      ! tolerance for the residual test for different math type/precision setups
      real(kind=rk), parameter       :: tol_res_real_double      = 5e-12_rk
      real(kind=rk), parameter       :: tol_res_real_single      = 3e-3_rk
      real(kind=rk), parameter       :: tol_res_complex_double   = 5e-12_rk
      real(kind=rk), parameter       :: tol_res_complex_single   = 3e-3_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
88
89
90
91
92
93
94
95
96
      ! tolerance for the orthogonality test for different math type/precision setups
      real(kind=rk), parameter       :: tol_orth_real_double     = 5e-12_rk
      real(kind=rk), parameter       :: tol_orth_real_single     = 9e-4_rk
      real(kind=rk), parameter       :: tol_orth_complex_double  = 5e-12_rk
      real(kind=rk), parameter       :: tol_orth_complex_single  = 9e-4_rk
      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
      if (nev .ge. 2) then
240
        if (errmax .gt. tol_orth .or. errmax .eq. 0.0_rk) then
241
242
243
          status = 1
        endif
      else
244
        if (errmax .gt. tol_orth) then
245
246
247
          status = 1
        endif
      endif
248
249

      endif  ! skiping test of orthogonality for generalized eigenproblem
250
251
252
253
254
255
    end function


#if REALCASE == 1

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

#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
269
    !c> int check_correctness_evp_numeric_residuals_complex_double_f(int na, int nev, int na_rows, int na_cols,
270
271
272
    !c>                                         complex double *as, complex double *z, double *ev,
    !c>                                         int sc_desc[9], int myid);
#else
273
    !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
274
    !c>                                         complex float *as, complex float *z, float *ev,
275
276
277
278
    !c>                                         int sc_desc[9], int myid);
#endif
#endif /* COMPLEXCASE */

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

      use iso_c_binding

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

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

Pavel Kus's avatar
Pavel Kus committed
301
302
303
      ! 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..)
304
      status = check_correctness_evp_numeric_residuals_&
305
306
307
      &MATH_DATATYPE&
      &_&
      &PRECISION&
Pavel Kus's avatar
Pavel Kus committed
308
      & (na, nev, as, z, ev, sc_desc, 0, myid, 0, 0, 0, 0)
309
310
311

    end function

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

      integer               :: status, ii, j, myid
      integer, intent(in)   :: na
323
324
325
      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
326
327

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

     ! analytic solution
     do ii=1, na
338
339
340
       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
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
369
     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)
370
     if (maxerr .gt. 8.e-13_c_double .or. maxerr .eq. 0.0_c_double) then
Andreas Marek's avatar
Andreas Marek committed
371
#else
372
     if (maxerr .gt. 8.e-4_c_float .or. maxerr .eq. 0.0_c_float) then
Andreas Marek's avatar
Andreas Marek committed
373
374
375
376
377
378
379
#endif
       status = 1
       if (myid .eq. 0) then
         print *,"Eigenvalues differ from analytic solution: maxerr = ",maxerr
       endif
     endif
    end function
380

381
382
383
384
385
386
387
388
389
390
    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

391
392
393
      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
394
395

#ifdef WITH_MPI
396
397
398
      real(kind=rck)                   :: p&
                                           &BLAS_CHAR&
                                           &lange
399
#else /* WITH_MPI */
400
401
      real(kind=rck)                   :: BLAS_CHAR&
                                          &lange
402
403
404
405
406
407
408
#endif /* WITH_MPI */

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

      status = 0
409
      tmp1(:,:) = 0.0_rck
410
411
412
413
414


#if REALCASE == 1
      ! tmp1 = a**T
#ifdef WITH_MPI
415
416
417
      call p&
          &BLAS_CHAR&
          &tran(na, na, 1.0_rck, a, 1, 1, sc_desc, 0.0_rck, tmp1, 1, 1, sc_desc)
418
419
420
421
422
423
424
425
#else /* WITH_MPI */
      tmp1 = transpose(a)
#endif /* WITH_MPI */
#endif /* REALCASE == 1 */

#if COMPLEXCASE == 1
      ! tmp1 = a**H
#ifdef WITH_MPI
426
427
428
      call p&
            &BLAS_CHAR&
            &tranc(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
429
430
431
#else /* WITH_MPI */
      tmp1 = transpose(conjg(a))
#endif /* WITH_MPI */
432
#endif /* COMPLEXCASE == 1 */
433

Pavel Kus's avatar
Pavel Kus committed
434
      ! tmp2 = a**T * a
435
#ifdef WITH_MPI
436
437
      call p&
            &BLAS_CHAR&
Pavel Kus's avatar
Pavel Kus committed
438
            &gemm("N","N", na, na, na, ONE, tmp1, 1, 1, sc_desc, a, 1, 1, &
439
               sc_desc, ZERO, tmp2, 1, 1, sc_desc)
440
#else /* WITH_MPI */
441
      call BLAS_CHAR&
Pavel Kus's avatar
Pavel Kus committed
442
                    &gemm("N","N", na, na, na, ONE, tmp1, na, a, na, ZERO, tmp2, na)
443
444
445
446
#endif /* WITH_MPI */

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

448
#ifdef WITH_MPI
449
450
451
      norm = p&
              &BLAS_CHAR&
              &lange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
452
#else /* WITH_MPI */
453
454
      norm = BLAS_CHAR&
             &lange("M", na, na, tmp2, na_rows, tmp1)
455
456
457
458
#endif /* WITH_MPI */


#ifdef WITH_MPI
459
      call mpi_allreduce(norm,normmax,1,MPI_REAL_PRECISION,MPI_MAX,MPI_COMM_WORLD,mpierr)
460
461
462
463
464
465
466
467
#else /* WITH_MPI */
      normmax = norm
#endif /* WITH_MPI */

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

468
#if REALCASE == 1
469
#ifdef DOUBLE_PRECISION_REAL
470
!      if (normmax .gt. 5e-12_rk8 .or. normmax .eq. 0.0_rk8) then
471
472
473
474
      if (normmax .gt. 5e-12_rk8) then
        status = 1
      endif
#else
475
476
!      if (normmax .gt. 5e-4_rk4 .or. normmax .eq. 0.0_rk4) then
      if (normmax .gt. 5e-4_rk4 ) then
477
478
479
        status = 1
      endif
#endif
480
#endif
481

482
483
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
484
485
!      if (normmax .gt. 5e-11_rk8 .or. normmax .eq. 0.0_rk8) then
      if (normmax .gt. 5e-11_rk8 ) then
486
487
488
        status = 1
      endif
#else
489
!      if (normmax .gt. 5e-3_rk4 .or. normmax .eq. 0.0_rk4) then
490
491
492
493
494
      if (normmax .gt. 5e-3_rk4) then
        status = 1
      endif
#endif
#endif
495
496
    end function

497
498
499
500
501
502
503
504
505
    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
506
507
      MATH_DATATYPE(kind=rck), intent(in)       :: a(:,:), b(:,:), c(:,:)
      MATH_DATATYPE(kind=rck), dimension(size(a,dim=1),size(a,dim=2)) :: tmp1, tmp2
508
509
510
      real(kind=rck)                   :: norm, normmax

#ifdef WITH_MPI
511
512
513
      real(kind=rck)                   :: p&
                                           &BLAS_CHAR&
                                           &lange
514
#else /* WITH_MPI */
515
516
      real(kind=rck)                   :: BLAS_CHAR&
                                          &lange
517
518
519
520
521
522
523
#endif /* WITH_MPI */

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

      status = 0
524
      tmp1(:,:) = ZERO
525
526
527
528

#if REALCASE == 1
      ! tmp1 = a**T
#ifdef WITH_MPI
529
530
531
      call p&
            &BLAS_CHAR&
            &tran(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
532
533
534
535
536
537
538
539
540
#else /* WITH_MPI */
      tmp1 = transpose(a)
#endif /* WITH_MPI */

#endif /* REALCASE == 1 */

#if COMPLEXCASE == 1
      ! tmp1 = a**H
#ifdef WITH_MPI
541
542
543
      call p&
            &BLAS_CHAR&
            &tranc(na, na, ONE, a, 1, 1, sc_desc, ZERO, tmp1, 1, 1, sc_desc)
544
545
546
#else /* WITH_MPI */
      tmp1 = transpose(conjg(a))
#endif /* WITH_MPI */
547
#endif /* COMPLEXCASE == 1 */
548
549
550

   ! tmp2 = tmp1 * b
#ifdef WITH_MPI
551
552
553
554
   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)
555
#else
556
557
   call BLAS_CHAR&
        &gemm("N","N", na, na, na, ONE, tmp1, na, b, na, ZERO, tmp2, na)
558
559
560
561
562
563
#endif

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

#ifdef WITH_MPI
564
565
566
      norm = p&
              &BLAS_CHAR&
              &lange("M",na, na, tmp2, 1, 1, sc_desc, tmp1)
567
#else /* WITH_MPI */
568
569
      norm = BLAS_CHAR&
             &lange("M", na, na, tmp2, na_rows, tmp1)
570
571
572
#endif /* WITH_MPI */

#ifdef WITH_MPI
573
      call mpi_allreduce(norm,normmax,1,MPI_REAL_PRECISION,MPI_MAX,MPI_COMM_WORLD,mpierr)
574
575
576
577
578
579
580
581
582
#else /* WITH_MPI */
      normmax = norm
#endif /* WITH_MPI */

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

#ifdef DOUBLE_PRECISION_REAL
583
      if (normmax .gt. 5e-11_rk8 ) then
584
585
586
        status = 1
      endif
#else
587
      if (normmax .gt. 5e-3_rk4 ) then
588
589
590
591
        status = 1
      endif
#endif

592
593
594
595
596
597
598
599
600
#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
601
602
603
    end function


604
605


606
! vim: syntax=fortran