test_analytic.F90 8.21 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.

44
#include "../Fortran/assert.h"  
Andreas Marek's avatar
Andreas Marek committed
45
46
#include "config-f90.h"

Pavel Kus's avatar
Pavel Kus committed
47
48
49
50
51
52
53
54
55
56
57
58
59
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
60

61
62
#include "../../src/general/prow_pcol.X90"
#include "../../src/general/map_global_to_local.X90"
63

Pavel Kus's avatar
Pavel Kus committed
64
65
66
67
68
69
  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
70

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

Andreas Marek's avatar
Andreas Marek committed
73
    if(.not. decompose(na, levels)) then
Pavel Kus's avatar
Pavel Kus committed
74
75
      print *, "can not decomopse matrix size"
      stop 1
Andreas Marek's avatar
Andreas Marek committed
76
77
    end if

Pavel Kus's avatar
Pavel Kus committed
78
79
80
81
82
83
84
85
86
    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
87
  end subroutine
Pavel Kus's avatar
Pavel Kus committed
88
89
90
91
92
93

  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
Andreas Marek's avatar
Andreas Marek committed
94
    integer(kind=ik)                :: status, mpierr
Pavel Kus's avatar
Pavel Kus committed
95
96
    real(kind=rk8), intent(inout)   :: z(:,:)
    real(kind=rk8), intent(inout)   :: ev(:)
Andreas Marek's avatar
Andreas Marek committed
97

Pavel Kus's avatar
Pavel Kus committed
98
99
100
101
    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
102
    if(.not. decompose(na, levels)) then
Pavel Kus's avatar
Pavel Kus committed
103
104
      print *, "can not decomopse matrix size"
      stop 1
Andreas Marek's avatar
Andreas Marek committed
105
106
    end if

Pavel Kus's avatar
Pavel Kus committed
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
    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
134
135
    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
136
137
    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
138
    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
139
  end function
Pavel Kus's avatar
Pavel Kus committed
140
141
142
143
144
145
146

  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
147

Pavel Kus's avatar
Pavel Kus committed
148
    decomposition = 0
Andreas Marek's avatar
Andreas Marek committed
149
    possible = .true.
Pavel Kus's avatar
Pavel Kus committed
150
151
152
153
154
155
156
157
    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
158
    end do
Pavel Kus's avatar
Pavel Kus committed
159
  end function
Andreas Marek's avatar
Andreas Marek committed
160

Pavel Kus's avatar
Pavel Kus committed
161
162
163
#define ANALYTIC_MATRIX 0
#define ANALYTIC_EIGENVECTORS 1
#define ANALYTIC_EIGENVALUES 2
Andreas Marek's avatar
Andreas Marek committed
164

Pavel Kus's avatar
Pavel Kus committed
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
190
191
  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
192

Pavel Kus's avatar
Pavel Kus committed
193
194
195
196
197
198

  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
199
    integer(kind=ik)               :: ii, jj, m
Pavel Kus's avatar
Pavel Kus committed
200
201
202
203

!    assert(i < 2**n)
!    assert(j < 2**n)
!    assert(i >= 0)
Andreas Marek's avatar
Andreas Marek committed
204
!    assert(j >= 0)
Pavel Kus's avatar
Pavel Kus committed
205
206
207
208
209
210
211
212
213
    ! 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
214
      if(what == ANALYTIC_MATRIX) then
Pavel Kus's avatar
Pavel Kus committed
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
        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
231
      element = element * mat(mod(ii,2) + 1, mod(jj,2) + 1)
Pavel Kus's avatar
Pavel Kus committed
232
233
234
235
236
      ii = ii / 2
      jj = jj / 2
    end do
    !write(*,*) "returning value ", element
  end function
Andreas Marek's avatar
Andreas Marek committed
237

Pavel Kus's avatar
Pavel Kus committed
238
239
240
241
242
243
244
245
246
247
248
249
  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