aboutsummaryrefslogtreecommitdiff
path: root/src/Coefficients_2_1.F90
blob: 8053c1c48b8b063219f49780105d114f1ab7fb65 (plain)
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
56
57
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
113
114
115
116
117
118
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
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
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
#include "cctk.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"

subroutine set_coeff_2_1 ( nsize, loc_order, bb, gsize, imin, imax, dd )

  use All_Coeffs_mod

  implicit none

  DECLARE_CCTK_FUNCTIONS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT, intent(IN) :: nsize, loc_order
  CCTK_INT, dimension(2), intent(IN) :: bb
  CCTK_INT, intent(IN) :: gsize
  CCTK_INT, dimension(nsize), intent(OUT) :: imin, imax
  CCTK_REAL, dimension(nsize,nsize), intent(OUT) :: dd

  CCTK_REAL, dimension(1), save :: a
  CCTK_REAL, dimension(3,2), save :: q 

  CCTK_INT :: i, il, ir

  logical, save :: first = .true.

  if ( first ) then
    call coeffs_1_2_1 ( a, q )
    first = .false.
  end if

  dd = zero
  imin = 0
  imax = -1

  if ( bb(1) == 0 ) then
    il = 1 + gsize
  else
    dd(1,1) = q(1,1); dd(2,1) = q(2,1)
    imin(1) = 1; imax(1) = 2
    dd(1,2) = q(1,2); dd(3,2) = q(3,2)
    imin(2) = 1; imax(2) = 3
    il = 3
  end if
  if ( bb(2) == 0 ) then
    ir = nsize - gsize
  else
    dd(nsize-2,nsize-1) = -q(3,2); dd(nsize,nsize-1) = -q(1,2)
    imin(nsize-1) = nsize-2; imax(nsize-1) = nsize
    dd(nsize-1,nsize) = -q(2,1); dd(nsize,nsize) = -q(1,1)
    imin(nsize) = nsize-1; imax(nsize) = nsize
    ir = nsize - 2
  end if
!$omp parallel do private(i)
  do i = il, ir
    dd(i-1,i) = -a(1); dd(i+1,i) = a(1)
    imin(i) = i-1; imax(i) = i+1
  end do

end subroutine set_coeff_2_1

subroutine set_coeff_up_2_1 ( nsize, dir, loc_order, bb, gsize, imin, imax, dd )

  use All_Coeffs_mod

  implicit none

  DECLARE_CCTK_FUNCTIONS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT, intent(IN) :: nsize, dir, loc_order
  CCTK_INT, dimension(2), intent(IN) :: bb
  CCTK_INT, intent(IN) :: gsize
  CCTK_INT, dimension(nsize), intent(OUT) :: imin, imax
  CCTK_REAL, dimension(nsize,nsize), intent(OUT) :: dd

  CCTK_REAL, dimension(-1:1), save :: a1, a2
  CCTK_REAL, dimension(3,2), save :: q1, q2

  CCTK_INT :: i, il, ir

  logical, save :: first = .true.

  if ( .not. ( dir == -1 .or. dir == 1 ) ) then
    call CCTK_WARN(0, 'Internal error: set_coeff_up_2_1 called with invalid dir value')
  end if

  if ( first ) then
    call coeffs_up_2_1 ( a1, q1, a2, q2 )
    first = .false.
  end if

  dd = zero
  imin = 0
  imax = -1

  if ( bb(1) == 0 ) then
    il = 1 + gsize
  else
    if ( dir == -1 ) then
      dd(1,1) = q1(1,1); dd(2,1) = q1(2,1)
      imin(1) = 1; imax(1) = 2
      dd(1,2) = q1(1,2); dd(2,2) = q1(2,2)
      imin(2) = 1; imax(2) = 2
    else
      dd(1,1) = q2(1,1); dd(2,1) = q2(2,1)
      imin(1) = 1; imax(1) = 2
      dd(2,2) = q2(2,2); dd(3,2) = q2(3,2)
      imin(2) = 2; imax(2) = 3
    end if
    il = 3
  end if
  if ( bb(2) == 0 ) then
    ir = nsize - gsize
  else
    if ( dir == -1 ) then
      dd(nsize-2,nsize-1) = -q2(3,2); dd(nsize-1,nsize-1) = -q2(2,2)
      imin(nsize-1) = nsize-2; imax(nsize-1) = nsize-1
      dd(nsize-1,nsize) = -q2(2,1); dd(nsize,nsize) = -q2(1,1)
      imin(nsize) = nsize-1; imax(nsize) = nsize
    else
      dd(nsize-1,nsize-1) = -q1(2,2); dd(nsize,nsize-1) = -q1(1,2)
      imin(nsize-1) = nsize-1; imax(nsize-1) = nsize
      dd(nsize-1,nsize) = -q1(2,1); dd(nsize,nsize) = -q1(1,1)
      imin(nsize) = nsize-1; imax(nsize) = nsize
    end if
    ir = nsize - 2
  end if
  do i = il, ir
    if ( dir == -1 ) then
      dd(i-1:i,i) = a1(-1:0)
      imin(i) = i-1; imax(i) = i
    else
      dd(i:i+1,i) = a2(0:1)
      imin(i) = i; imax(i) = i+1
    end if
  end do
