aboutsummaryrefslogtreecommitdiff
path: root/src/Derivatives_6_3_min_err_coeff.F90
blob: b2d00e709e332bc3010ae30052761f84d24a510f (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
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
#include "cctk.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"

subroutine deriv_gf_6_3_opt ( 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(3), save :: a
  CCTK_REAL, dimension(9,6), save :: q 
  CCTK_REAL :: idel

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

  logical, save :: first = .true.

  if ( first ) then
    a(1) = 3.0_wp/4.0_wp; a(2) = -3.0_wp/20.0_wp; a(3) = 1.0_wp/60.0_wp

    q(1,1) = -1.582533518939116418785258993332844897062_wp
    q(2,1) = 2.033426786468126253898161347360808173712_wp
    q(3,1) = -0.1417052898146741610733887894481170575600_wp
    q(4,1) = -0.4501096599735708523162117824920488989702_wp
    q(5,1) = 0.1042956382142412661862395105494407610836_wp
    q(6,1) = 0.03662604404499391209045870736276191879693_wp
    q(7,1) = 0.0_wp
    q(8,1) = 0.0_wp
    q(9,1) = 0.0_wp
    q(1,2) = -0.4620701275035953590186631853846278325646_wp
    q(2,2) = 0.0_wp
    q(3,2) = 0.2873679417026202568532985205129449923126_wp
    q(4,2) = 0.2585974499280928196267362923074433487080_wp
    q(5,2) = -0.06894808744606961472005221923058251153103_wp
    q(6,2) = -0.01494717668104810274131940820517799692506_wp
    q(7,2) = 0.0_wp
    q(8,2) = 0.0_wp
    q(9,2) = 0.0_wp
    q(1,3) = 0.07134398748360337973038301686379010397038_wp
    q(2,3) = -0.6366933020423417826592908754928085932593_wp
    q(3,3) = 0.0_wp
    q(4,3) = 0.6067199374180168986519150843189505198519_wp
    q(5,3) = -0.02338660408468356531858175098561718651857_wp
    q(6,3) = -0.01798401877459493040442547470431484404443_wp
    q(7,3) = 0.0_wp
    q(8,3) = 0.0_wp
    q(9,3) = 0.0_wp
    q(1,4) = 0.1146397975178068401430112823144985150596_wp
    q(2,4) = -0.2898424301162697370942324201800071793273_wp
    q(3,4) = -0.3069262456316931913128086944558079603132_wp
    q(4,4) = 0.0_wp
    q(5,4) = 0.5203848121857539166740071338174418292578_wp
    q(6,4) = -0.05169127637022742348368508279860701098408_wp
    q(7,4) = 0.01343534241462959507370778130248180630715_wp
    q(8,4) = 0.0_wp
    q(9,4) = 0.0_wp
    q(1,5) = -0.03614399304268576976452921364705641609825_wp
    q(2,5) = 0.1051508663818248421520867474440761344449_wp
    q(3,5) = 0.01609777419666805778308369351834662756172_wp
    q(4,5) = -0.7080721616106272031118456849378369336023_wp
    q(5,5) = 0.0_wp
    q(6,5) = 0.7692160858661111736140494493705980473867_wp
    q(7,5) = -0.1645296432652024882569506157166433921544_wp
    q(8,5) = 0.01828107147391138758410562396851593246160_wp
    q(9,5) = 0.0_wp
    q(1,6) = -0.01141318406360863692889821914555232596651_wp
    q(2,6) = 0.02049729840293952857599941220163960606616_wp
    q(3,6) = 0.01113095018331244864875173213474522093204_wp
    q(4,6) = 0.06324365883611076515355091406993789453750_wp
    q(5,6) = -0.6916640154753724474963890679085181638850_wp
    q(6,6) = 0.0_wp
    q(7,6) = 0.7397091390607520376247117645715851236273_wp
    q(8,6) = -0.1479418278121504075249423529143170247255_wp
    q(9,6) = 0.01643798086801671194721581699047966941394_wp

    first = .false. 
  end if

  idel = 1.0_wp / delta

  if (gsize < 3) 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,:,:) + &
                      q(5,1) * var(5,:,:) + q(6,1) * var(6,:,:) ) * idel
      dvar(2,:,:) = ( q(1,2) * var(1,:,:) + q(3,2) * var(3,:,:) + &
                      q(4,2) * var(4,:,:) + q(5,2) * var(5,:,:) + &
                      q(6,2) * var(6,:,:) ) * idel
      dvar(3,:,:) = ( q(1,3) * var(1,:,:) + q(2,3) * var(2,:,:) + &
                      q(4,3) * var(4,:,:) + q(5,3) * var(5,:,:) + &
                      q(6,3) * var(6,:,:) ) * idel
      dvar(4,:,:) = ( q(1,4) * var(1,:,:) + q(2,4) * var(2,:,:) + &
                      q(3,4) * var(3,:,:) + q(5,4) * var(5,:,:) + &
                      q(6,4) * var(6,:,:) + q(7,4) * var(7,:,:) ) * idel
      dvar(5,:,:) = ( q(1,5) * var(1,:,:) + q(2,5) * var(2,:,:) + &
                      q(3,5) * var(3,:,:) + q(4,5) * var(4,:,:) + &
                      q(6,5) * var(6,:,:) + q(7,5) * var(7,:,:) + &
                      q(8,5) * var(8,:,:) ) * idel
      dvar(6,:,:) = ( q(1,6) * var(1,:,:) + q(2,6) * var(2,:,:) + &
                      q(3,6) * var(3,:,:) + q(4,6) * var(4,:,:) + &
                      q(5,6) * var(5,:,:) + q(7,6) * var(7,:,:) + &
                      q(8,6) * var(8,:,:) + q(9,6) * var(9,:,:) ) * idel

      il = 7
    end if
    if ( bb(2) == 0 ) then
      ir = ni - gsize
    else
      dvar(ni-5,:,:) = - ( q(1,6) * var(ni,:,:) + q(2,6) * var(ni-1,:,:) + &
                           q(3,6) * var(ni-2,:,:) + q(4,6) * var(ni-3,:,:) + &
                           q(5,6) * var(ni-4,:,:) + q(7,6) * var(ni-6,:,:) + &
                           q(8,6) * var(ni-7,:,:) + &
                           q(9,6) * var(ni-8,:,:) ) * idel
      dvar(ni-4,:,:) = - ( q(1,5) * var(ni,:,:) + q(2,5) * var(ni-1,:,:) + &
                           q(3,5) * var(ni-2,:,:) + q(4,5) * var(ni-3,:,:) + &
                           q(6,5) * var(ni-5,:,:) + q(7,5) * var(ni-6,:,:) + &
                           q(8,5) * var(ni-7,:,:) ) * idel
      dvar(ni-3,:,:) = - ( q(1,4) * var(ni,:,:) + q(2,4) * var(ni-1,:,:) + &
                           q(3,4) * var(ni-2,:,:) + q(5,4) * var(ni-4,:,:) + &
                           q(6,4) * var(ni-5,:,:) + &
                           q(7,4) * var(ni-6,:,:) ) * 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,:,:) + &
                           q(6,3) * var(ni-5,:,:) ) * idel
      dvar(ni-1,:,:) = - ( q(1,2) * var(ni,:,:) + q(3,2) * var(ni-2,:,:) + &
                           q(4,2) * var(ni-3,:,:) + q(5,2) * var(ni-4,:,:) + &
                           q(6,2) * var(ni-5,:,:) ) * 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,:,:) + &
                           q(5,1) * var(ni-4,:,:) + &
                           q(6,1) * var(ni-5,:,:) ) * idel

      ir = ni - 6
    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,:,:) ) + &
                        a(3) * ( var(il+3:ir+3,:,:) - &
                                 var(il-3:ir-3,:,:) ) ) * 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,:) + &
                        q(5,1) * var(:,5,:) + q(6,1) * var(:,6,:) ) * idel
        dvar(:,2,:) = ( q(1,2) * var(:,1,:) + q(3,2) * var(:,3,:) + &
                        q(4,2) * var(:,4,:) + q(5,2) * var(:,5,:) + &
                        q(6,2) * var(:,6,:) ) * idel
        dvar(:,3,:) = ( q(1,3) * var(:,1,:) + q(2,3) * var(:,2,:) + &
                        q(4,3) * var(:,4,:) + q(5,3) * var(:,5,:) + &
                        q(6,3) * var(:,6,:) ) * idel
        dvar(:,4,:) = ( q(1,4) * var(:,1,:) + q(2,4) * var(:,2,:) + &
                        q(3,4) * var(:,3,:) + q(5,4) * var(:,5,:) + &
                        q(6,4) * var(:,6,:) + q(7,4) * var(:,7,:) ) * idel
        dvar(:,5,:) = ( q(1,5) * var(:,1,:) + q(2,5) * var(:,2,:) + &
                        q(3,5) * var(:,3,:) + q(4,5) * var(:,4,:) + &
                        q(6,5) * var(:,6,:) + q(7,5) * var(:,7,:) + &
                        q(8,5) * var(:,8,:) ) * idel
        dvar(:,6,:) = ( q(1,6) * var(:,1,:) + q(2,6) * var(:,2,:) + &
                        q(3,6) * var(:,3,:) + q(4,6) * var(:,4,:) + &
                        q(5,6) * var(:,5,:) + q(7,6) * var(:,7,:) + &
                        q(8,6) * var(:,8,:) + q(9,6) * var(:,9,:) ) * idel
  
        jl = 7
      end if
      if ( bb(2) == 0 ) then
        jr = nj - gsize
      else
        dvar(:,nj-5,:) = - ( q(1,6) * var(:,nj,:) + q(2,6) * var(:,nj-1,:) + &
                             q(3,6) * var(:,nj-2,:) + q(4,6) * var(:,nj-3,:) + &
                             q(5,6) * var(:,nj-4,:) + q(7,6) * var(:,nj-6,:) + &
                             q(8,6) * var(:,nj-7,:) + &
                             q(9,6) * var(:,nj-8,:) ) * idel
        dvar(:,nj-4,:) = - ( q(1,5) * var(:,nj,:) + q(2,5) * var(:,nj-1,:) + &
                             q(3,5) * var(:,nj-2,:) + q(4,5) * var(:,nj-3,:) + &
                             q(6,5) * var(:,nj-5,:) + q(7,5) * var(:,nj-6,:) + &
                             q(8,5) * var(:,nj-7,:) ) * idel
        dvar(:,nj-3,:) = - ( q(1,4) * var(:,nj,:) + q(2,4) * var(:,nj-1,:) + &
                             q(3,4) * var(:,nj-2,:) + q(5,4) * var(:,nj-4,:) + &
                             q(6,4) * var(:,nj-5,:) + &
                             q(7,4) * var(:,nj-6,:) ) * 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,:) + &
                             q(6,3) * var(:,nj-5,:) ) * idel
        dvar(:,nj-1,:) = - ( q(1,2) * var(:,nj,:) + q(3,2) * var(:,nj-2,:) + &
                             q(4,2) * var(:,nj-3,:) + q(5,2) * var(:,nj-4,:) + &
                             q(6,2) * var(:,nj-5,:) ) * 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,:) + &
                             q(5,1) * var(:,nj-4,:) + &
                             q(6,1) * var(:,nj-5,:) ) * idel
  
        jr = nj - 6
      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,:) ) + &
                          a(3) * ( var(:,jl+3:jr+3,:) - &
                                   var(:,jl-3:jr-3,:) ) ) * 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) + &
                        q(5,1) * var(:,:,5) + q(6,1) * var(:,:,6) ) * idel
        dvar(:,:,2) = ( q(1,2) * var(:,:,1) + q(3,2) * var(:,:,3) + &
                        q(4,2) * var(:,:,4) + q(5,2) * var(:,:,5) + &
                        q(6,2) * var(:,:,6) ) * idel
        dvar(:,:,3) = ( q(1,3) * var(:,:,1) + q(2,3) * var(:,:,2) + &
                        q(4,3) * var(:,:,4) + q(5,3) * var(:,:,5) + &
                        q(6,3) * var(:,:,6) ) * idel
        dvar(:,:,4) = ( q(1,4) * var(:,:,1) + q(2,4) * var(:,:,2) + &
                        q(3,4) * var(:,:,3) + q(5,4) * var(:,:,5) + &
                        q(6,4) * var(:,:,6) + q(7,4) * var(:,:,7) ) * idel
        dvar(:,:,5) = ( q(1,5) * var(:,:,1) + q(2,5) * var(:,:,2) + &
                        q(3,5) * var(:,:,3) + q(4,5) * var(:,:,4) + &
                        q(6,5) * var(:,:,6) + q(7,5) * var(:,:,7) + &
                        q(8,5) * var(:,:,8) ) * idel
        dvar(:,:,6) = ( q(1,6) * var(:,:,1) + q(2,6) * var(:,:,2) + &
                        q(3,6) * var(:,:,3) + q(4,6) * var(:,:,4) + &
                        q(5,6) * var(:,:,5) + q(7,6) * var(:,:,7) + &
                        q(8,6) * var(:,:,8) + q(9,6) * var(:,:,9) ) * idel
  
        kl = 7
      end if
      if ( bb(2) == 0 ) then
        kr = nk - gsize
      else 
        dvar(:,:,nk-5) = - ( q(1,6) * var(:,:,nk) + q(2,6) * var(:,:,nk-1) + &
                             q(3,6) * var(:,:,nk-2) + q(4,6) * var(:,:,nk-3) + &
                             q(5,6) * var(:,:,nk-4) + q(7,6) * var(:,:,nk-6) + &
                             q(8,6) * var(:,:,nk-7) + &
                             q(9,6) * var(:,:,nk-8) ) * idel
        dvar(:,:,nk-4) = - ( q(1,5) * var(:,:,nk) + q(2,5) * var(:,:,nk-1) + &
                             q(3,5) * var(:,:,nk-2) + q(4,5) * var(:,:,nk-3) + &
                             q(6,5) * var(:,:,nk-5) + q(7,5) * var(:,:,nk-6) + &
                             q(8,5) * var(:,:,nk-7) ) * idel
        dvar(:,:,nk-3) = - ( q(1,4) * var(:,:,nk) + q(2,4) * var(:,:,nk-1) + &
                             q(3,4) * var(:,:,nk-2) + q(5,4) * var(:,:,nk-4) + &
                             q(6,4) * var(:,:,nk-5) + &
                             q(7,4) * var(:,:,nk-6) ) * 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) + &
                             q(6,3) * var(:,:,nk-5) ) * idel
        dvar(:,:,nk-1) = - ( q(1,2) * var(:,:,nk) + q(3,2) * var(:,:,nk-2) + &
                             q(4,2) * var(:,:,nk-3) + q(5,2) * var(:,:,nk-4) + &
                             q(6,2) * var(:,:,nk-5) ) * 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) + &
                             q(5,1) * var(:,:,nk-4) + &
                             q(6,1) * var(:,:,nk-5) ) * idel
  
        kr = nk - 6
      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) ) + &
                          a(3) * ( var(:,:,kl+3:kr+3) - &
                                   var(:,:,kl-3:kr-3) ) ) * idel
    end if
  end select direction
