aboutsummaryrefslogtreecommitdiff
path: root/src/Derivatives_4_2.F90
blob: fb6888fb0c27717cf158aa74bd738e3dad8e0581 (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
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
#include "cctk.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"

subroutine deriv_gf_4_2 ( 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(2), save :: a
  CCTK_REAL, dimension(6,4), save :: q 
  CCTK_REAL :: idel

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

  logical, save :: first = .true.

  if ( first ) then
    a(1) = 2.0_wp/3.0_wp; a(2) = -1.0_wp/12.0_wp

    q(1,1) = -24.0_wp/17.0_wp; q(2,1) = 59.0_wp/34.0_wp
    q(3,1) = -4.0_wp/17.0_wp; q(4,1) = -3.0_wp/34.0_wp
    q(5,1) = zero; q(6,1) = zero
    q(1,2) = -1.0_wp/2.0_wp; q(2,2) = zero
    q(3,2) = 1.0_wp/2.0_wp; q(4,2) = zero
    q(5,2) = zero; q(6,2) = zero
    q(1,3) = 4.0_wp/43.0_wp; q(2,3) = -59.0_wp/86.0_wp
    q(3,3) = zero; q(4,3) = 59.0_wp/86.0_wp
    q(5,3) = -4.0_wp/43.0_wp; q(6,3) = zero
    q(1,4) = 3.0_wp/98.0_wp; q(2,4) = zero
    q(3,4) = -59.0_wp/98.0_wp; q(4,4) = zero
    q(5,4) = 32.0_wp/49.0_wp; q(6,4) = -4.0_wp/49.0_wp
    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
      dvar(1,:,:) = ( q(1,1) * var(1,:,:) + q(2,1) * var(2,:,:)  + &
                      q(3,1) * var(3,:,:) + q(4,1) * var(4,:,:) ) * idel
      dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(3,2) * var(3,:,:) ) * idel
      dvar(3,:,:) = ( q(1,3) * var(1,:,:) + q(2,3) * var(2,:,:) + &
                      q(4,3) * var(4,:,:) + q(5,3) * var(5,:,:) ) * idel
      dvar(4,:,:) = ( q(1,4) * var(1,:,:) + q(3,4) * var(3,:,:) + &
                      q(5,4) * var(5,:,:) + q(6,4) * var(6,:,:) ) * idel
      il = 5
    end if
    if ( bb(2) == 0 ) then
      ir = ni - gsize
    else
      dvar(ni-3,:,:) = - ( q(1,4) * var(ni,:,:) + &
                           q(3,4) * var(ni-2,:,:) + &
                           q(5,4) * var(ni-4,:,:) + &
                           q(6,4) * var(ni-5,:,:) ) * idel
      dvar(ni-2,:,:) = - ( q(1,3) * var(ni,:,:) + &
                           q(2,3) * var(ni-1,:,:) + &
                           q(4,3) * var(ni-3,:,:) + &
                           q(5,3) * var(ni-4,:,:) ) * idel
      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,:,:) + &
                           q(3,1) * var(ni-2,:,:) + &
                           q(4,1) * var(ni-3,:,:) ) * idel
      ir = ni - 4
    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,:,:) ) + &
                        a(2) * ( var(il+2:ir+2,:,:) - &
                                 var(il-2:ir-2,:,:) ) ) * 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,:)  + &
                        q(3,1) * var(:,3,:) + q(4,1) * var(:,4,:) ) * idel
        dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(3,2) * var(:,3,:) ) * idel
        dvar(:,3,:) = ( q(1,3) * var(:,1,:) + q(2,3) * var(:,2,:) + &
                        q(4,3) * var(:,4,:) + q(5,3) * var(:,5,:) ) * idel
        dvar(:,4,:) = ( q(1,4) * var(:,1,:) + q(3,4) * var(:,3,:) + &
                        q(5,4) * var(:,5,:) + q(6,4) * var(:,6,:) ) * idel
        jl = 5
      end if
      if ( bb(2) == 0 ) then
        jr = nj - gsize
      else
        dvar(:,nj-3,:) = - ( q(1,4) * var(:,nj,:) + &
                             q(3,4) * var(:,nj-2,:) + &
                             q(5,4) * var(:,nj-4,:) + &
                             q(6,4) * var(:,nj-5,:) ) * idel
        dvar(:,nj-2,:) = - ( q(1,3) * var(:,nj,:) + &
                             q(2,3) * var(:,nj-1,:) + &
                             q(4,3) * var(:,nj-3,:) + &
                             q(5,3) * var(:,nj-4,:) ) * idel
        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,:) + &
                             q(3,1) * var(:,nj-2,:) + &
                             q(4,1) * var(:,nj-3,:) ) * idel
        jr = nj - 4
      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,:) ) + &
                          a(2) * ( var(:,jl+2:jr+2,:) - &
                                   var(:,jl-2:jr-2,:) ) ) * 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)  + &
                        q(3,1) * var(:,:,3) + q(4,1) * var(:,:,4) ) * idel
        dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(3,2) * var(:,:,3) ) * idel
        dvar(:,:,3) = ( q(1,3) * var(:,:,1) + q(2,3) * var(:,:,2) + &
                        q(4,3) * var(:,:,4) + q(5,3) * var(:,:,5) ) * idel
        dvar(:,:,4) = ( q(1,4) * var(:,:,1) + q(3,4) * var(:,:,3) + &
                        q(5,4) * var(:,:,5) + q(6,4) * var(:,:,6) ) * idel
        kl = 5
      end if
      if ( bb(2) == 0 ) then
        kr = nk - gsize
      else 
        dvar(:,:,nk-3) = - ( q(1,4) * var(:,:,nk) + &
                             q(3,4) * var(:,:,nk-2) + &
                             q(5,4) * var(:,:,nk-4) + &
                             q(6,4) * var(:,:,nk-5) ) * idel
        dvar(:,:,nk-2) = - ( q(1,3) * var(:,:,nk) + &
                             q(2,3) * var(:,:,nk-1) + &
                             q(4,3) * var(:,:,nk-3) + &
                             q(5,3) * var(:,:,nk-4) ) * idel
        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) + &
                             q(3,1) * var(:,:,nk-2) + &
                             q(4,1) * var(:,:,nk-3) ) * idel
        kr = nk - 4
      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) ) + &
                          a(2) * ( var(:,:,kl+2:kr+2) - &
                                   var(:,:,kl-2:kr-2) ) ) * idel
    end if
  end select direction