end subroutine set_coeff_up_2_1

subroutine set_coeff_up_2 ( nsize, dir, loc_order, bb, gsize, imin, imax, dd )

  use All_Coeffs_mod

  implicit none

  DECLARE_CCTK_FUNCTIONS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT, intent(IN) :: nsize, dir, loc_order
  CCTK_INT, dimension(2), intent(IN) :: bb
  CCTK_INT, intent(IN) :: gsize
  CCTK_INT, dimension(nsize), intent(OUT) :: imin, imax
  CCTK_REAL, dimension(nsize,nsize), intent(OUT) :: dd

  CCTK_REAL, dimension(-2:2), save :: a1, a2
  CCTK_REAL, dimension(3,2), save :: q1, q2

  CCTK_INT :: i, il, ir

  logical, save :: first = .true.

  if ( .not. ( dir == -1 .or. dir == 1 ) ) then
    call CCTK_WARN(0, 'Internal error: set_coeff_up_2 called with invalid dir value')
  end if

  if ( first ) then
    call coeffs_up_2 ( a1, q1, a2, q2 )
    first = .false.
  end if

  dd = zero
  imin = 0
  imax = -1

  if ( bb(1) == 0 ) then
    il = 1 + gsize
  else
    if ( dir == -1 ) then
      dd(1,1) = q1(1,1);
      imin(1) = 1; imax(1) = 1

      dd(1:2,2) = q1(1:2,2);
      imin(2) = 1; imax(2) = 2
    else
      dd(1,1) = q2(1,1);
      imin(1) = 1; imax(1) = 1

      dd(2:4,2) = a2(0:2);
      imin(2) = 2; imax(2) = 4
    end if
    il = 3
  end if
  if ( bb(2) == 0 ) then
    ir = nsize - gsize
  else
    if ( dir == -1 ) then
      dd(nsize,nsize) = q1(1,1)
      imin(nsize) = nsize; imax(nsize) = nsize

      dd(nsize-3:nsize-1,nsize-1) = a1(-2:0)
      imin(nsize-1) = nsize-3; imax(nsize-1) = nsize-1
    else
      dd(nsize,nsize) = q2(1,1)
      imin(nsize) = nsize; imax(nsize) = nsize

      dd(nsize-1:nsize,nsize-1) = q2(2:3,2)
      imin(nsize-1) = nsize-1; imax(nsize-1) = nsize
    end if
    ir = nsize - 2
  end if
  do i = il, ir
    if ( dir == -1 ) then
      dd(i-2:i,i) = a1(-2:0)
      imin(i) = i-2; imax(i) = i
    else
      dd(i:i+2,i) = a2(0:2)
      imin(i) = i; imax(i) = i+2
    end if
  end do
end subroutine set_coeff_up_2