test_check_correctness_template.X90 13.4 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
49
50
51
52
53
54
55
!    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

    function check_correctness_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, nev, as, z, ev, sc_desc, myid) result(status)
      implicit none
      integer(kind=ik)                 :: status
      integer(kind=ik), intent(in)     :: na, nev, myid
#if REALCASE == 1
      real(kind=C_DATATYPE_KIND), intent(in)     :: as(:,:), z(:,:)
      real(kind=C_DATATYPE_KIND)                 :: ev(:)
      real(kind=C_DATATYPE_KIND), dimension(size(as,dim=1),size(as,dim=2)) :: tmp1, tmp2
Pavel Kus's avatar
Pavel Kus committed
56
      real(kind=C_DATATYPE_KIND)                :: xc
57

58
59
60
61
62
63
#ifdef DOUBLE_PRECISION_REAL
      real(kind=C_DATATYPE_KIND), parameter      :: ZERO=0.0_rk8, ONE = 1.0_rk8
#else
      real(kind=C_DATATYPE_KIND), parameter      :: ZERO=0.0_rk4, ONE = 1.0_rk4
#endif

64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
#ifndef WITH_MPI

#ifdef DOUBLE_PRECISION_REAL
      real(kind=C_DATATYPE_KIND)                 :: dnrm2
#else
      real(kind=C_DATATYPE_KIND)                 :: snrm2
#endif

#endif
#endif /* REALCASE */

#if COMPLEXCASE == 1
      complex(kind=C_DATATYPE_KIND), intent(in)    :: as(:,:), z(:,:)
      real(kind=C_DATATYPE_KIND)                   :: ev(:)
      complex(kind=C_DATATYPE_KIND), dimension(size(as,dim=1),size(as,dim=2)) :: tmp1, tmp2
      complex(kind=C_DATATYPE_KIND)                :: xc
#ifdef DOUBLE_PRECISION_COMPLEX
81
      complex(kind=C_DATATYPE_KIND), parameter     :: ZERO = (0.0_rk8,0.0_rk8), ONE = (1.0_rk8,0.0_rk8)
82
83
84
85
86
#ifndef WITH_MPI
      complex(kind=C_DATATYPE_KIND)                :: zdotc, cdotc
#endif

#else /* DOUBLE_PRECISION_COMPLEX */
87
      complex(kind=C_DATATYPE_KIND), parameter     :: ZERO = (0.0_rk4,0.0_rk4), ONE = (1.0_rk4,0.0_rk4)
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
#ifndef WITH_MPI
      complex(kind=C_DATATYPE_KIND)                :: zdotc, cdotc
#endif

#endif /* DOUBLE_PRECISION_COMPLEX */

#endif /* COMPLEXCASE */

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

      integer(kind=ik)                 :: i
      real(kind=C_DATATYPE_KIND)                   :: err, errmax

      integer :: mpierr

      status = 0

      ! 1. Residual (maximum of || A*Zi - Zi*EVi ||)
      ! tmp1 =  A * Z
      ! as is original stored matrix, Z are the EVs

#ifdef WITH_MPI
110
111
      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)
112
#else /* WITH_MPI */
113
      call PRECISION_GEMM('N','N',na,nev,na,ONE,as,na,z,na,ZERO,tmp1,na)
114
115
116
117
118
119
#endif /* WITH_MPI */


      ! tmp2 = Zi*EVi
      tmp2(:,:) = z(:,:)
      do i=1,nev
Pavel Kus's avatar
Pavel Kus committed
120
        xc = ev(i)
121
122
123
124
#if REALCASE == 1
#ifdef WITH_MPI

#ifdef DOUBLE_PRECISION_REAL
Pavel Kus's avatar
Pavel Kus committed
125
        call pdscal(na, xc, tmp2, 1, i, sc_desc, 1)
126
#else
Pavel Kus's avatar
Pavel Kus committed
127
        call psscal(na, xc, tmp2, 1, i, sc_desc, 1)
128
129
130
131
132
#endif

#else /* WITH_MPI */

#ifdef DOUBLE_PRECISION_REAL
Pavel Kus's avatar
Pavel Kus committed
133
        call dscal(na,xc,tmp2(:,i),1)
134
#else
Pavel Kus's avatar
Pavel Kus committed
135
        call sscal(na,xc,tmp2(:,i),1)
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
#endif

#endif /* WITH_MPI */
#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef WITH_MPI

#ifdef DOUBLE_PRECISION_COMPLEX
        call pzscal(na, xc, tmp2, 1, i, sc_desc, 1)
#else
        call pcscal(na, xc, tmp2, 1, i, sc_desc, 1)
#endif

#else /* WITH_MPI */

