aboutsummaryrefslogtreecommitdiff
path: root/src/Derivatives_2_1.F90
blob: 788d3252af945748b18a969c36bbbf393ffb3e12 (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
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
#include "cctk.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"

subroutine deriv_gf_2_1 ( var, ni, nj, nk, dir, bb, gsize, delta, dvar )

  implicit none

  DECLARE_CCTK_FUNCTIONS
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL, parameter :: zero = 0.0
  integer, parameter :: wp = kind(zero)
  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_REAL, intent(IN) :: delta
  CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar

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

  CCTK_INT :: il, ir, jl, jr, kl, kr

  logical, save :: first = .true.

  if ( first ) then
    a(1) = 1.0_wp/2.0_wp

    q(1,1) = -1.0_wp; q(2,1) = 1.0_wp; q(3,1) = zero
    q(1,2) = -1.0_wp/2.0_wp; q(2,2) = zero; q(3,2) = 1.0_wp/2.0_wp
    first = .false.
  end if

  idel = 1.0_wp / delta

  if (gsize < 1) call CCTK_WARN (0, "not enough ghostzones")

  direction: select case (dir)
  case (0) direction
    if ( bb(1) == 0 ) then
      il = 1 + gsize
    else
      dvar(1,:,:) = ( q(1,1) * var(1,:,:) + q(2,1) * var(2,:,:) ) * idel
      dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(3,2) * var(3,:,:) ) * idel
      il = 3
    end if
    if ( bb(2) == 0 ) then
      ir = ni - gsize
    else
      dvar(ni-1,:,:) = - ( q(1,2) * var(ni,:,:) + &
                           q(3,2) * var(ni-2,:,:) ) * idel
      dvar(ni,:,:) = - ( q(1,1) * var(ni,:,:) + &
                         q(2,1) * var(ni-1,:,:) ) * idel
      ir = ni - 2
    end if
    if (il > ir+1) call CCTK_WARN (0, "domain too small")
    dvar(il:ir,:,:) = a(1) * ( var(il+1:ir+1,:,:) - var(il-1:ir-1,:,:) ) * idel
  case (1) direction
    if ( zero_derivs_y /= 0 ) then
      dvar = zero
    else
      if ( bb(1) == 0 ) then
        jl = 1 + gsize
      else
        dvar(:,1,:) = ( q(1,1) * var(:,1,:) + q(2,1) * var(:,2,:) ) * idel
        dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(3,2) * var(:,3,:) ) * idel
        jl = 3
      end if
      if ( bb(2) == 0 ) then
        jr = nj - gsize
      else
        dvar(:,nj-1,:) = - ( q(1,2) * var(:,nj,:) + &
                             q(3,2) * var(:,nj-2,:) ) * idel
        dvar(:,nj,:) = - ( q(1,1) * var(:,nj,:) + &
                           q(2,1) * var(:,nj-1,:) ) * idel
        jr = nj - 2
      end if
      if (jl > jr+1) call CCTK_WARN (0, "domain too small")
      dvar(:,jl:jr,:) = a(1) * ( var(:,jl+1:jr+1,:) - &
                                 var(:,jl-1:jr-1,:) ) * idel
    end if
  case (2) direction
    if ( zero_derivs_z /= 0 ) then
      dvar = zero
    else
      if ( bb(1) == 0 ) then
        kl = 1 + gsize
      else
        dvar(:,:,1) = ( q(1,1) * var(:,:,1) + q(2,1) * var(:,:,2) ) * idel
        dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(3,2) * var(:,:,3) ) * idel
        kl = 3
      end if
      if ( bb(2) == 0 ) then
        kr = nk - gsize
      else 
        dvar(:,:,nk-1) = - ( q(1,2) * var(:,:,nk) + &
                             q(3,2) * var(:,:,nk-2) ) * idel
        dvar(:,:,nk) = - ( q(1,1) * var(:,:,nk) + &
                           q(2,1) * var(:,:,nk-1) ) * idel
        kr = nk - 2
      end if
      if (kl > kr+1) call CCTK_WARN (0, "domain too small")
      dvar(:,:,kl:kr) = a(1) * ( var(:,:,kl+1:kr+1) -  &
                                 var(:,:,kl-1:kr-1) ) * idel
    end if
  end select direction
end subroutine deriv_gf_2_1

