test_analytic.F90 9.77 KB
Newer Older
Andreas Marek's avatar
Andreas Marek committed
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
! (c) Copyright Pavel Kus, 2017, MPCDF
!
!    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.

Pavel Kus's avatar
Pavel Kus committed
44
45
46
47
48
49
50
51
52
53
54
55
56
57
!#include "assert.h"  ! why complains?
module test_analytic

  use test_util

  interface prepare_matrix_analytic
    module procedure prepare_matrix_analytic_real_double
  end interface

  interface check_correctness_analytic
    module procedure check_correctness_analytic_real_double
  end interface

  contains
Andreas Marek's avatar
Andreas Marek committed
58

59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
    !Processor col for global col number
    pure function pcol(global_col, nblk, np_cols) result(local_col)
      implicit none
      integer(kind=c_int), intent(in) :: global_col, nblk, np_cols
      integer(kind=c_int)             :: local_col
      local_col = MOD((global_col-1)/nblk,np_cols)
    end function

    !Processor row for global row number
    pure function prow(global_row, nblk, np_rows) result(local_row)
      implicit none
      integer(kind=c_int), intent(in) :: global_row, nblk, np_rows
      integer(kind=c_int)             :: local_row
      local_row = MOD((global_row-1)/nblk,np_rows)
    end function

    function map_global_array_index_to_local_index(iGLobal, jGlobal, iLocal, jLocal , nblk, np_rows, np_cols, my_prow, my_pcol) &
      result(possible)
      implicit none

      integer(kind=c_int)              :: pi, pj, li, lj, xi, xj
      integer(kind=c_int), intent(in)  :: iGlobal, jGlobal, nblk, np_rows, np_cols, my_prow, my_pcol
      integer(kind=c_int), intent(out) :: iLocal, jLocal
      logical                       :: possible

      possible = .true.
      iLocal = 0
      jLocal = 0

      pi = prow(iGlobal, nblk, np_rows)

      if (my_prow .ne. pi) then
        possible = .false.
        return
      endif

      pj = pcol(jGlobal, nblk, np_cols)

      if (my_pcol .ne. pj) then
        possible = .false.
        return
      endif
      li = (iGlobal-1)/(np_rows*nblk) ! block number for rows
      lj = (jGlobal-1)/(np_cols*nblk) ! block number for columns

      xi = mod( (iGlobal-1),nblk)+1   ! offset in block li
      xj = mod( (jGlobal-1),nblk)+1   ! offset in block lj

      iLocal = li * nblk + xi
      jLocal = lj * nblk + xj

    end function


Pavel Kus's avatar
Pavel Kus committed
113
114
115
116
117
118
  subroutine prepare_matrix_analytic_real_double (na, a, nblk, myid, np_rows, &
                            np_cols, my_prow, my_pcol)
    implicit none

    integer(kind=ik), intent(in)    :: na, nblk, myid, np_rows, np_cols, my_prow, my_pcol
    real(kind=rk8), intent(inout)   :: a(:,:)
Andreas Marek's avatar
Andreas Marek committed
119

Pavel Kus's avatar
Pavel Kus committed
120
121
    integer(kind=ik) :: globI, globJ, locI, locJ, levels

Andreas Marek's avatar
Andreas Marek committed
122
    if(.not. decompose(na, levels)) then
Pavel Kus's avatar
Pavel Kus committed
123
124
      print *, "can not decomopse matrix size"
      stop 1
Andreas Marek's avatar
Andreas Marek committed
125
126
    end if

Pavel Kus's avatar
Pavel Kus committed
127
128
129
130
131
132
133
134
135
    do globI = 1, na
      do globJ = 1, na
        if(map_global_array_index_to_local_index(globI, globJ, locI, locJ, &
                 nblk, np_rows, np_cols, my_prow, my_pcol)) then
           a(locI, locJ) = analytic_matrix(levels, globI, globJ)
        end if
      end do
    end do

Andreas Marek's avatar
Andreas Marek committed
136
  end subroutine