#ifdef DOUBLE_PRECISION_COMPLEX
        call zscal(na,xc,tmp2(1,i),1)
#else
        call cscal(na,xc,tmp2(1,i),1)
#endif

#endif /* WITH_MPI */
#endif /* COMPLEXCASE */
      enddo

      !  tmp1 = A*Zi - Zi*EVi
      tmp1(:,:) =  tmp1(:,:) - tmp2(:,:)

      ! Get maximum norm of columns of tmp1
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
      errmax = 0.0_rk8
#else
      errmax = 0.0_rk4
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
      errmax = 0.0_rk8
#else
      errmax = 0.0_rk4
#endif
#endif

      do i=1,nev
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
        err = 0.0_rk8
#else
        err = 0.0_rk4
#endif

#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
190
        call scal_PRECISION_NRM2(na, err, tmp1, 1, i, sc_desc, 1)
191
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
192
        err = PRECISION_NRM2(na,tmp1(1,i),1)
193
194
195
196
197
198
199
#endif /* WITH_MPI */
        errmax = max(errmax, err)
#endif /* REALCASE */

#if COMPLEXCASE == 1
        xc = 0
#ifdef WITH_MPI
Pavel Kus's avatar
Pavel Kus committed
200
        call scal_PRECISION_DOTC(na, xc, tmp1, 1, i, sc_desc, 1, tmp1, 1, i, sc_desc, 1)
201
#else /* WITH_MPI */
Pavel Kus's avatar
Pavel Kus committed
202
        xc = PRECISION_DOTC(na,tmp1,1,tmp1,1)
203
204
205
206
207
208
209
210
211
212
213
214
215
216
#endif /* WITH_MPI */


#ifdef DOUBLE_PRECISION_COMPLEX
        errmax = max(errmax, sqrt(real(xc,kind=rk8)))
#else
        errmax = max(errmax, sqrt(real(xc,kind=rk4)))
#endif
#endif /* COMPLEXCASE */
      enddo

      ! Get maximum error norm over all processors
      err = errmax
#ifdef WITH_MPI
217
      call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
#else /* WITH_MPI */
      errmax = err
#endif /* WITH_MPI */

      if (myid==0) print *,'Error Residual     :',errmax
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
      if (errmax .gt. 5e-12_rk8 .or. errmax .eq. 0.0_rk8) then
#else
      if (errmax .gt. 3e-3_rk4 .or. errmax .eq. 0.0_rk4) then
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
      if (errmax .gt. 5e-12_rk8 .or. errmax .eq. 0.0_rk8) then
#else
      if (errmax .gt. 3e-3_rk4 .or. errmax .eq. 0.0_rk4) then
#endif
#endif

        status = 1
      endif

      ! 2. Eigenvector orthogonality

      ! tmp1 = Z**T * Z
      tmp1 = 0
#ifdef WITH_MPI
246
247
      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)
248
#else /* WITH_MPI */
249
      call PRECISION_GEMM(BLAS_TRANS_OR_CONJ,'N',nev,nev,na,ONE,z,na,z,na,ZERO,tmp1,na)
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
#endif /* WITH_MPI */

      ! Initialize tmp2 to unit matrix
      tmp2 = 0
#if REALCASE == 1

#ifdef WITH_MPI

#ifdef DOUBLE_PRECISION_REAL
      call pdlaset('A', nev, nev, 0.0_rk8, 1.0_rk8, tmp2, 1, 1, sc_desc)
#else
      call pslaset('A', nev, nev, 0.0_rk4, 1.0_rk4, tmp2, 1, 1, sc_desc)
#endif

#else /* WITH_MPI */

#ifdef DOUBLE_PRECISION_REAL
      call dlaset('A',nev,nev,0.0_rk8,1.0_rk8,tmp2,na)
#else
      call slaset('A',nev,nev,0.0_rk4,1.0_rk4,tmp2,na)
#endif

#endif /* WITH_MPI */

#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef WITH_MPI

#ifdef DOUBLE_PRECISION_COMPLEX
280
      call pzlaset('A', nev, nev, ZERO, ONE, tmp2, 1, 1, sc_desc)
281
#else
282
      call pclaset('A', nev, nev, ZERO, ONE, tmp2, 1, 1, sc_desc)
283
284
285
286
287
#endif

#else /* WITH_MPI */

#ifdef DOUBLE_PRECISION_COMPLEX
288
      call zlaset('A',nev,nev,ZERO,ONE,tmp2,na)
289
#else
290
      call claset('A',nev,nev,ZERO,ONE,tmp2,na)
291
292
293
294
295
296
297
298
299
300
301
302
#endif

#endif /* WITH_MPI */
#endif /* COMPLEXCASE */


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

      ! Get maximum error (max abs value in tmp1)
      err = maxval(abs(tmp1))