subroutine up_deriv_gf_2_1 ( var, ni, nj, nk, dir, bb, gsize, delta, up, dvar )

  implicit none

  DECLARE_CCTK_FUNCTIONS
  DECLARE_CCTK_PARAMETERS

  CCTK_REAL, parameter :: zero = 0.0
  integer, parameter :: wp = kind(zero)
  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_REAL, intent(IN) :: delta
  CCTK_REAL, dimension(ni,nj,nk), intent(IN) :: up
  CCTK_REAL, dimension(ni,nj,nk), intent(OUT) :: dvar

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

  CCTK_INT :: il, ir, jl, jr, kl, kr

  logical, save :: first = .true.

  if ( first ) then
    a1(-1) = -1.000000000000000000000000000000000000000_wp
    a1(0) = 1.000000000000000000000000000000000000000_wp
    a1(1) = 0.0_wp

    q1(1,1) = 0.0_wp
    q1(2,1) = 0.0_wp
    q1(3,1) = 0.0_wp
    q1(1,2) = -1.000000000000000000000000000000000000000_wp
    q1(2,2) = 1.000000000000000000000000000000000000000_wp
    q1(3,2) = 0.0_wp

    a2(-1) = 0.0_wp
    a2(0) = -1.000000000000000000000000000000000000000_wp
    a2(1) = 1.000000000000000000000000000000000000000_wp

    q2(1,1) = -2.000000000000000000000000000000000000000_wp
    q2(2,1) = 2.000000000000000000000000000000000000000_wp
    q2(3,1) = 0.0_wp
    q2(1,2) = 0.0_wp
    q2(2,2) = -1.000000000000000000000000000000000000000_wp
    q2(3,2) = 1.000000000000000000000000000000000000000_wp

    first = .false.
  end if

  idel = 1.0_wp / delta

  if (gsize < 1) call CCTK_WARN (0, "not enough ghostzones")

  direction: select case (dir)
  case (0) direction
    if ( bb(1) == 0 ) then
      il = 1 + gsize
    else
      where ( up(1,:,:) < zero )
        dvar(1,:,:) = ( q1(1,1) * var(1,:,:) + q1(2,1) * var(2,:,:) ) * idel
      elsewhere
        dvar(1,:,:) = ( q2(1,1) * var(1,:,:) + q2(2,1) * var(2,:,:) ) * idel
      end where
      where ( up(2,:,:) < zero ) 
        dvar(2,:,:) = ( q1(1,2) * var(1,:,:) + q1(2,2) * var(2,:,:) + &
                        q1(3,2) * var(3,:,:) ) * idel
      elsewhere
        dvar(2,:,:) = ( q2(1,2) * var(1,:,:) + q2(2,2) * var(2,:,:) + &
                        q2(3,2) * var(3,:,:) ) * idel
      end where
      il = 3
    end if
    if ( bb(2) == 0 ) then
      ir = ni - gsize
    else
      where ( up(ni-1,:,:) < zero )
        dvar(ni-1,:,:) = - ( q2(1,2) * var(ni,:,:) + &
                             q2(2,2) * var(ni-1,:,:) + &
                             q2(3,2) * var(ni-2,:,:) ) * idel
      elsewhere
        dvar(ni-1,:,:) = - ( q1(1,2) * var(ni,:,:) + &
                             q1(2,2) * var(ni-1,:,:) + &
                             q1(3,2) * var(ni-2,:,:) ) * idel
      end where
      where ( up(ni,:,:) < zero )
        dvar(ni,:,:) = - ( q2(1,1) * var(ni,:,:) + &
                           q2(2,1) * var(ni-1,:,:) ) * idel
      elsewhere
        dvar(ni,:,:) = - ( q1(1,1) * var(ni,:,:) + &
                           q1(2,1) * var(ni-1,:,:) ) * idel
      end where
      ir = ni - 2
    end if
    if (il > ir+1) call CCTK_WARN (0, "domain too small")
    where ( up(il:ir,:,:) < zero ) 
      dvar(il:ir,:,:) = ( a1(-1) * var(il-1:ir-1,:,:) + &
                          a1(0)  * var(il:ir,:,:) + &
                          a1(1)  * var(il+1:ir+1,:,:) ) * idel
    elsewhere
      dvar(il:ir,:,:) = ( a2(-1) * var(il-1:ir-1,:,:) + &
                          a2(0)  * var(il:ir,:,:) + &
                          a2(1)  * var(il+1:ir+1,:,:) ) * idel
    end where
  case (1) direction
    if ( zero_derivs_y /= 0 ) then
      dvar = zero
    else
      if ( bb(1) == 0 ) then
        jl = 1 + gsize
      else
        where ( up(:,1,:) < zero )
          dvar(:,1,:) = ( q1(1,1) * var(:,1,:) + q1(2,1) * var(:,2,:) ) * idel
        elsewhere
          dvar(:,1,:) = ( q2(1,1) * var(:,1,:) + q2(2,1) * var(:,2,:) ) * idel
        end where
        where ( up(:,2,:) < zero )
          dvar(:,2,:) = ( q1(1,2) * var(:,1,:) + q1(2,2) * var(:,2,:) + &
                          q1(3,2) * var(:,3,:) ) * idel
        elsewhere
          dvar(:,2,:) = ( q2(1,2) * var(:,1,:) + q2(2,2) * var(:,2,:) + &
                          q2(3,2) * var(:,3,:) ) * idel
        end where
        jl = 3
      end if
      if ( bb(2) == 0 ) then
        jr = nj - gsize
      else
        where ( up(:,nj-1,:) < zero )
          dvar(:,nj-1,:) = - ( q2(1,2) * var(:,nj,:) + &
                               q2(2,2) * var(:,nj-1,:) + &
                               q2(3,2) * var(:,nj-2,:) ) * idel
        elsewhere
          dvar(:,nj-1,:) = - ( q1(1,2) * var(:,nj,:) + &
                               q1(2,2) * var(:,nj-1,:) + &
                               q1(3,2) * var(:,nj-2,:) ) * idel
        end where
        where ( up(:,nj,:) < zero )
          dvar(:,nj,:) = - ( q2(1,1) * var(:,nj,:) + &
                             q2(2,1) * var(:,nj-1,:) ) * idel
        elsewhere
          dvar(:,nj,:) = - ( q1(1,1) * var(:,nj,:) + &
                             q1(2,1) * var(:,nj-1,:) ) * idel
        end where
        jr = nj - 2
      end if
      if (jl > jr+1) call CCTK_WARN (0, "domain too small")
      where ( up(:,jl:jr,:) < zero ) 
        dvar(:,jl:jr,:) = ( a1(-1) * var(:,jl-1:jr-1,:) + &
                            a1(0)  * var(:,jl:jr,:) + &
                            a1(1)  * var(:,jl+1:jr+1,:) ) * idel
      elsewhere
        dvar(:,jl:jr,:) = ( a2(-1) * var(:,jl-1:jr-1,:) + &
                            a2(0)  * var(:,jl:jr,:) + &
                            a2(1)  * var(:,jl+1:jr+1,:) ) * idel
      end where
    end if
  case (2) direction
    if ( zero_derivs_z /= 0 ) then
      dvar = zero
    else
      if ( bb(1) == 0 ) then
        kl = 1 + gsize
      else
        where ( up(:,:,1) < zero )
          dvar(:,:,1) = ( q1(1,1) * var(:,:,1) + q1(2,1) * var(:,:,2) ) * idel
        elsewhere
          dvar(:,:,1) = ( q2(1,1) * var(:,:,1) + q2(2,1) * var(:,:,2) ) * idel
        end where
        where ( up(:,:,2) < zero )  
          dvar(:,:,2) = ( q1(1,2) * var(:,:,1) + q1(2,2) * var(:,:,2) + &
                          q1(3,2) * var(:,:,3) ) * idel
        elsewhere
          dvar(:,:,2) = ( q2(1,2) * var(:,:,1) + q2(2,2) * var(:,:,2) + &
                          q2(3,2) * var(:,:,3) ) * idel
        end where
        kl = 3
      end if
      if ( bb(2) == 0 ) then
        kr = nk - gsize
      else 
        where ( up(:,:,nk-1) < zero )
          dvar(:,:,nk-1) = - ( q2(1,2) * var(:,:,nk) + &
                               q2(2,2) * var(:,:,nk-1) + &
                               q2(3,2) * var(:,:,nk-2) ) * idel
        elsewhere
          dvar(:,:,nk-1) = - ( q1(1,2) * var(:,:,nk) + &
                               q1(2,2) * var(:,:,nk-1) + &
                               q1(3,2) * var(:,:,nk-2) ) * idel
        end where
        where ( up(:,:,nk) < zero )
          dvar(:,:,nk) = - ( q2(1,1) * var(:,:,nk) + &
                             q2(2,1) * var(:,:,nk-1) ) * idel
        elsewhere
          dvar(:,:,nk) = - ( q1(1,1) * var(:,:,nk) + &
                             q1(2,1) * var(:,:,nk-1) ) * idel
        end where
        kr = nk - 2
      end if
      if (kl > kr+1) call CCTK_WARN (0, "domain too small")
      where ( up(:,:,kl:kr) < zero )
        dvar(:,:,kl:kr) = ( a1(-1) * var(:,:,kl-1:kr-1) + &
                            a1(0)  * var(:,:,kl:kr) + &
                            a1(1)  * var(:,:,kl+1:kr+1) ) * idel
      elsewhere
        dvar(:,:,kl:kr) = ( a2(-1) * var(:,:,kl-1:kr-1) + &
                            a2(0)  * var(:,:,kl:kr) + &
                            a2(1)  * var(:,:,kl+1:kr+1) ) * idel
      end where
    end if
  end select direction
end subroutine up_deriv_gf_2_1