end subroutine deriv_gf_4_2

subroutine up_deriv_gf_4_2 ( 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(-2:2), save :: a1, a2
  CCTK_REAL, dimension(6,4), save :: q1, q2
  CCTK_REAL :: idel

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

  logical, save :: first = .true.

  if ( first ) then
    a1(-2) = 0.1666666666666666666666666666666666666667_wp
    a1(-1) = -1.000000000000000000000000000000000000000_wp
    a1(0) = 0.5000000000000000000000000000000000000000_wp
    a1(1) = 0.3333333333333333333333333333333333333333_wp
    a1(2) = 0.0_wp

    q1(1,1) = -1.176470588235294117647058823529411764706_wp
    q1(2,1) = 1.264705882352941176470588235294117647059_wp
    q1(3,1) = 0.0_wp
    q1(4,1) = -0.08823529411764705882352941176470588235294_wp
    q1(5,1) = 0.0_wp
    q1(6,1) = 0.0_wp
    q1(1,2) = -0.6355932203389830508474576271186440677966_wp
    q1(2,2) = 0.3389830508474576271186440677966101694915_wp
    q1(3,2) = 0.2288135593220338983050847457627118644068_wp
    q1(4,2) = 0.06779661016949152542372881355932203389831_wp
    q1(5,2) = 0.0_wp
    q1(6,2) = 0.0_wp
    q1(1,3) = 0.1860465116279069767441860465116279069767_wp
    q1(2,3) = -1.058139534883720930232558139534883720930_wp
    q1(3,3) = 0.5581395348837209302325581395348837209302_wp
    q1(4,3) = 0.3139534883720930232558139534883720930233_wp
    q1(5,3) = 0.0_wp
    q1(6,3) = 0.0_wp
    q1(1,4) = 0.03061224489795918367346938775510204081633_wp
    q1(2,4) = 0.08163265306122448979591836734693877551020_wp
    q1(3,4) = -0.9285714285714285714285714285714285714286_wp
    q1(4,4) = 0.4897959183673469387755102040816326530612_wp
    q1(5,4) = 0.3265306122448979591836734693877551020408_wp
    q1(6,4) = 0.0_wp

    a2(-2) = 0.0_wp
    a2(-1) = -0.3333333333333333333333333333333333333333_wp
    a2(0) = -0.5000000000000000000000000000000000000000_wp
    a2(1) = 1.000000000000000000000000000000000000000_wp
    a2(2) = -0.1666666666666666666666666666666666666667_wp

    q2(1,1) = -1.647058823529411764705882352941176470588_wp
    q2(2,1) = 2.205882352941176470588235294117647058824_wp
    q2(3,1) = -0.4705882352941176470588235294117647058824_wp
    q2(4,1) = -0.08823529411764705882352941176470588235294_wp
    q2(5,1) = 0.0_wp
    q2(6,1) = 0.0_wp
    q2(1,2) = -0.3644067796610169491525423728813559322034_wp
    q2(2,2) = -0.3389830508474576271186440677966101694915_wp
    q2(3,2) = 0.7711864406779661016949152542372881355932_wp
    q2(4,2) = -0.06779661016949152542372881355932203389831_wp
    q2(5,2) = 0.0_wp
    q2(6,2) = 0.0_wp
    q2(1,3) = 0.0_wp
    q2(2,3) = -0.3139534883720930232558139534883720930233_wp
    q2(3,3) = -0.5581395348837209302325581395348837209302_wp
    q2(4,3) = 1.058139534883720930232558139534883720930_wp
    q2(5,3) = -0.1860465116279069767441860465116279069767_wp
    q2(6,3) = 0.0_wp
    q2(1,4) = 0.03061224489795918367346938775510204081633_wp
    q2(2,4) = -0.08163265306122448979591836734693877551020_wp
    q2(3,4) = -0.2755102040816326530612244897959183673469_wp
    q2(4,4) = -0.4897959183673469387755102040816326530612_wp
    q2(5,4) = 0.9795918367346938775510204081632653061224_wp
    q2(6,4) = -0.1632653061224489795918367346938775510204_wp

    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
      where ( up(1,:,:) < zero )
        dvar(1,:,:) = ( q1(1,1) * var(1,:,:) + q1(2,1) * var(2,:,:) + &
                        q1(3,1) * var(3,:,:) + q1(4,1) * var(4,:,:) ) * idel
      elsewhere
        dvar(1,:,:) = ( q2(1,1) * var(1,:,:) + q2(2,1) * var(2,:,:) + &
                        q2(3,1) * var(3,:,:) + q2(4,1) * var(4,:,:) ) * idel
      end where
      where ( up(2,:,:) < zero )
        dvar(2,:,:) = ( q1(1,2) * var(1,:,:) + q1(2,2) * var(2,:,:) + &
                        q1(3,2) * var(3,:,:) + q1(4,2) * var(4,:,:) ) * idel
      elsewhere
        dvar(2,:,:) = ( q2(1,2) * var(1,:,:) + q2(2,2) * var(2,:,:) + &
                        q2(3,2) * var(3,:,:) + q2(4,2) * var(4,:,:) ) * idel
      end where
      where ( up(3,:,:) < zero ) 
        dvar(3,:,:) = ( q1(1,3) * var(1,:,:) + q1(2,3) * var(2,:,:) + &
                        q1(3,3) * var(3,:,:) + q1(4,3) * var(4,:,:) + &
                        q1(5,3) * var(5,:,:) ) * idel
      elsewhere
        dvar(3,:,:) = ( q2(1,3) * var(1,:,:) + q2(2,3) * var(2,:,:) + &
                        q2(3,3) * var(3,:,:) + q2(4,3) * var(4,:,:) + &
                        q2(5,3) * var(5,:,:) ) * idel
      end where
      where ( up(4,:,:) < zero ) 
        dvar(4,:,:) = ( q1(1,4) * var(1,:,:) + q1(2,4) * var(2,:,:) + &
                        q1(3,4) * var(3,:,:) + q1(4,4) * var(4,:,:) + &
                        q1(5,4) * var(5,:,:) + q1(6,4) * var(6,:,:) ) * idel
      elsewhere
        dvar(4,:,:) = ( q2(1,4) * var(1,:,:) + q2(2,4) * var(2,:,:) + &
                        q2(3,4) * var(3,:,:) + q2(4,4) * var(4,:,:) + &
                        q2(5,4) * var(5,:,:) + q2(6,4) * var(6,:,:) ) * idel
      end where
      il = 5
    end if
    if ( bb(2) == 0 ) then
      ir = ni - gsize
    else
      where ( up(ni-3,:,:) < zero ) 
        dvar(ni-3,:,:) = - ( q2(1,4) * var(ni,:,:) + &
                             q2(2,4) * var(ni-1,:,:) + &
                             q2(3,4) * var(ni-2,:,:) + &
                             q2(4,4) * var(ni-3,:,:) + &
                             q2(5,4) * var(ni-4,:,:) + &
                             q2(6,4) * var(ni-5,:,:) ) * idel
      elsewhere
        dvar(ni-3,:,:) = - ( q1(1,4) * var(ni,:,:) + &
                             q1(2,4) * var(ni-1,:,:) + &
                             q1(3,4) * var(ni-2,:,:) + &
                             q1(4,4) * var(ni-3,:,:) + &
                             q1(5,4) * var(ni-4,:,:) + &
                             q1(6,4) * var(ni-5,:,:) ) * idel
      end where
      where ( up(ni-2,:,:) < zero ) 
        dvar(ni-2,:,:) = - ( q2(1,3) * var(ni,:,:) + &
                             q2(2,3) * var(ni-1,:,:) + &
                             q2(3,3) * var(ni-2,:,:) + &
                             q2(4,3) * var(ni-3,:,:) + &
                             q2(5,3) * var(ni-4,:,:) ) * idel
      elsewhere
        dvar(ni-2,:,:) = - ( q1(1,3) * var(ni,:,:) + &
                             q1(2,3) * var(ni-1,:,:) + &
                             q1(3,3) * var(ni-2,:,:) + &
                             q1(4,3) * var(ni-3,:,:) + &
                             q1(5,3) * var(ni-4,:,:) ) * idel
      end where
      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,:,:) + &
                             q2(4,2) * var(ni-3,:,:) ) * idel
      elsewhere
        dvar(ni-1,:,:) = - ( q1(1,2) * var(ni,:,:) + &
                             q1(2,2) * var(ni-1,:,:) + &
                             q1(3,2) * var(ni-2,:,:) + &
                             q1(4,2) * var(ni-3,:,:) ) * idel
      end where
      where ( up(ni,:,:) < zero ) 
        dvar(ni,:,:) = - ( q2(1,1) * var(ni,:,:) + &
                           q2(2,1) * var(ni-1,:,:) + &
                           q2(3,1) * var(ni-2,:,:) + &
                           q2(4,1) * var(ni-3,:,:) ) * idel
      elsewhere
        dvar(ni,:,:) = - ( q1(1,1) * var(ni,:,:) + &
                           q1(2,1) * var(ni-1,:,:) + &
                           q1(3,1) * var(ni-2,:,:) + &
                           q1(4,1) * var(ni-3,:,:) ) * idel
      end where
      ir = ni - 4
    end if
    if (il > ir+1) call CCTK_WARN (0, "domain too small")
    where ( up(il:ir,:,:) < zero ) 
      dvar(il:ir,:,:) = ( a1(-2) * var(il-2:ir-2,:,:) + &
                          a1(-1) * var(il-1:ir-1,:,:) + &
                          a1(0)  * var(il:ir,:,:) + &
                          a1(1)  * var(il+1:ir+1,:,:) + &
                          a1(2)  * var(il+2:ir+2,:,:) ) * idel
    elsewhere
      dvar(il:ir,:,:) = ( a2(-2) * var(il-2:ir-2,:,:) + &
                          a2(-1) * var(il-1:ir-1,:,:) + &
                          a2(0)  * var(il:ir,:,:) + &
                          a2(1)  * var(il+1:ir+1,:,:) + &
                          a2(2)  * var(il+2:ir+2,:,:) ) * 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,:) + &
                          q1(3,1) * var(:,3,:) + q1(4,1) * var(:,4,:) ) * idel
        elsewhere
          dvar(:,1,:) = ( q2(1,1) * var(:,1,:) + q2(2,1) * var(:,2,:) + &
                          q2(3,1) * var(:,3,:) + q2(4,1) * var(:,4,:) ) * idel
        end where
        where ( up(:,2,:) < zero )
          dvar(:,2,:) = ( q1(1,2) * var(:,1,:) + q1(2,2) * var(:,2,:) + &
                          q1(3,2) * var(:,3,:) + q1(4,2) * var(:,4,:) ) * idel
        elsewhere
          dvar(:,2,:) = ( q2(1,2) * var(:,1,:) + q2(2,2) * var(:,2,:) + &
                          q2(3,2) * var(:,3,:) + q2(4,2) * var(:,4,:) ) * idel
        end where
        where ( up(:,3,:) < zero )
          dvar(:,3,:) = ( q1(1,3) * var(:,1,:) + q1(2,3) * var(:,2,:) + &
                          q1(3,3) * var(:,3,:) + q1(4,3) * var(:,4,:) + &
                          q1(5,3) * var(:,5,:) ) * idel
        elsewhere
          dvar(:,3,:) = ( q2(1,3) * var(:,1,:) + q2(2,3) * var(:,2,:) + &
                          q2(3,3) * var(:,3,:) + q2(4,3) * var(:,4,:) + &
                          q2(5,3) * var(:,5,:) ) * idel
        end where
        where ( up(:,4,:) < zero )
          dvar(:,4,:) = ( q1(1,4) * var(:,1,:) + q1(2,4) * var(:,2,:) + &
                          q1(3,4) * var(:,3,:) + q1(4,4) * var(:,4,:) + &
                          q1(5,4) * var(:,5,:) + q1(6,4) * var(:,6,:) ) * idel
        elsewhere
          dvar(:,4,:) = ( q2(1,4) * var(:,1,:) + q2(2,4) * var(:,2,:) + &
                          q2(3,4) * var(:,3,:) + q2(4,4) * var(:,4,:) + &
                          q2(5,4) * var(:,5,:) + q2(6,4) * var(:,6,:) ) * idel
        end where
        jl = 5
      end if
      if ( bb(2) == 0 ) then
        jr = nj - gsize
      else
        if (jl > jr+1) call CCTK_WARN (0, "domain too small")
        where ( up(:,nj-3,:) < zero )
          dvar(:,nj-3,:) = - ( q2(1,4) * var(:,nj,:) + &
                               q2(2,4) * var(:,nj-1,:) + &
                               q2(3,4) * var(:,nj-2,:) + &
                               q2(4,4) * var(:,nj-3,:) + &
                               q2(5,4) * var(:,nj-4,:) + &
                               q2(6,4) * var(:,nj-5,:) ) * idel
        elsewhere
          dvar(:,nj-3,:) = - ( q1(1,4) * var(:,nj,:) + &
                               q1(2,4) * var(:,nj-1,:) + &
                               q1(3,4) * var(:,nj-2,:) + &
                               q1(4,4) * var(:,nj-3,:) + &
                               q1(5,4) * var(:,nj-4,:) + &
                               q1(6,4) * var(:,nj-5,:) ) * idel
        end where
        where ( up(:,nj-2,:) < zero )
          dvar(:,nj-2,:) = - ( q2(1,3) * var(:,nj,:) + &
                               q2(2,3) * var(:,nj-1,:) + &
                               q2(3,3) * var(:,nj-2,:) + &
                               q2(4,3) * var(:,nj-3,:) + &
                               q2(5,3) * var(:,nj-4,:) ) * idel
        elsewhere
          dvar(:,nj-2,:) = - ( q1(1,3) * var(:,nj,:) + &
                               q1(2,3) * var(:,nj-1,:) + &
                               q1(3,3) * var(:,nj-2,:) + &
                               q1(4,3) * var(:,nj-3,:) + &
                               q1(5,3) * var(:,nj-4,:) ) * idel
        end where
        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,:) + &
                               q2(4,2) * var(:,nj-3,:) ) * idel
        elsewhere
          dvar(:,nj-1,:) = - ( q1(1,2) * var(:,nj,:) + &
                               q1(2,2) * var(:,nj-1,:) + &
                               q1(3,2) * var(:,nj-2,:) + &
                               q1(4,2) * var(:,nj-3,:) ) * idel
        end where
        where ( up(:,nj,:) < zero )
          dvar(:,nj,:) = - ( q2(1,1) * var(:,nj,:) + &
                             q2(2,1) * var(:,nj-1,:) + &
                             q2(3,1) * var(:,nj-2,:) + &
                             q2(4,1) * var(:,nj-3,:) ) * idel
        elsewhere
          dvar(:,nj,:) = - ( q1(1,1) * var(:,nj,:) + &
                             q1(2,1) * var(:,nj-1,:) + &
                             q1(3,1) * var(:,nj-2,:) + &
                             q1(4,1) * var(:,nj-3,:) ) * idel
        end where
        jr = nj - 4
      end if
      where ( up(:,jl:jr,:) < zero )
        dvar(:,jl:jr,:) = ( a1(-2) * var(:,jl-2:jr-2,:) + &
                            a1(-1) * var(:,jl-1:jr-1,:) + &
                            a1(0)  * var(:,jl:jr,:) + &
                            a1(1)  * var(:,jl+1:jr+1,:) + &
                            a1(2)  * var(:,jl+2:jr+2,:) ) * idel
      elsewhere
        dvar(:,jl:jr,:) = ( a2(-2) * var(:,jl-2:jr-2,:) + &
                            a2(-1) * var(:,jl-1:jr-1,:) + &
                            a2(0)  * var(:,jl:jr,:) + &
                            a2(1)  * var(:,jl+1:jr+1,:) + &
                            a2(2)  * var(:,jl+2:jr+2,:) ) * 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) + &
                          q1(3,1) * var(:,:,3) + q1(4,1) * var(:,:,4) ) * idel
        elsewhere
          dvar(:,:,1) = ( q2(1,1) * var(:,:,1) + q2(2,1) * var(:,:,2) + &
                          q2(3,1) * var(:,:,3) + q2(4,1) * var(:,:,4) ) * idel
        end where
        where ( up(:,:,2) < zero ) 
          dvar(:,:,2) = ( q1(1,2) * var(:,:,1) + q1(2,2) * var(:,:,2) + &
                          q1(3,2) * var(:,:,3) + q1(4,2) * var(:,:,4) ) * idel
        elsewhere
          dvar(:,:,2) = ( q2(1,2) * var(:,:,1) + q2(2,2) * var(:,:,2) + &
                          q2(3,2) * var(:,:,3) + q2(4,2) * var(:,:,4) ) * idel
        end where
        where ( up(:,:,3) < zero ) 
          dvar(:,:,3) = ( q1(1,3) * var(:,:,1) + q1(2,3) * var(:,:,2) + &
                          q1(3,3) * var(:,:,3) + q1(4,3) * var(:,:,4) + &
                          q1(5,3) * var(:,:,5) ) * idel
        elsewhere
          dvar(:,:,3) = ( q2(1,3) * var(:,:,1) + q2(2,3) * var(:,:,2) + &
                          q2(3,3) * var(:,:,3) + q2(4,3) * var(:,:,4) + &
                          q2(5,3) * var(:,:,5) ) * idel
        end where
        where ( up(:,:,4) < zero ) 
          dvar(:,:,4) = ( q1(1,4) * var(:,:,1) + q1(2,4) * var(:,:,2) + &
                          q1(3,4) * var(:,:,3) + q1(4,4) * var(:,:,4) + &
                          q1(5,4) * var(:,:,5) + q1(6,4) * var(:,:,6) ) * idel
        elsewhere
          dvar(:,:,4) = ( q2(1,4) * var(:,:,1) + q2(2,4) * var(:,:,2) + &
                          q2(3,4) * var(:,:,3) + q2(4,4) * var(:,:,4) + &
                          q2(5,4) * var(:,:,5) + q2(6,4) * var(:,:,6) ) * idel
        end where
        kl = 5
      end if
      if ( bb(2) == 0 ) then
        kr = nk - gsize
      else
        where ( up(:,:,nk-3) < zero ) 
          dvar(:,:,nk-3) = - ( q2(1,4) * var(:,:,nk) + &
                               q2(2,4) * var(:,:,nk-1) + &
                               q2(3,4) * var(:,:,nk-2) + &
                               q2(4,4) * var(:,:,nk-3) + &
                               q2(5,4) * var(:,:,nk-4) + &
                               q2(6,4) * var(:,:,nk-5) ) * idel
        elsewhere
          dvar(:,:,nk-3) = - ( q1(1,4) * var(:,:,nk) + &
                               q1(2,4) * var(:,:,nk-1) + &
                               q1(3,4) * var(:,:,nk-2) + &
                               q1(4,4) * var(:,:,nk-3) + &
                               q1(5,4) * var(:,:,nk-4) + &
                               q1(6,4) * var(:,:,nk-5) ) * idel
        end where
        where ( up(:,:,nk-2) < zero ) 
          dvar(:,:,nk-2) = - ( q2(1,3) * var(:,:,nk) + &
                               q2(2,3) * var(:,:,nk-1) + &
                               q2(3,3) * var(:,:,nk-2) + &
                               q2(4,3) * var(:,:,nk-3) + &
                               q2(5,3) * var(:,:,nk-4) ) * idel
        elsewhere
          dvar(:,:,nk-2) = - ( q1(1,3) * var(:,:,nk) + &
                               q1(2,3) * var(:,:,nk-1) + &
                               q1(3,3) * var(:,:,nk-2) + &
                               q1(4,3) * var(:,:,nk-3) + &
                               q1(5,3) * var(:,:,nk-4) ) * idel
        end where
        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) + &
                               q2(4,2) * var(:,:,nk-3) ) * idel
        elsewhere
          dvar(:,:,nk-1) = - ( q1(1,2) * var(:,:,nk) + &
                               q1(2,2) * var(:,:,nk-1) + &
                               q1(3,2) * var(:,:,nk-2) + &
                               q1(4,2) * var(:,:,nk-3) ) * idel
        end where
        where ( up(:,:,nk) < zero ) 
          dvar(:,:,nk) = - ( q2(1,1) * var(:,:,nk) + &
                             q2(2,1) * var(:,:,nk-1) + &
                             q2(3,1) * var(:,:,nk-2) + &
                             q2(4,1) * var(:,:,nk-3) ) * idel
        elsewhere
          dvar(:,:,nk) = - ( q1(1,1) * var(:,:,nk) + &
                             q1(2,1) * var(:,:,nk-1) + &
                             q1(3,1) * var(:,:,nk-2) + &
                             q1(4,1) * var(:,:,nk-3) ) * idel
        end where
        kr = nk - 4
      end if
      if (kl > kr+1) call CCTK_WARN (0, "domain too small")
      where ( up(:,:,kl:kr) < zero )
        dvar(:,:,kl:kr) = ( a1(-2) * var(:,:,kl-2:kr-2) + &
                            a1(-1) * var(:,:,kl-1:kr-1) + &
                            a1(0)  * var(:,:,kl:kr) + &
                            a1(1)  * var(:,:,kl+1:kr+1) + &
                            a1(2)  * var(:,:,kl+2:kr+2) ) * idel
      elsewhere
        dvar(:,:,kl:kr) = ( a2(-2) * var(:,:,kl-2:kr-2) + &
                            a2(-1) * var(:,:,kl-1:kr-1) + &
                            a2(0)  * var(:,:,kl:kr) + &
                            a2(1)  * var(:,:,kl+1:kr+1) + &
                            a2(2)  * var(:,:,kl+2:kr+2) ) * idel
      end where
    end if
  end select direction
end subroutine up_deriv_gf_4_2