#ifdef WITH_MPI
303
      call mpi_allreduce(err, errmax, 1, MPI_REAL_PRECISION, MPI_MAX, MPI_COMM_WORLD, mpierr)
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
#else /* WITH_MPI */
      errmax = err
#endif /* WITH_MPI */

      if (myid==0) print *,'Error Orthogonality:',errmax
#if REALCASE == 1
#ifdef DOUBLE_PRECISION_REAL
      if (errmax .gt. 5e-12_rk8 .or. errmax .eq. 0.0_rk8) then
#else
      if (errmax .gt. 9e-4_rk4 .or. errmax .eq. 0.0_rk4) then
#endif
#endif
#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
      if (errmax .gt. 5e-12_rk8 .or. errmax .eq. 0.0_rk8) then
#else
      if (errmax .gt. 9e-4_rk4 .or. errmax .eq. 0.0_rk4) then
#endif
#endif

        status = 1
      endif
    end function


#if REALCASE == 1

#ifdef DOUBLE_PRECISION_REAL
    !c> int check_correctness_real_double_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         double *as, double *z, double *ev,
    !c>                                         int sc_desc[9], int myid);
#else
    !c> int check_correctness_real_single_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         float *as, float *z, float *ev,
    !c>                                         int sc_desc[9], int myid);
#endif

#endif /* REALCASE */

#if COMPLEXCASE == 1
#ifdef DOUBLE_PRECISION_COMPLEX
    !c> int check_correctness_complex_double_f(int na, int nev, int na_rows, int na_cols,
    !c>                                         complex double *as, complex double *z, double *ev,
    !c>                                         int sc_desc[9], int myid);
#else
    !c> int check_correctness_complex_single_f(int na, int nev, int na_rows, int na_cols,
Andreas Marek's avatar
Andreas Marek committed
350
    !c>                                         complex float *as, complex float *z, float *ev,
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
    !c>                                         int sc_desc[9], int myid);
#endif
#endif /* COMPLEXCASE */

function check_correctness_&
&MATH_DATATYPE&
&_&
&PRECISION&
&_f (na, nev, na_rows, na_cols, as, z, ev, sc_desc, myid) result(status) &
      bind(C,name="check_correctness_&
      &MATH_DATATYPE&
      &_&
      &PRECISION&
      &_f")

      use iso_c_binding

      implicit none

      integer(kind=c_int)            :: status
      integer(kind=c_int), value     :: na, nev, myid, na_rows, na_cols
#if REALCASE == 1
      real(kind=C_DATATYPE_KIND)     :: as(1:na_rows,1:na_cols), z(1:na_rows,1:na_cols)
#endif

#if COMPLEXCASE == 1
      complex(kind=C_DATATYPE_KIND) :: as(1:na_rows,1:na_cols), z(1:na_rows,1:na_cols)
#endif
      real(kind=C_DATATYPE_KIND)    :: ev(1:na)
      integer(kind=c_int)            :: sc_desc(1:9)

      status = check_correctness_&
      &MATH_DATATYPE&
      &_&
      &PRECISION&
      & (na, nev, as, z, ev, sc_desc, myid)

    end function

Andreas Marek's avatar
Andreas Marek committed
390
391
392
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
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
    function check_correctness_eigenvalues_toeplitz_&
    &MATH_DATATYPE&
    &_&
    &PRECISION&
    & (na, diagonalElement, subdiagonalElement, ev, z, myid) result(status)
      use iso_c_binding
      implicit none

      integer               :: status, ii, j, myid
      integer, intent(in)   :: na
      real(kind=C_DATATYPE_KIND) :: diagonalElement, subdiagonalElement
      real(kind=C_DATATYPE_KIND) :: ev_analytic(na), ev(na)
#if REALCASE == 1
      real(kind=C_DATATYPE_KIND) :: z(:,:)
#endif
#if COMPLEXCASE == 1
      complex(kind=C_DATATYPE_KIND) :: z(:,:)
#endif

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

     ! analytic solution
     do ii=1, na
#if defined(DOUBLE_PRECISION_REAL) || defined(DOUBLE_PRECISION_COMPLEX)
       ev_analytic(ii) = diagonalElement + 2.0_c_double * &
                         subdiagonalElement *cos( pi*real(ii,kind=c_double)/ &
			 real(na+1,kind=c_double) )
#else
       ev_analytic(ii) = diagonalElement + 2.0_c_float * &
                         subdiagonalElement *cos( pi*real(ii,kind=c_float)/ &
                         real(na+1,kind=c_float) )
#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

     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)
     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 *,"Eigenvalues differ from analytic solution: maxerr = ",maxerr
       endif
     endif
    end function
468
469

! vim: syntax=fortran