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
221
222
223
224
225
226
227
|
#include "cctk.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"
subroutine deriv_gf_4_3 ( var, ni, nj, nk, dir, bb, gsize, offset, delta, dvar )
use All_Coeffs_mod
implicit none
DECLARE_CCTK_FUNCTIONS
DECLARE_CCTK_PARAMETERS
CCTK_INT, intent(IN) :: ni, nj, nk
CCTK_REAL, dimension(ni,nj,nk), intent(IN) :: var
CCTK_INT, intent(IN) :: dir
CCTK_INT, intent(IN) :: bb(2)
CCTK_INT, intent(IN) :: gsize
CCTK_INT, intent(IN) :: offset(2)
CCTK_REAL, intent(IN) :: delta
CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar
CCTK_REAL, dimension(2), save :: a
CCTK_REAL, dimension(7,5), save :: q
CCTK_REAL :: idel
CCTK_INT :: il, ir, jl, jr, kl, kr, ol, or
logical, save :: first = .true.
if ( first ) then
call coeffs_1_4_3 ( a, q )
first = .false.
end if
idel = 1.0_wp / delta
if (gsize < 2) call CCTK_WARN (0, "not enough ghostzones")
direction: select case (dir)
case (0) direction
if ( bb(1) == 0 ) then
il = 1 + gsize
else
ol = offset(1)
!$omp parallel workshare
dvar(1+ol,:,:) = ( q(1,1) * var(1+ol,:,:) + q(2,1) * var(2+ol,:,:) + &
q(3,1) * var(3+ol,:,:) + q(4,1) * var(4+ol,:,:) ) * idel
dvar(2+ol,:,:) = ( q(1,2) * var(1+ol,:,:) + q(2,2) * var(2+ol,:,:) + &
q(3,2) * var(3+ol,:,:) + q(4,2) * var(4+ol,:,:) + &
q(5,2) * var(5+ol,:,:) + q(6,2) * var(6+ol,:,:) ) * idel
dvar(3+ol,:,:) = ( q(1,3) * var(1+ol,:,:) + q(2,3) * var(2+ol,:,:) + &
q(3,3) * var(3+ol,:,:) + q(4,3) * var(4+ol,:,:) + &
q(5,3) * var(5+ol,:,:) + q(6,3) * var(6+ol,:,:) ) * idel
dvar(4+ol,:,:) = ( q(1,4) * var(1+ol,:,:) + q(2,4) * var(2+ol,:,:) + &
q(3,4) * var(3+ol,:,:) + q(4,4) * var(4+ol,:,:) + &
q(5,4) * var(5+ol,:,:) + q(6,4) * var(6+ol,:,:) + &
q(7,4) * var(7+ol,:,:) ) * idel
dvar(5+ol,:,:) = ( q(1,5) * var(1+ol,:,:) + q(2,5) * var(2+ol,:,:) + &
q(3,5) * var(3+ol,:,:) + q(4,5) * var(4+ol,:,:) + &
q(5,5) * var(5+ol,:,:) + q(6,5) * var(6+ol,:,:) + &
q(7,5) * var(7+ol,:,:) ) * idel
!$omp end parallel workshare
il = 6 + ol
end if
if ( bb(2) == 0 ) then
ir = ni - gsize
else
or = ni - offset(2)
!$omp parallel workshare
dvar(or-4,:,:) = - ( q(1,5) * var(or,:,:) + q(2,5) * var(or-1,:,:) + &
q(3,5) * var(or-2,:,:) + q(4,5) * var(or-3,:,:) + &
q(5,5) * var(or-4,:,:) + q(6,5) * var(or-5,:,:) + &
q(7,5) * var(or-6,:,:) ) * idel
dvar(or-3,:,:) = - ( q(1,4) * var(or,:,:) + q(2,4) * var(or-1,:,:) + &
q(3,4) * var(or-2,:,:) + q(4,4) * var(or-3,:,:) + &
q(5,4) * var(or-4,:,:) + q(6,4) * var(or-5,:,:) + &
q(7,4) * var(or-6,:,:) ) * idel
dvar(or-2,:,:) = - ( q(1,3) * var(or,:,:) + q(2,3) * var(or-1,:,:) + &
q(3,3) * var(or-2,:,:) + q(4,3) * var(or-3,:,:) + &
q(5,3) * var(or-4,:,:) + &
q(6,3) * var(or-5,:,:) ) * idel
dvar(or-1,:,:) = - ( q(1,2) * var(or,:,:) + q(2,2) * var(or-1,:,:) + &
q(3,2) * var(or-2,:,:) + q(4,2) * var(or-3,:,:) + &
q(5,2) * var(or-4,:,:) + &
q(6,2) * var(or-5,:,:) ) * idel
dvar(or,:,:) = - ( q(1,1) * var(or,:,:) + q(2,1) * var(or-1,:,:) + &
q(3,1) * var(or-2,:,:) + &
q(4,1) * var(or-3,:,:) ) * idel
!$omp end parallel workshare
ir = or - 5
end if
if (il > ir+1) call CCTK_WARN (0, "domain too small")
!$omp parallel workshare
dvar(il:ir,:,:) = ( a(1) * ( var(il+1:ir+1,:,:) - &
var(il-1:ir-1,:,:) ) + &
a(2) * ( var(il+2:ir+2,:,:) - &
var(il-2:ir-2,:,:) ) ) * idel
!$omp end parallel workshare
case (1) direction
if ( zero_derivs_y /= 0 ) then
dvar = zero
else
if ( bb(1) == 0 ) then
jl = 1 + gsize
else
ol = offset(1)
!$omp parallel workshare
dvar(:,1+ol,:) = ( q(1,1) * var(:,1+ol,:) + q(2,1) * var(:,2+ol,:) + &
q(3,1) * var(:,3+ol,:) + q(4,1) * var(:,4+ol,:) ) * idel
dvar(:,2+ol,:) = ( q(1,2) * var(:,1+ol,:) + q(2,2) * var(:,2+ol,:) + &
q(3,2) * var(:,3+ol,:) + q(4,2) * var(:,4+ol,:) + &
q(5,2) * var(:,5+ol,:) + q(6,2) * var(:,6+ol,:) ) * idel
dvar(:,3+ol,:) = ( q(1,3) * var(:,1+ol,:) + q(2,3) * var(:,2+ol,:) + &
q(3,3) * var(:,3+ol,:) + q(4,3) * var(:,4+ol,:) + &
q(5,3) * var(:,5+ol,:) + q(6,3) * var(:,6+ol,:) ) * idel
dvar(:,4+ol,:) = ( q(1,4) * var(:,1+ol,:) + q(2,4) * var(:,2+ol,:) + &
q(3,4) * var(:,3+ol,:) + q(4,4) * var(:,4+ol,:) + &
q(5,4) * var(:,5+ol,:) + q(6,4) * var(:,6+ol,:) + &
q(7,4) * var(:,7+ol,:) ) * idel
dvar(:,5+ol,:) = ( q(1,5) * var(:,1+ol,:) + q(2,5) * var(:,2+ol,:) + &
q(3,5) * var(:,3+ol,:) + q(4,5) * var(:,4+ol,:) + &
q(5,5) * var(:,5+ol,:) + q(6,5) * var(:,6+ol,:) + &
q(7,5) * var(:,7+ol,:) ) * idel
!$omp end parallel workshare
jl = 6 + ol
end if
if ( bb(2) == 0 ) then
jr = nj - gsize
else
or = nj - offset(2)
!$omp parallel workshare
dvar(:,or-4,:) = - ( q(1,5) * var(:,or,:) + q(2,5) * var(:,or-1,:) + &
q(3,5) * var(:,or-2,:) + q(4,5) * var(:,or-3,:) + &
q(5,5) * var(:,or-4,:) + q(6,5) * var(:,or-5,:) + &
q(7,5) * var(:,or-6,:) ) * idel
dvar(:,or-3,:) = - ( q(1,4) * var(:,or,:) + q(2,4) * var(:,or-1,:) + &
q(3,4) * var(:,or-2,:) + q(4,4) * var(:,or-3,:) + &
q(5,4) * var(:,or-4,:) + q(6,4) * var(:,or-5,:) + &
q(7,4) * var(:,or-6,:) ) * idel
dvar(:,or-2,:) = - ( q(1,3) * var(:,or,:) + q(2,3) * var(:,or-1,:) + &
q(3,3) * var(:,or-2,:) + q(4,3) * var(:,or-3,:) + &
q(5,3) * var(:,or-4,:) + &
q(6,3) * var(:,or-5,:) ) * idel
dvar(:,or-1,:) = - ( q(1,2) * var(:,or,:) + q(2,2) * var(:,or-1,:) + &
q(3,2) * var(:,or-2,:) + q(4,2) * var(:,or-3,:) + &
q(5,2) * var(:,or-4,:) + &
q(6,2) * var(:,or-5,:) ) * idel
dvar(:,or,:) = - ( q(1,1) * var(:,or,:) + q(2,1) * var(:,or-1,:) + &
q(3,1) * var(:,or-2,:) + &
q(4,1) * var(:,or-3,:) ) * idel
!$omp end parallel workshare
jr = or - 5
end if
if (jl > jr+1) call CCTK_WARN (0, "domain too small")
!$omp parallel workshare
dvar(:,jl:jr,:) = ( a(1) * ( var(:,jl+1:jr+1,:) - &
var(:,jl-1:jr-1,:) ) + &
a(2) * ( var(:,jl+2:jr+2,:) - &
var(:,jl-2:jr-2,:) ) ) * idel
!$omp end parallel workshare
end if
case (2) direction
if ( zero_derivs_z /= 0 ) then
dvar = zero
else
if ( bb(1) == 0 ) then
kl = 1 + gsize
else
ol = offset(1)
!$omp parallel workshare
dvar(:,:,1+ol) = ( q(1,1) * var(:,:,1+ol) + q(2,1) * var(:,:,2+ol) + &
q(3,1) * var(:,:,3+ol) + q(4,1) * var(:,:,4+ol) ) * idel
dvar(:,:,2+ol) = ( q(1,2) * var(:,:,1+ol) + q(2,2) * var(:,:,2+ol) + &
q(3,2) * var(:,:,3+ol) + q(4,2) * var(:,:,4+ol) + &
q(5,2) * var(:,:,5+ol) + q(6,2) * var(:,:,6+ol) ) * idel
dvar(:,:,3+ol) = ( q(1,3) * var(:,:,1+ol) + q(2,3) * var(:,:,2+ol) + &
q(3,3) * var(:,:,3+ol) + q(4,3) * var(:,:,4+ol) + &
q(5,3) * var(:,:,5+ol) + q(6,3) * var(:,:,6+ol) ) * idel
dvar(:,:,4+ol) = ( q(1,4) * var(:,:,1+ol) + q(2,4) * var(:,:,2+ol) + &
q(3,4) * var(:,:,3+ol) + q(4,4) * var(:,:,4+ol) + &
q(5,4) * var(:,:,5+ol) + q(6,4) * var(:,:,6+ol) + &
q(7,4) * var(:,:,7+ol) ) * idel
dvar(:,:,5+ol) = ( q(1,5) * var(:,:,1+ol) + q(2,5) * var(:,:,2+ol) + &
q(3,5) * var(:,:,3+ol) + q(4,5) * var(:,:,4+ol) + &
q(5,5) * var(:,:,5+ol) + q(6,5) * var(:,:,6+ol) + &
q(7,5) * var(:,:,7+ol) ) * idel
!$omp end parallel workshare
kl = 6 + ol
end if
if ( bb(2) == 0 ) then
kr = nk - gsize
else
or = nk - offset(2)
!$omp parallel workshare
dvar(:,:,or-4) = - ( q(1,5) * var(:,:,or) + q(2,5) * var(:,:,or-1) + &
q(3,5) * var(:,:,or-2) + q(4,5) * var(:,:,or-3) + &
q(5,5) * var(:,:,or-4) + q(6,5) * var(:,:,or-5) + &
q(7,5) * var(:,:,or-6) ) * idel
dvar(:,:,or-3) = - ( q(1,4) * var(:,:,or) + q(2,4) * var(:,:,or-1) + &
q(3,4) * var(:,:,or-2) + q(4,4) * var(:,:,or-3) + &
q(5,4) * var(:,:,or-4) + q(6,4) * var(:,:,or-5) + &
q(7,4) * var(:,:,or-6) ) * idel
dvar(:,:,or-2) = - ( q(1,3) * var(:,:,or) + q(2,3) * var(:,:,or-1) + &
q(3,3) * var(:,:,or-2) + q(4,3) * var(:,:,or-3) + &
q(5,3) * var(:,:,or-4) + &
q(6,3) * var(:,:,or-5) ) * idel
dvar(:,:,or-1) = - ( q(1,2) * var(:,:,or) + q(2,2) * var(:,:,or-1) + &
q(3,2) * var(:,:,or-2) + q(4,2) * var(:,:,or-3) + &
q(5,2) * var(:,:,or-4) + &
q(6,2) * var(:,:,or-5) ) * idel
dvar(:,:,or) = - ( q(1,1) * var(:,:,or) + q(2,1) * var(:,:,or-1) + &
q(3,1) * var(:,:,or-2) + &
q(4,1) * var(:,:,or-3) ) * idel
!$omp end parallel workshare
kr = or - 5
end if
if (kl > kr+1) call CCTK_WARN (0, "domain too small")
!$omp parallel workshare
dvar(:,:,kl:kr) = ( a(1) * ( var(:,:,kl+1:kr+1) - &
var(:,:,kl-1:kr-1) ) + &
a(2) * ( var(:,:,kl+2:kr+2) - &
var(:,:,kl-2:kr-2) ) ) * idel
!$omp end parallel workshare
end if
end select direction
end subroutine deriv_gf_4_3
|