end subroutine deriv_gf_6_3_opt

subroutine up_deriv_gf_6_3_opt ( 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(-3:3), save :: a1, a2
  CCTK_REAL, dimension(9,6), save :: q1, q2 
  CCTK_REAL :: idel

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

  logical, save :: first = .true.

  if ( first ) then
    a1(-3) = -0.03333333333333333333333333333333333333333_wp
    a1(-2) = 0.2500000000000000000000000000000000000000_wp
    a1(-1) = -1.000000000000000000000000000000000000000_wp
    a1(0) = 0.3333333333333333333333333333333333333333_wp
    a1(1) = 0.5000000000000000000000000000000000000000_wp
    a1(2) = -0.05000000000000000000000000000000000000000_wp
    a1(3) = 0.0_wp

    q1(1,1) = -1.529782401641145871492417026888416733827_wp
    q1(2,1) = 1.875173434574214612019635448027523684006_wp
    q1(3,1) = 0.01654806207923748080513710988516743214616_wp
    q1(4,1) = -0.5028607772715413996090537489364770622056_wp
    q1(5,1) = 0.1042956382142412661862395105494407610836_wp
    q1(6,1) = 0.03662604404499391209045870736276191879693_wp
    q1(7,1) = 0.0_wp
    q1(8,1) = 0.0_wp
    q1(9,1) = 0.0_wp
    q1(1,2) = -0.4980311697078740570957463452947252270539_wp
    q1(2,2) = 0.1198701406809289935902771997003246482977_wp
    q1(3,2) = 0.1435237728855054645449658808725554143554_wp
    q1(4,2) = 0.3305195343366502157809026121276381376866_wp
    q1(5,2) = -0.08093510151416251407907993920061497636080_wp
    q1(6,2) = -0.01494717668104810274131940820517799692506_wp
    q1(7,2) = 0.0_wp
    q1(8,2) = 0.0_wp
    q1(9,2) = 0.0_wp
    q1(1,3) = 0.1510193840162481602541749755506215314879_wp
    q1(2,3) = -0.9553948881729209047544587102401343033294_wp
    q1(3,3) = 0.5046108447067502766506824050165990409443_wp
    q1(4,3) = 0.2083429547547929960329552908847933822643_wp
    q1(5,3) = 0.1359641889806059957290021663880456685165_wp
    q1(6,3) = -0.04454248428547652391235612759992531988360_wp
    q1(7,3) = 0.0_wp
    q1(8,3) = 0.0_wp
    q1(9,3) = 0.0_wp
    q1(1,4) = 0.1012044551031772450693035010120167087525_wp
    q1(2,4) = -0.2092303756284921666519857323651163414844_wp
    q1(3,4) = -0.5084563818511371174184254139930350549204_wp
    q1(4,4) = 0.2687068482925919014741556260496361261429_wp
    q1(5,4) = 0.3188546759663099905683904142802147346506_wp
    q1(6,4) = 0.02892077811755014695856160501628382685880_wp
    q1(7,4) = 0.0_wp
    q1(8,4) = 0.0_wp
    q1(9,4) = 0.0_wp
    q1(1,5) = -0.03614399304268576976452921364705641609825_wp
    q1(2,5) = 0.08686979490791345456798112347556020198327_wp
    q1(3,5) = 0.1257842030401363832877174373294422223313_wp
    q1(4,5) = -0.9822882337192980168734300444655759205262_wp
    q1(5,5) = 0.3656214294782277516821124793703186492319_wp
    q1(6,5) = 0.4950000137574403598524650898428590604628_wp
    q1(7,5) = -0.05484321442173416275231687190554779738479_wp
    q1(8,5) = 0.0_wp
    q1(9,5) = 0.0_wp
    q1(1,6) = -0.01141318406360863692889821914555232596651_wp
    q1(2,6) = 0.02049729840293952857599941220163960606616_wp
    q1(3,6) = -0.005307030684704263298464084855734448481896_wp
    q1(4,6) = 0.1618715440442110368368458160128159110211_wp
    q1(5,6) = -0.9382337284956231267046263227657132050941_wp
    q1(6,6) = 0.3287596173603342389443163398095933882788_wp
    q1(7,6) = 0.4931394260405013584164745097143900824182_wp
    q1(8,6) = -0.04931394260405013584164745097143900824182_wp
    q1(9,6) = 0.0_wp

    a2(-3) = 0.0_wp
    a2(-2) = 0.05000000000000000000000000000000000000000_wp
    a2(-1) = -0.5000000000000000000000000000000000000000_wp
    a2(0) = -0.3333333333333333333333333333333333333333_wp
    a2(1) = 1.000000000000000000000000000000000000000_wp
    a2(2) = -0.2500000000000000000000000000000000000000_wp
    a2(3) = 0.03333333333333333333333333333333333333333_wp

    q2(1,1) = -1.635284636237086966078100959777273060297_wp
    q2(2,1) = 2.191680138362037895776687246694092663418_wp
    q2(3,1) = -0.2999586417085858029519146887814015472662_wp
    q2(4,1) = -0.3973585426756003050233698160476207357348_wp
    q2(5,1) = 0.1042956382142412661862395105494407610836_wp
    q2(6,1) = 0.03662604404499391209045870736276191879693_wp
    q2(7,1) = 0.0_wp
    q2(8,1) = 0.0_wp
    q2(9,1) = 0.0_wp
    q2(1,2) = -0.4261090852993166609415800254745304380753_wp
    q2(2,2) = -0.1198701406809289935902771997003246482977_wp
    q2(3,2) = 0.4312121105197350491616311601533345702699_wp
    q2(4,2) = 0.1866753655195354234725699724872485597294_wp
    q2(5,2) = -0.05696107337797671536102449926055004670126_wp
    q2(6,2) = -0.01494717668104810274131940820517799692506_wp
    q2(7,2) = 0.0_wp
    q2(8,2) = 0.0_wp
    q2(9,2) = 0.0_wp
    q2(1,3) = -0.008331409049041400793408941823041323547141_wp
    q2(2,3) = -0.3179917159117626605641230407454828831892_wp
    q2(3,3) = -0.5046108447067502766506824050165990409443_wp
    q2(4,3) = 1.005096920081240801270874877753107657440_wp
    q2(5,3) = -0.1827373971499731263661656683592800415536_wp
    q2(6,3) = 0.008574446736286663103505178191295631794744_wp
    q2(7,3) = 0.0_wp
    q2(8,3) = 0.0_wp
    q2(9,3) = 0.0_wp
    q2(1,4) = 0.1280751399324364352167190636169803213668_wp
    q2(2,4) = -0.3704544846040473075364791079948980171701_wp
    q2(3,4) = -0.1053961094122492652071919749185808657060_wp
    q2(4,4) = -0.2687068482925919014741556260496361261429_wp
    q2(5,4) = 0.7219149484051978427796238533546689238650_wp
    q2(6,4) = -0.1323033308580049939259317706134978488270_wp
    q2(7,4) = 0.02687068482925919014741556260496361261429_wp
    q2(8,4) = 0.0_wp
    q2(9,4) = 0.0_wp
    q2(1,5) = -0.03614399304268576976452921364705641609825_wp
    q2(2,5) = 0.1234319378557362297361923714125920669065_wp
    q2(3,5) = -0.09358865464680026772155005029274896720786_wp
    q2(4,5) = -0.4338560895019563893502613254100979466783_wp
    q2(5,5) = -0.3656214294782277516821124793703186492319_wp
    q2(6,5) = 1.043432157974781987375633808898337034311_wp
    q2(7,5) = -0.2742160721086708137615843595277389869240_wp
    q2(8,5) = 0.03656214294782277516821124793703186492319_wp
    q2(9,5) = 0.0_wp
    q2(1,6) = -0.01141318406360863692889821914555232596651_wp
    q2(2,6) = 0.02049729840293952857599941220163960606616_wp
    q2(3,6) = 0.02756893105132916059596754912522489034598_wp
    q2(4,6) = -0.03538422637198950652974398787294012194614_wp
    q2(5,6) = -0.4450943024551217682881518130513231226759_wp
    q2(6,6) = -0.3287596173603342389443163398095933882788_wp
    q2(7,6) = 0.9862788520810027168329490194287801648364_wp
    q2(8,6) = -0.2465697130202506792082372548571950412091_wp
    q2(9,6) = 0.03287596173603342389443163398095933882788_wp

    first = .false. 
  end if

  idel = 1.0_wp / delta

  if (gsize < 3) 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,:,:) + &
                        q1(5,1) * var(5,:,:) + q1(6,1) * var(6,:,:) ) * idel
      elsewhere
        dvar(1,:,:) = ( q2(1,1) * var(1,:,:) + q2(2,1) * var(2,:,:) + &
                        q2(3,1) * var(3,:,:) + q2(4,1) * var(4,:,:) + &
                        q2(5,1) * var(5,:,:) + q2(6,1) * var(6,:,:) ) * 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,:,:) + &
                        q1(5,2) * var(5,:,:) + q1(6,2) * var(6,:,:) ) * idel
      elsewhere
        dvar(2,:,:) = ( q2(1,2) * var(1,:,:) + q2(2,2) * var(2,:,:) + &
                        q2(3,2) * var(3,:,:) + q2(4,2) * var(4,:,:) + &
                        q2(5,2) * var(5,:,:) + q2(6,2) * var(6,:,:) ) * 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,:,:) + q1(6,3) * var(6,:,:) ) * 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,:,:) + q2(6,3) * var(6,:,:) ) * 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,:,:) + &
                        q1(7,4) * var(7,:,:) ) * 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,:,:) + &
                        q2(7,4) * var(7,:,:) ) * idel
      end where
      where ( up(5,:,:) < zero )
        dvar(5,:,:) = ( q1(1,5) * var(1,:,:) + q1(2,5) * var(2,:,:) + &
                        q1(3,5) * var(3,:,:) + q1(4,5) * var(4,:,:) + &
                        q1(5,5) * var(5,:,:) + q1(6,5) * var(6,:,:) + &
                        q1(7,5) * var(7,:,:) + q1(8,5) * var(8,:,:) ) * idel
      elsewhere
        dvar(5,:,:) = ( q2(1,5) * var(1,:,:) + q2(2,5) * var(2,:,:) + &
                        q2(3,5) * var(3,:,:) + q2(4,5) * var(4,:,:) + &
                        q2(5,5) * var(5,:,:) + q2(6,5) * var(6,:,:) + &
                        q2(7,5) * var(7,:,:) + q2(8,5) * var(8,:,:) ) * idel
      end where
      where ( up(6,:,:) < zero )
        dvar(6,:,:) = ( q1(1,6) * var(1,:,:) + q1(2,6) * var(2,:,:) + &
                        q1(3,6) * var(3,:,:) + q1(4,6) * var(4,:,:) + &
                        q1(5,6) * var(5,:,:) + q1(6,6) * var(6,:,:) + &
                        q1(7,6) * var(7,:,:) + q1(8,6) * var(8,:,:) + &
                        q1(9,6) * var(9,:,:) ) * idel
      elsewhere
        dvar(6,:,:) = ( q2(1,6) * var(1,:,:) + q2(2,6) * var(2,:,:) + &
                        q2(3,6) * var(3,:,:) + q2(4,6) * var(4,:,:) + &
                        q2(5,6) * var(5,:,:) + q2(6,6) * var(6,:,:) + &
                        q2(7,6) * var(7,:,:) + q2(8,6) * var(8,:,:) + &
                        q2(9,6) * var(9,:,:) ) * idel
      end where

      il = 7
    end if
    if ( bb(2) == 0 ) then
      ir = ni - gsize
    else
      where ( up(ni-5,:,:) < zero )
        dvar(ni-5,:,:) = - ( q2(1,6) * var(ni,:,:) + &
                             q2(2,6) * var(ni-1,:,:) + &
                             q2(3,6) * var(ni-2,:,:) + &
                             q2(4,6) * var(ni-3,:,:) + &
                             q2(5,6) * var(ni-4,:,:) + &
                             q2(6,6) * var(ni-5,:,:) + &
                             q2(7,6) * var(ni-6,:,:) + &
                             q2(8,6) * var(ni-7,:,:) + &
                             q2(9,6) * var(ni-8,:,:) ) * idel
      elsewhere
        dvar(ni-5,:,:) = - ( q1(1,6) * var(ni,:,:) + &
                             q1(2,6) * var(ni-1,:,:) + &
                             q1(3,6) * var(ni-2,:,:) + &
                             q1(4,6) * var(ni-3,:,:) + &
                             q1(5,6) * var(ni-4,:,:) + &
                             q1(6,6) * var(ni-5,:,:) + &
                             q1(7,6) * var(ni-6,:,:) + &
                             q1(8,6) * var(ni-7,:,:) + &
                             q1(9,6) * var(ni-8,:,:) ) * idel
      end where
      where ( up(ni-4,:,:) < zero )
        dvar(ni-4,:,:) = - ( q2(1,5) * var(ni,:,:) + &
                             q2(2,5) * var(ni-1,:,:) + &
                             q2(3,5) * var(ni-2,:,:) + &
                             q2(4,5) * var(ni-3,:,:) + &
                             q2(5,5) * var(ni-4,:,:) + &
                             q2(6,5) * var(ni-5,:,:) + &
                             q2(7,5) * var(ni-6,:,:) + &
                             q2(8,5) * var(ni-7,:,:) ) * idel
      elsewhere
        dvar(ni-4,:,:) = - ( q1(1,5) * var(ni,:,:) + &
                             q1(2,5) * var(ni-1,:,:) + &
                             q1(3,5) * var(ni-2,:,:) + &
                             q1(4,5) * var(ni-3,:,:) + &
                             q1(5,5) * var(ni-4,:,:) + &
                             q1(6,5) * var(ni-5,:,:) + &
                             q1(7,5) * var(ni-6,:,:) + &
                             q1(8,5) * var(ni-7,:,:) ) * idel
      end where
      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,:,:) + &
                             q2(7,4) * var(ni-6,:,:) ) * 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,:,:) + &
                             q1(7,4) * var(ni-6,:,:) ) * 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,:,:) + &
                             q2(6,3) * var(ni-5,:,:) ) * 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,:,:) + &
                             q1(6,3) * var(ni-5,:,:) ) * 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,:,:) + &
                             q2(5,2) * var(ni-4,:,:) + &
                             q2(6,2) * var(ni-5,:,:) ) * 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,:,:) + &
                             q1(5,2) * var(ni-4,:,:) + &
                             q1(6,2) * var(ni-5,:,:) ) * 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,:,:) + &
                             q2(5,1) * var(ni-4,:,:) + &
                             q2(6,1) * var(ni-5,:,:) ) * 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,:,:) + &
                             q1(5,1) * var(ni-4,:,:) + &
                             q1(6,1) * var(ni-5,:,:) ) * idel
      end where

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

    end if
  end select direction
end subroutine up_deriv_gf_6_3_opt