Pavel Kus's avatar
Pavel Kus committed
137
138
139
140
141
142
143
144
145

  function check_correctness_analytic_real_double (na, nev, ev, z, nblk, myid, np_rows, &
                            np_cols, my_prow, my_pcol) result(status)
    implicit none

    integer(kind=ik), intent(in)    :: na, nev, nblk, myid, np_rows, np_cols, my_prow, my_pcol
    integer(kind=ik)                :: status
    real(kind=rk8), intent(inout)   :: z(:,:)
    real(kind=rk8), intent(inout)   :: ev(:)
Andreas Marek's avatar
Andreas Marek committed
146

Pavel Kus's avatar
Pavel Kus committed
147
148
149
150
    integer(kind=ik) :: globI, globJ, locI, locJ, levels
    real(kind=rk8)   :: diff, max_z_diff, max_ev_diff, glob_max_z_diff, max_curr_z_diff_minus, max_curr_z_diff_plus
    real(kind=rk8)   :: computed, expected

Andreas Marek's avatar
Andreas Marek committed
151
    if(.not. decompose(na, levels)) then
Pavel Kus's avatar
Pavel Kus committed
152
153
      print *, "can not decomopse matrix size"
      stop 1
Andreas Marek's avatar
Andreas Marek committed
154
155
    end if

Pavel Kus's avatar
Pavel Kus committed
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
    max_z_diff = 0.0_rk8
    max_ev_diff = 0.0_rk8
    do globJ = 1, nev
      diff = abs(ev(globJ) - analytic_eigenvalues(levels, globJ))
      max_ev_diff = max(diff, max_ev_diff)

      ! calculated eigenvector can be in opposite direction
      max_curr_z_diff_minus = 0.0_rk8
      max_curr_z_diff_plus  = 0.0_rk8
      do globI = 1, na
        if(map_global_array_index_to_local_index(globI, globJ, locI, locJ, &
                 nblk, np_rows, np_cols, my_prow, my_pcol)) then
           computed = z(locI, locJ)
           expected = analytic_eigenvectors(levels, globI, globJ)
           max_curr_z_diff_minus = max(abs(computed - expected), max_curr_z_diff_minus)
           max_curr_z_diff_plus = max(abs(computed + expected), max_curr_z_diff_plus)
        end if
      end do
      ! we have max difference of one of the eigenvectors, update global
      max_z_diff = max(max_z_diff, min(max_curr_z_diff_minus, max_curr_z_diff_plus))
    end do

#ifdef WITH_MPI
    call mpi_allreduce(max_z_diff, glob_max_z_diff, 1, MPI_REAL8, MPI_MAX, MPI_COMM_WORLD, mpierr)
#else
    glob_max_z_diff = max_z_diff
#endif
Andreas Marek's avatar
Andreas Marek committed
183
184
    if(myid == 0) print *, 'Maximal error in eigenvalues      :', max_ev_diff
    if(myid == 0) print *, 'Maximal error in eigenvectors     :', glob_max_z_diff
Pavel Kus's avatar
Pavel Kus committed
185
186
    status = 0
    if (max_ev_diff .gt. 5e-14_rk8 .or. max_ev_diff .eq. 0.0_rk8) status = 1
Pavel Kus's avatar
Pavel Kus committed
187
    if (glob_max_z_diff .gt. 1e-12_rk8 .or. glob_max_z_diff .eq. 0.0_rk8) status = 1
Andreas Marek's avatar
Andreas Marek committed
188
  end function
Pavel Kus's avatar
Pavel Kus committed
189
190
191
192
193
194
195

  function decompose(num, decomposition) result(possible)
    implicit none
    integer(kind=ik), intent(in)   :: num
    integer(kind=ik), intent(out)  :: decomposition
    logical                        :: possible
    integer(kind=ik)               :: reminder
Andreas Marek's avatar
Andreas Marek committed
196

Pavel Kus's avatar
Pavel Kus committed
197
    decomposition = 0
Andreas Marek's avatar
Andreas Marek committed
198
    possible = .true.
Pavel Kus's avatar
Pavel Kus committed
199
200
201
202
203
204
205
206
    reminder = num
    do while (reminder > 1)
      if (MOD(reminder, 2) == 0) then
        decomposition = decomposition + 1
        reminder = reminder / 2
      else
        possible = .false.
      end if
Andreas Marek's avatar
Andreas Marek committed
207
    end do
