test_analytic.F90 6.3 KB
Newer Older
Pavel Kus's avatar
Pavel Kus committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
!#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
15

Pavel Kus's avatar
Pavel Kus committed
16
17
18
19
20
21
22
  subroutine prepare_matrix_analytic_real_double (na, a, nblk, myid, np_rows, &
                            np_cols, my_prow, my_pcol)
    use elpa_utilities
    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
23

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

Andreas Marek's avatar
Andreas Marek committed
26
    if(.not. decompose(na, levels)) then
Pavel Kus's avatar
Pavel Kus committed
27
28
      print *, "can not decomopse matrix size"
      stop 1
Andreas Marek's avatar
Andreas Marek committed
29
30
    end if

Pavel Kus's avatar
Pavel Kus committed
31
32
33
34
35
36
37
38
39
    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
40
  end subroutine
Pavel Kus's avatar
Pavel Kus committed
41
42
43
44
45
46
47
48
49
50

  function check_correctness_analytic_real_double (na, nev, ev, z, nblk, myid, np_rows, &
                            np_cols, my_prow, my_pcol) result(status)
    use elpa_utilities
    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
51

Pavel Kus's avatar
Pavel Kus committed
52
53
54
55
    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
56
    if(.not. decompose(na, levels)) then
Pavel Kus's avatar
Pavel Kus committed
57
58
      print *, "can not decomopse matrix size"
      stop 1
Andreas Marek's avatar
Andreas Marek committed
59
60
    end if

Pavel Kus's avatar
Pavel Kus committed
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
    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
88
89
    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
90
91
    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
92
    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
93
  end function
Pavel Kus's avatar
Pavel Kus committed
94
95
96
97
98
99
100

  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
101

Pavel Kus's avatar
Pavel Kus committed
102
    decomposition = 0
Andreas Marek's avatar
Andreas Marek committed
103
    possible = .true.
Pavel Kus's avatar
Pavel Kus committed
104
105
106
107
108
109
110
111
    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
112
    end do
Pavel Kus's avatar
Pavel Kus committed
113
  end function
Andreas Marek's avatar
Andreas Marek committed
114

Pavel Kus's avatar
Pavel Kus committed
115
116
117
#define ANALYTIC_MATRIX 0
#define ANALYTIC_EIGENVECTORS 1
#define ANALYTIC_EIGENVALUES 2
Andreas Marek's avatar
Andreas Marek committed
118

Pavel Kus's avatar
Pavel Kus committed
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
  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
146

Pavel Kus's avatar
Pavel Kus committed
147
148
149
150
151
152

  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
153
    integer(kind=ik)               :: ii, jj, m
Pavel Kus's avatar
Pavel Kus committed
154
155
156
157

!    assert(i < 2**n)
!    assert(j < 2**n)
!    assert(i >= 0)
Andreas Marek's avatar
Andreas Marek committed
158
!    assert(j >= 0)
Pavel Kus's avatar
Pavel Kus committed
159
160
161
162
163
164
165
166
167
    ! 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
168
      if(what == ANALYTIC_MATRIX) then
Pavel Kus's avatar
Pavel Kus committed
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
        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
185
      element = element * mat(mod(ii,2) + 1, mod(jj,2) + 1)
Pavel Kus's avatar
Pavel Kus committed
186
187
188
189
190
      ii = ii / 2
      jj = jj / 2
    end do
    !write(*,*) "returning value ", element
  end function
Andreas Marek's avatar
Andreas Marek committed
191

Pavel Kus's avatar
Pavel Kus committed
192
193
194
195
196
197
198
199
200
201
202
203
  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