Pavel Kus's avatar
Pavel Kus committed
208
  end function
Andreas Marek's avatar
Andreas Marek committed
209

Pavel Kus's avatar
Pavel Kus committed
210
211
212
#define ANALYTIC_MATRIX 0
#define ANALYTIC_EIGENVECTORS 1
#define ANALYTIC_EIGENVALUES 2
Andreas Marek's avatar
Andreas Marek committed
213

Pavel Kus's avatar
Pavel Kus committed
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
  function analytic_matrix(levels, i, j) result(element)
    implicit none
    integer(kind=ik), intent(in) :: levels, i, j
    real(kind=rk8)               :: element

    element = analytic(levels, i, j, ANALYTIC_MATRIX)

  end function

  function analytic_eigenvectors(levels, i, j) result(element)
    implicit none
    integer(kind=ik), intent(in) :: levels, i, j
    real(kind=rk8)               :: element

    element = analytic(levels, i, j, ANALYTIC_EIGENVECTORS)

  end function

  function analytic_eigenvalues(levels, i) result(element)
    implicit none
    integer(kind=ik), intent(in) :: levels, i
    real(kind=rk8)               :: element

    element = analytic(levels, i, i, ANALYTIC_EIGENVALUES)

  end function

Andreas Marek's avatar
Andreas Marek committed
241

Pavel Kus's avatar
Pavel Kus committed
242
243
244
245
246
247

  function analytic(n, i, j, what) result(element)
    implicit none
    integer(kind=ik), intent(in)   :: n, i, j, what
    real(kind=rk8)                 :: element, am
    real(kind=rk8)                 :: a, s, c, mat(2,2)
Andreas Marek's avatar
Andreas Marek committed
248
    integer(kind=ik)               :: ii, jj, m
Pavel Kus's avatar
Pavel Kus committed
249
250
251
252

!    assert(i < 2**n)
!    assert(j < 2**n)
!    assert(i >= 0)
Andreas Marek's avatar
Andreas Marek committed
253
!    assert(j >= 0)
Pavel Kus's avatar
Pavel Kus committed
254
255
256
257
258
259
260
261
262
    ! go to zero-based indexing
    ii = i - 1
    jj = j - 1
    a = get_a(n)
    s = 0.5_rk8
    c = sqrt(3.0_rk8)/2.0_rk8
    element = 1.0_rk8
    do  m = 1, n
      am = a**(2**(m-1))
Andreas Marek's avatar
Andreas Marek committed
263
      if(what == ANALYTIC_MATRIX) then
Pavel Kus's avatar
Pavel Kus committed
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
        mat = reshape((/ c*c + am * s*s, (1.0_rk8-am) * s*c,  &
                         (1.0_rk8-am) * s*c, s*s + am * c*c  /), &
                                    (/2, 2/))
      else if(what == ANALYTIC_EIGENVECTORS) then
        mat = reshape((/ c, s,  &
                         -s,  c  /), &
                              (/2, 2/))
      else if(what == ANALYTIC_EIGENVALUES) then
        mat = reshape((/ 1.0_rk8, 0.0_rk8,  &
                         0.0_rk8, am  /), &
                               (/2, 2/))
      else
        !assert(0)
      end if
!      write(*,*) "calc value, elem: ", element, ", mat: ", mod(ii,2), mod(jj,2),  mat(mod(ii,2), mod(jj,2)), "am ", am
!      write(*,*) " matrix mat", mat
Andreas Marek's avatar
Andreas Marek committed
280
      element = element * mat(mod(ii,2) + 1, mod(jj,2) + 1)
Pavel Kus's avatar
Pavel Kus committed
281
282
283
284
285
      ii = ii / 2
      jj = jj / 2
    end do
    !write(*,*) "returning value ", element
  end function
Andreas Marek's avatar
Andreas Marek committed
286

Pavel Kus's avatar
Pavel Kus committed
287
288
289
290
291
292
293
294
295
296
297
298
  function get_a(n) result (a)
    implicit none
    integer(kind=ik), intent(in)   :: n
    real(kind=rk8)                 :: a
    real(kind=rk8), parameter      :: largest_ev = 2.0

    a = exp(log(largest_ev)/(2**n-1))
  end function



end module