aboutsummaryrefslogtreecommitdiff
path: root/src/All_Coeffs_mod.F90
blob: 69e58a9257d00d0f1ca5222569bf91aa30acb7f7 (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
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
! $Header$

#include "cctk.h"

module All_Coeffs_mod

  implicit none

  CCTK_REAL, parameter :: zero = 0.0
  integer, parameter :: wp = kind(zero)

  contains

    subroutine coeffs_1_2_1 ( a, q )

      CCTK_REAL, dimension(1), intent(OUT) :: a
      CCTK_REAL, dimension(3,2), intent(OUT) :: q

      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

    end subroutine coeffs_1_2_1

    subroutine coeffs_up_2_1 ( a1, q1, a2, q2 )

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

      a1(-1) = -1.0_wp
      a1(0) = 1.0_wp
      a1(1) = zero

      q1(1,1) = zero
      q1(2,1) = zero
      q1(3,1) = zero

      q1(1,2) = -1.0_wp
      q1(2,2) = 1.0_wp
      q1(3,2) = zero

      a2(-1) = zero
      a2(0) = -1.0_wp
      a2(1) = 1.0_wp

      q2(1,1) = -2.0_wp
      q2(2,1) = 2.0_wp
      q2(3,1) = zero

      q2(1,2) = zero
      q2(2,2) = -1.0_wp
      q2(3,2) = 1.0_wp

    end subroutine coeffs_up_2_1

    subroutine coeffs_2_2_1 ( a, q )

      CCTK_REAL, dimension(2), intent(OUT) :: a
      CCTK_REAL, dimension(3,2), intent(OUT) :: q

      a(1) = -2.0_wp
      a(2) = 1.0_wp

      q(1,1) = 1.0_wp
      q(2,1) = -2.0_wp
      q(3,1) = 1.0_wp
      q(1,2) = 1.0_wp
      q(2,2) = -2.0_wp
      q(3,2) = 1.0_wp

    end subroutine coeffs_2_2_1

    subroutine coeffs_1_4_2 ( a, q )

      CCTK_REAL, dimension(2), intent(OUT) :: a
      CCTK_REAL, dimension(6,4), intent(OUT) :: q

      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

    end subroutine coeffs_1_4_2

    subroutine coeffs_up_4_2 ( a1, q1, a2, q2 )

      CCTK_REAL, dimension(-2:2), intent(OUT) :: a1, a2
      CCTK_REAL, dimension(6,4), intent(OUT) :: q1, q2

      a1(-2) = 0.1666666666666666666666666666666666666667_wp
      a1(-1) = -1.0_wp
      a1(0) = 0.5_wp
      a1(1) = 0.3333333333333333333333333333333333333333_wp
      a1(2) = zero

      q1(1,1) = -1.176470588235294117647058823529411764706_wp
      q1(2,1) = 1.264705882352941176470588235294117647059_wp
      q1(3,1) = zero
      q1(4,1) = -0.08823529411764705882352941176470588235294_wp
      q1(5,1) = zero
      q1(6,1) = zero

      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) = zero
      q1(6,2) = zero

      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) = zero
      q1(6,3) = zero

      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) = zero

      a2(-2) = zero
      a2(-1) = -0.3333333333333333333333333333333333333333_wp
      a2(0) = -0.5_wp
      a2(1) = 1.0_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) = zero
      q2(6,1) = zero

      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) = zero
      q2(6,2) = zero

      q2(1,3) = zero
      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) = zero

      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

    end subroutine coeffs_up_4_2

    subroutine coeffs_2_4_2 ( a, q )

      CCTK_REAL, dimension(3), intent(OUT) :: a
      CCTK_REAL, dimension(6,4), intent(OUT) :: q

      a(1) = -2.500000000000000000000000000000000000000_wp
      a(2) = 1.333333333333333333333333333333333333333_wp
      a(3) = -0.08333333333333333333333333333333333333333_wp

      q(1,1) = 2.000000000000000000000000000000000000000_wp
      q(2,1) = -5.000000000000000000000000000000000000000_wp
      q(3,1) = 4.000000000000000000000000000000000000000_wp
      q(4,1) = -1.000000000000000000000000000000000000000_wp
      q(5,1) = 0_wp
      q(6,1) = 0_wp

      q(1,2) = 1.000000000000000000000000000000000000000_wp
      q(2,2) = -2.000000000000000000000000000000000000000_wp
      q(3,2) = 1.000000000000000000000000000000000000000_wp
      q(4,2) = 0_wp
      q(5,2) = 0_wp
      q(6,2) = 0_wp

      q(1,3) = -0.09302325581395348837209302325581395348837_wp
      q(2,3) = 1.372093023255813953488372093023255813953_wp
      q(3,3) = -2.558139534883720930232558139534883720930_wp
      q(4,3) = 1.372093023255813953488372093023255813953_wp
      q(5,3) = -0.09302325581395348837209302325581395348837_wp
      q(6,3) = 0_wp

      q(1,4) = -0.02040816326530612244897959183673469387755_wp
      q(2,4) = 0_wp
      q(3,4) = 1.204081632653061224489795918367346938776_wp
      q(4,4) = -2.408163265306122448979591836734693877551_wp
      q(5,4) = 1.306122448979591836734693877551020408163_wp
      q(6,4) = -0.08163265306122448979591836734693877551020_wp

    end subroutine coeffs_2_4_2

    subroutine coeffs_2_4_2_opt ( a, q )

      CCTK_REAL, dimension(3), intent(OUT) :: a
      CCTK_REAL, dimension(6,4), intent(OUT) :: q

      a(1) = -2.5_wp
      a(2) = 1.333333333333333333333333333333333333333_wp
      a(3) = -0.08333333333333333333333333333333333333333_wp

      q(1,1) = 2.282982974062943480262699765139082828341_wp
      q(2,1) = -6.061186152736038050985124119271560606278_wp
      q(3,1) = 5.414914870314717401313498825695414141704_wp
      q(4,1) = -1.707457435157358700656749412847707070852_wp
      q(5,1) = zero
      q(6,1) = 0.07074574351573587006567494128477070708519_wp

      q(1,2) = 0.6942344983641924259873371181759910117504_wp
      q(2,2) = -0.7761428473759173089500505257318397680987_wp
      q(3,2) = -0.8377735941382550240731690788325510451097_wp
      q(4,2) = 1.227832883028344666046439209128781626417_wp
      q(5,2) = -0.3089460859592171540098546697125061038620_wp
      q(6,2) = 0.0007951460808523949992979469721242789031010_wp

      q(1,3) = 0.4663616929151208330774297682981869862550_wp
      q(2,3) = -1.149503303585047591170162224444663061895_wp
      q(3,3) = 1.934396285188982033906351214796782385029_wp
      q(4,3) = -2.569785963207868885472377980704238646268_wp
      q(5,3) = 1.602587820613377868519202373305847453754_wp
      q(6,3) = -0.2840565319245642588604431512519151168747_wp

      q(1,4) = -0.2658525795443897532890763269063473511119_wp
      q(2,4) = 1.478411022421884393811018639563226856298_wp
      q(3,4) = -2.255118294243640042353311289189433914072_wp
      q(4,4) = 1.553414543643511297084585299252414115549_wp
      q(5,4) = -0.9258553965216912759079296546576971585129_wp
      q(6,4) = 0.4150007042443253806547133319378374518503_wp

    end subroutine coeffs_2_4_2_opt

    subroutine coeffs_1_6_3 ( a, q )

      CCTK_REAL, dimension(3), intent(OUT) :: a
      CCTK_REAL, dimension(9,6), intent(OUT) :: q

      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) = -21600.0_wp/13649.0_wp
      q(2,1) = 81763.0_wp/40947.0_wp
      q(3,1) = 131.0_wp/27298.0_wp
      q(4,1) = -9143.0_wp/13649.0_wp
      q(5,1) = 20539.0_wp/81894.0_wp
      q(6,1) = zero
      q(7,1) = zero
      q(8,1) = zero
      q(9,1) = zero

      q(1,2) = -81763.0_wp/180195.0_wp
      q(2,2) = zero
      q(3,2) = 7357.0_wp/36039.0_wp
      q(4,2) = 30637.0_wp/72078.0_wp
      q(5,2) = -2328.0_wp/12013.0_wp
      q(6,2) = 6611.0_wp/360390.0_wp
      q(7,2) = zero
      q(8,2) = zero
      q(9,2) = zero

      q(1,3) = -131.0_wp/54220.0_wp
      q(2,3) = -7357.0_wp/16266.0_wp
      q(3,3) = zero
      q(4,3) = 645.0_wp/2711.0_wp
      q(5,3) = 11237.0_wp/32532.0_wp
      q(6,3) = -3487.0_wp/27110.0_wp
      q(7,3) = zero
      q(8,3) = zero
      q(9,3) = zero

      q(1,4) = 9143.0_wp/53590.0_wp
      q(2,4) = -30637.0_wp/64308.0_wp
      q(3,4) =  -645.0_wp/5359.0_wp
      q(4,4) = zero
      q(5,4) = 13733.0_wp/32154.0_wp
      q(6,4) = -67.0_wp/4660.0_wp
      q(7,4) = 72.0_wp/5359.0_wp
      q(8,4) = zero
      q(9,4) = zero

      q(1,5) = -20539.0_wp/236310.0_wp
      q(2,5) = 2328.0_wp/7877.0_wp
      q(3,5) = -11237.0_wp/47262.0_wp
      q(4,5) = -13733.0_wp/23631.0_wp
      q(5,5) = zero
      q(6,5) = 89387.0_wp/118155.0_wp
      q(7,5) = -1296.0_wp/7877.0_wp
      q(8,5) = 144.0_wp/7877.0_wp
      q(9,5) = zero

      q(1,6) = zero
      q(2,6) = -6611.0_wp/262806.0_wp
      q(3,6) = 3487.0_wp/43801.0_wp
      q(4,6) = 1541.0_wp/87602.0_wp
      q(5,6) = -89387.0_wp/131403.0_wp
      q(6,6) = zero
      q(7,6) = 32400.0_wp/43801.0_wp
      q(8,6) = -6480.0_wp/43801.0_wp
      q(9,6) = 720.0_wp/43801.0_wp

    end subroutine coeffs_1_6_3 

    subroutine coeffs_1_6_3_opt ( a, q )

      CCTK_REAL, dimension(3), intent(OUT) :: a
      CCTK_REAL, dimension(9,6), intent(OUT) :: q

      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) = zero
      q(8,1) = zero
      q(9,1) = zero

      q(1,2) = -0.4620701275035953590186631853846278325646_wp
      q(2,2) = zero
      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) = zero
      q(8,2) = zero
      q(9,2) = zero

      q(1,3) = 0.07134398748360337973038301686379010397038_wp
      q(2,3) = -0.6366933020423417826592908754928085932593_wp
      q(3,3) = zero
      q(4,3) = 0.6067199374180168986519150843189505198519_wp
      q(5,3) = -0.02338660408468356531858175098561718651857_wp
      q(6,3) = -0.01798401877459493040442547470431484404443_wp
      q(7,3) = zero
      q(8,3) = zero
      q(9,3) = zero

      q(1,4) = 0.1146397975178068401430112823144985150596_wp
      q(2,4) = -0.2898424301162697370942324201800071793273_wp
      q(3,4) = -0.3069262456316931913128086944558079603132_wp
      q(4,4) = zero
      q(5,4) = 0.5203848121857539166740071338174418292578_wp
      q(6,4) = -0.05169127637022742348368508279860701098408_wp
      q(7,4) = 0.01343534241462959507370778130248180630715_wp
      q(8,4) = zero
      q(9,4) = zero

      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) = zero
      q(6,5) = 0.7692160858661111736140494493705980473867_wp
      q(7,5) = -0.1645296432652024882569506157166433921544_wp
      q(8,5) = 0.01828107147391138758410562396851593246160_wp
      q(9,5) = zero

      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) = zero
      q(7,6) = 0.7397091390607520376247117645715851236273_wp
      q(8,6) = -0.1479418278121504075249423529143170247255_wp
      q(9,6) = 0.01643798086801671194721581699047966941394_wp

    end subroutine coeffs_1_6_3_opt

    subroutine coeffs_up_6_3_opt ( a1, q1, a2, q2 )

      CCTK_REAL, dimension(-3:3), intent(OUT) :: a1, a2
      CCTK_REAL, dimension(9,6), intent(OUT) :: q1, q2

      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) = zero

      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) = zero
      q1(8,1) = zero
      q1(9,1) = zero

      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) = zero
      q1(8,2) = zero
      q1(9,2) = zero

      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) = zero
      q1(8,3) = zero
      q1(9,3) = zero

      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) = zero
      q1(8,4) = zero
      q1(9,4) = zero

      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) = zero
      q1(9,5) = zero

      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) = zero

      a2(-3) = zero
      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) = zero
      q2(8,1) = zero
      q2(9,1) = zero

      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) = zero
      q2(8,2) = zero
      q2(9,2) = zero

      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) = zero
      q2(8,3) = zero
      q2(9,3) = zero

      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) = zero
      q2(9,4) = zero

      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) = zero

      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

    end subroutine coeffs_up_6_3_opt

    subroutine coeffs_2_6_3 ( a, q )

      CCTK_REAL, dimension(4), intent(OUT) :: a
      CCTK_REAL, dimension(9,6), intent(OUT) :: q

      a(1) = -2.722222222222222222222222222222222222222_wp
      a(2) = 1.500000000000000000000000000000000000000_wp
      a(3) = -0.1500000000000000000000000000000000000000_wp
      a(4) = 0.01111111111111111111111111111111111111111_wp

      q(1,1) = 2.788238454587637677973966346740908979901_wp
      q(2,1) = -8.024525606271521723203165067037878232838_wp
      q(3,1) = 8.215717879209710113072996800742423132342_wp
      q(4,1) = -3.382384545876376779739663467409089799008_wp
      q(5,1) = 0.2745256062715217232031650670378782328376_wp
      q(6,1) = 0.1284282120790289886927003199257576867658_wp
      q(7,1) = 0_wp
      q(8,1) = 0_wp
      q(9,1) = 0_wp

      q(1,2) = 1.053412969283276450511945392491467576792_wp
      q(2,2) = -2.350398179749715585893060295790671217292_wp
      q(3,2) = 1.867463026166097838452787258248009101251_wp
      q(4,2) = -1.034129692832764505119453924914675767918_wp
      q(5,2) = 0.6003981797497155858930602957906712172924_wp
      q(6,2) = -0.1367463026166097838452787258248009101251_wp
      q(7,2) = 0_wp
      q(8,2) = 0_wp
      q(9,2) = 0_wp

      q(1,3) = -0.6441780400836099840157383499323742776343_wp
      q(2,3) = 4.137556867084716586745358416328538054838_wp
      q(3,3) = -8.108447067502766506824050165990409443010_wp
      q(4,3) = 6.941780400836099840157383499323742776343_wp
      q(5,3) = -2.887556867084716586745358416328538054838_wp
      q(6,3) = 0.5608447067502766506824050165990409443010_wp
      q(7,3) = 0_wp
      q(8,3) = 0_wp
      q(9,3) = 0_wp

      q(1,4) = 0.2133575915904708589911053057162405921503_wp
      q(2,4) = -1.159078186228774025004665049449524164956_wp
      q(3,4) = 3.511693723953473906823412328170678609193_wp
      q(4,4) = -4.723144865335572557069104932512284630217_wp
      q(5,4) = 2.489690240716551595446911737264415002799_wp
      q(6,4) = -0.3414753996392361759034645767245132798408_wp
      q(7,4) = 0.008956894943086396715805187534987870871431_wp
      q(8,4) = 0_wp
      q(9,4) = 0_wp

      q(1,5) = -0.1790783293131903008759680081249206550717_wp
      q(2,5) = 0.9156510515847827006897719097795268926410_wp
      q(3,5) = -1.987601032541999915365409843002835258770_wp
      q(4,5) = 3.387647581566586263806017519360162498413_wp
      q(5,5) = -3.719859506580339384706529558630612331260_wp
      q(6,5) = 1.735582497566755532986331513689644957894_wp
      q(7,5) = -0.1645296432652024882569506157166433921544_wp
      q(8,5) = 0.01218738098260759172273708264567728830773_wp
      q(9,5) = 0_wp

      q(1,6) = 0.04002001476374207590389869333272451922711_wp
      q(2,6) = -0.1875223548929628699496967344733377472356_wp
      q(3,6) = 0.3471267779274445788908929019885390744504_wp
      q(4,6) = -0.4177910702190969764769449708149737829426_wp
      q(5,6) = 1.560601736642238000654475164189554271972_wp
      q(6,6) = -2.684870208442729618045250108445012670944_wp
      q(7,6) = 1.479418278121504075249423529143170247255_wp
      q(8,6) = -0.1479418278121504075249423529143170247255_wp
      q(9,6) = 0.01095865391201114129814387799365311294263_wp

    end subroutine coeffs_2_6_3

    subroutine coeffs_1_8_4 ( a, q )

      CCTK_REAL, dimension(4), intent(OUT) :: a
      CCTK_REAL, dimension(12,8), intent(OUT) :: q
      CCTK_REAL, parameter :: x1 = 0.541_wp, x2 = -0.0675_wp, x3 = 0.748_wp

      a(1) = 4.0_wp/5.0_wp
      a(2) = -1.0_wp/5.0_wp
      a(3) = 4.0_wp/105.0_wp
      a(4) = -1.0_wp/280.0_wp

      q(1,1) = -2540160.0_wp/1498139.0_wp
      q(2,1) = 9.0_wp * ( 2257920.0_wp*x1 + 11289600.0_wp*x2 + &
                          22579200.0_wp*x3 - 15849163.0_wp ) / 5992556.0_wp
      q(3,1) = 3.0_wp * ( -33868800.0_wp*x1 - 162570240.0_wp*x2 - &
                           304819200.0_wp*x3 + 235236677.0_wp ) / 5992556.0_wp
      q(4,1) = ( 609638400.0_wp*x1 + 2743372800.0_wp*x2 + &
                 4572288000.0_wp * x3 - 3577778591.0_wp ) / 17977668.0_wp
      q(5,1) = 3.0_wp * ( -16934400*x1 - 67737600.0_wp*x2 - &
                           84672000.0_wp*x3 + 67906303.0_wp ) / 1498139.0_wp
      q(6,1) = 105.0_wp * ( 967680.0_wp*x1 + 2903040.0_wp*x2 - &
                            305821.0_wp ) / 5992556.0_wp
      q(7,1) = 49.0_wp * ( -1244160.0_wp*x1 + 18662400.0_wp * x3 - &
                            13322233.0_wp ) / 17977668.0_wp
      q(8,1) = 3.0_wp * ( -6773760.0_wp*x2 - 33868800.0_wp*x3 + &
                           24839327.0_wp ) / 5992556.0_wp
      q(9,1) = zero
      q(10,1) = zero
      q(11,1) = zero
      q(12,1) = zero

      q(1,2) = 9.0_wp * ( -2257920.0_wp*x1 - 11289600.0_wp*x2 - &
                           22579200.0_wp*x3 + 15849163.0_wp ) / 31004596.0_wp
      q(2,2) = zero
      q(3,2) = 3.0_wp * ( 7257600.0_wp*x1 + 33868800.0_wp*x2 + &
                          60963840.0_wp*x3 - 47167457.0_wp ) / 2214614.0_wp
      q(4,2) = 3.0_wp * ( -9676800.0_wp*x1 - 42336000.0_wp*x2 - &
                           67737600.0_wp*x3 + 53224573.0_wp ) / 1107307.0_wp
      q(5,2) = 7.0_wp * ( 55987200.0_wp*x1 + 217728000.0_wp*x2 + &
                          261273600.0_wp*x3 - 211102099.0_wp ) / 13287684.0_wp
      q(6,2) = 3.0_wp * ( -11612160.0_wp*x1 - 33868800.0_wp*x2 + &
                           3884117.0_wp ) / 2214614.0_wp
      q(7,2) = 150.0_wp * ( 24192.0_wp*x1 - 338688.0_wp*x3 + &
                            240463.0_wp ) / 1107307.0_wp
      q(8,2) = ( 152409600.0_wp*x2 + 731566080.0_wp*x3 - &
                 536324953.0_wp ) / 46506894.0_wp
      q(9,2) = zero
      q(10,2) = zero
      q(11,2) = zero
      q(12,2) = zero

      q(1,3) = ( 33868800.0_wp*x1 + 162570240.0_wp*x2 + &
               304819200.0_wp*x3 - 235236677.0_wp ) / 1743924.0_wp
      q(2,3) = ( -7257600.0_wp*x1 - 33868800.0_wp*x2 - &
                  60963840.0_wp*x3 + 47167457.0_wp ) / 124566.0_wp
      q(3,3) = zero
      q(4,3) = ( 24192000.0_wp*x1 + 101606400.0_wp*x2 + &
                 152409600.0_wp*x3 - 120219461.0_wp ) / 124566.0_wp
      q(5,3) = ( -72576000.0_wp*x1 - 270950400.0_wp*x2 - &
                  304819200.0_wp*x3 + 249289259.0_wp ) / 249132.0_wp
      q(6,3) = 9.0_wp * ( 806400.0_wp*x1 + 2257920.0_wp*x2 - &
                          290167.0_wp ) / 41522.0_wp
      q(7,3) = 6.0_wp * ( -134400.0_wp*x1 + 1693440.0_wp*x3 - &
                           1191611.0_wp ) / 20761.0_wp
      q(8,3) = 5.0_wp * ( -2257920.0_wp*x2 - 10160640.0_wp*x3 + &
                           7439833.0_wp ) / 290654.0_wp
      q(9,3) = zero
      q(10,3) = zero
      q(11,3) = zero
      q(12,3) = zero

      q(1,4) = ( -609638400.0_wp*x1 - 2743372800.0_wp*x2 - &
                 4572288000.0_wp*x3 + 3577778591.0_wp ) / 109619916.0_wp
      q(2,4) = 3.0_wp * ( 9676800.0_wp*x1 + 42336000.0_wp*x2 + &
                          67737600.0_wp*x3 - 53224573.0_wp ) / 1304999.0_wp
      q(3,4) = 3.0_wp * ( -24192000.0_wp*x1 - 101606400.0_wp*x2 - &
                           152409600.0_wp*x3 + 120219461.0_wp ) / 2609998.0_wp
      q(4,4) = zero
      q(5,4) = 9.0_wp * ( 16128000.0_wp*x1 + 56448000.0_wp*x2 + &
                          56448000.0_wp*x3 - 47206049.0_wp ) / 5219996.0_wp
      q(6,4) = 3.0_wp * ( -19353600.0_wp*x1 - 50803200.0_wp*x2 + &
                           7628371.0_wp ) / 2609998.0_wp
      q(7,4) = 2.0_wp * ( 10886400.0_wp*x1 - 114307200.0_wp*x3 + &
                          79048289.0_wp ) / 3914997.0_wp
      q(8,4) = 75.0_wp * ( 1354752.0_wp*x2 + 5419008.0_wp*x3 - &
                           3952831.0_wp ) / 18269986.0_wp
      q(9,4) = zero
      q(10,4) = zero
      q(11,4) = zero
      q(12,4) = zero

      q(1,5) = 3.0_wp * ( 16934400.0_wp*x1 + 67737600.0_wp*x2 + &
                          84672000.0_wp*x3 - 67906303.0_wp ) / 2096689.0_wp
      q(2,5) = 7.0_wp * ( -55987200.0_wp*x1 - 217728000.0_wp*x2 - &
                           261273600.0_wp*x3 + 211102099.0_wp ) / 3594324.0_wp
      q(3,5) = 3.0_wp * ( 72576000.0_wp*x1 + 270950400.0_wp*x2 + &
                          304819200.0_wp*x3 - 249289259.0_wp ) / 1198108.0_wp
      q(4,5) = 9.0_wp * ( -16128000.0_wp*x1 - 56448000.0_wp*x2 - &
                           56448000.0_wp*x3 + 47206049.0_wp ) / 1198108.0_wp
      q(5,5) = zero
      q(6,5) = 105.0_wp * ( 414720.0_wp*x1 + 967680.0_wp*x2 - &
                            165527.0_wp ) / 1198108.0_wp
      q(7,5) = 15.0_wp * ( -967680.0_wp*x1 + 6773760.0_wp*x3 - &
                            4472029.0_wp ) / 1198108.0_wp
      q(8,5) = ( -304819200.0_wp*x2 - 914457600.0_wp*x3 + &
                 657798011.0_wp ) / 25160268.0_wp
      q(9,5) = -2592.0_wp/299527.0_wp;
      q(10,5) = zero
      q(11,5) = zero
      q(12,5) = zero

      q(1,6) = 5.0_wp * ( -967680.0_wp*x1 - 2903040.0_wp*x2 + &
                           305821.0_wp ) / 1237164.0_wp
      q(2,6) = ( 11612160.0_wp*x1 + 33868800.0_wp*x2 - &
                 3884117.0_wp ) / 618582.0_wp
      q(3,6) = 9.0_wp * ( -806400.0_wp*x1 - 2257920.0_wp*x2 + &
                           290167.0_wp ) / 206194.0_wp
      q(4,6) = ( 19353600.0_wp*x1 + 50803200.0_wp*x2 - &
                 7628371.0_wp ) / 618582.0_wp
      q(5,6) = 35.0_wp * ( -414720.0_wp*x1 - 967680.0_wp*x2 + &
                            165527.0_wp ) / 1237164.0_wp
      q(6,6) = zero
      q(7,6) = 80640.0_wp*x1 / 103097.0_wp
      q(8,6) = 80640.0_wp*x2 / 103097.0_wp
      q(9,6) = 3072.0_wp / 103097.0_wp
      q(10,6) = -288.0_wp/103097.0_wp
      q(11,6) = zero
      q(12,6) = zero

      q(1,7) = 7.0_wp * ( 1244160.0_wp*x1 - 18662400.0_wp*x3 + &
                          13322233.0_wp ) / 8041092.0_wp
      q(2,7) = 150.0_wp * ( -24192.0_wp*x1 + 338688.0_wp*x3 - &
                             240463.0_wp ) / 670091.0_wp
      q(3,7) = 54.0_wp * ( 134400.0_wp*x1 - 1693440.0_wp*x3 + &
                           1191611.0_wp ) / 670091.0_wp
      q(4,7) = 2.0_wp * ( -10886400.0_wp*x1 + 114307200.0_wp*x3 - &
                           79048289.0_wp ) / 2010273.0_wp
      q(5,7) = 15.0_wp * ( 967680.0_wp*x1 - 6773760.0_wp*x3 + &
                           4472029.0_wp ) / 2680364.0_wp
      q(6,7) = -725760.0_wp*x1 / 670091.0_wp
      q(7,7) = zero
      q(8,7) = 725760.0_wp*x3 / 670091.0_wp
      q(9,7) = -145152.0_wp/670091.0_wp
      q(10,7) = 27648.0_wp/670091.0_wp
      q(11,7) = -2592.0_wp/670091.0_wp
      q(12,7) = zero

      q(1,8) = 3.0_wp * ( 6773760.0_wp*x2 + 33868800.0_wp*x3 - &
                          24839327.0_wp ) / 20510956.0_wp
      q(2,8) = ( -152409600.0_wp*x2 - 731566080.0_wp*x3 + &
                  536324953.0_wp ) / 30766434.0_wp
      q(3,8) = 45.0_wp * ( 2257920.0_wp*x2 + 10160640.0_wp*x3 - &
                           7439833.0_wp ) / 10255478.0_wp
      q(4,8) = 75.0_wp * ( -1354752.0_wp*x2 - 5419008.0_wp*x3 + &
                            3952831.0_wp ) / 10255478.0_wp
      q(5,8) = ( 304819200.0_wp*x2 + 914457600.0_wp*x3 - &
                 657798011.0_wp ) / 61532868.0_wp
      q(6,8) = -5080320.0_wp*x2 / 5127739.0_wp
      q(7,8) = -5080320.0_wp*x3 / 5127739.0_wp
      q(8,8) = zero
      q(9,8) = 4064256.0_wp/5127739.0_wp
      q(10,8) = -1016064.0_wp/5127739.0_wp
      q(11,8) = 193536.0_wp/5127739.0_wp
      q(12,8) = -18144.0_wp/5127739.0_wp

    end subroutine coeffs_1_8_4

    subroutine coeffs_1_8_4_opt ( a, q )
  
      CCTK_REAL, dimension(4), intent(OUT) :: a
      CCTK_REAL, dimension(12,8), intent(OUT) :: q

      a(1) = 4.0_wp/5.0_wp
      a(2) = -1.0_wp/5.0_wp
      a(3) = 4.0_wp/105.0_wp
      a(4) = -1.0_wp/280.0_wp

      q(1,1) = -1.695543604431898508749855654248370812054_wp
      q(2,1) = 2.244525109969933844775276252674718713114_wp
      q(3,1) = -0.002163206607021876385102025416828692597556_wp
      q(4,1) = -0.8963677947307838167920190976769921451579_wp
      q(5,1) = 0.2272938673777916084251861567415263328923_wp
      q(6,1) = 0.1564962758424078514454013191006671546757_wp
      q(7,1) = 0.003067629343740613524557473056482141287158_wp
      q(8,1) = -0.03730827676416971624344442423120269215938_wp
      q(9,1) = zero
      q(10,1) = zero
      q(11,1) = zero
      q(12,1) = zero
  
      q(1,2) = -0.4338209217401506177055540526837828066711_wp
      q(2,2) = zero
      q(3,2) = 0.06631961192627293470305188523936353462439_wp
      q(4,2) = 0.6021342012739069461791324177355805603656_wp
      q(5,2) = -0.2213339967864188453859900127368926875520_wp
      q(6,2) = -0.001759980003751516130708214237001705688807_wp
      q(7,2) = -0.03586646217917856587352426166797879385310_wp
      q(8,2) = 0.02432754750931966421359223835071189877503_wp
      q(9,2) = zero
      q(10,2) = zero
      q(11,2) = zero
      q(12,2) = zero
  
      q(1,3) = 0.002477771724790106958560398469919805908522_wp
      q(2,3) = -0.3930241559935857537756812928554282732813_wp
      q(3,3) = zero
      q(4,3) = -0.6080263371305341035810612202824357742461_wp
      q(5,3) = 2.372533124692864412954923458457782022478_wp
      q(6,3) = -2.126112217612194556872952028715623246817_wp
      q(7,3) = 0.9075309435559160631495435310307745215131_wp
      q(8,3) = -0.1553791292372561688333328461049890555549_wp
      q(9,3) = zero
      q(10,3) = zero
      q(11,3) = zero
      q(12,3) = zero
  
      q(1,4) = 0.1470043328582935680963279007411074464266_wp
      q(2,4) = -0.5109179516689331400658365102850902594996_wp
      q(3,4) = 0.08705685833207777685654755900085198684605_wp
      q(4,4) = zero
      q(5,4) = -0.1306945333157231839535552435788086367322_wp
      q(6,4) = 0.5981933016362239542598518794410035592720_wp
      q(7,4) = -0.2031099149801233772210520666270517105660_wp
      q(8,4) = 0.01246790713818440202771648130798761425319_wp
      q(9,4) = zero
      q(10,4) = zero
      q(11,4) = zero
      q(12,4) = zero
  
      q(1,5) = -0.1624073990846984662267508265053107632238_wp
      q(2,5) = 0.8182390368133059538132603839842499379861_wp
      q(3,5) = -1.480018301574606037840376638130713134769_wp
      q(4,5) = 0.5694185675497882973361521309100738568707_wp
      q(5,5) = zero
      q(6,5) = 0.02929508293831017534661096502884570435894_wp
      q(7,5) = 0.2825910005984090924748504462528740123097_wp
      q(8,5) = -0.04846434332860791750685068920563907354793_wp
      q(9,5) = -0.008653643911901097396895772334380539984709_wp
      q(10,5) = zero
      q(11,5) = zero
      q(12,5) = zero
  
      q(1,6) = -0.03609686950604370828405582087760384548297_wp
      q(2,6) = 0.002100328577309696555612805397001149641337_wp
      q(3,6) = 0.4281425817419204360479874008765051769418_wp
      q(4,6) = -0.8413238238875046736839550701643761633826_wp
      q(5,6) = -0.009456755727630000971085851751473603919416_wp
      q(6,6) = zero
      q(7,6) = 0.5172389789041733771234394082779034212042_wp
      q(8,6) = -0.08760813565107717873901807278223424202774_wp
      q(9,6) = 0.02979718129528502284256573906127239395909_wp
      q(10,6) = -0.002793485746432970891490538036994286933664_wp
      q(11,6) = zero
      q(12,6) = zero
  
      q(1,7) = -0.0009797678134978722516935350416937185004514_wp
      q(2,7) = 0.05926834509975463070197111976550139351983_wp
      q(3,7) = -0.2530570463899371286637621744353665227114_wp
      q(4,7) = 0.3955555826583941989223787901885502942688_wp
      q(5,7) = -0.1263166266018192756531792392597193430222_wp
      q(6,7) = -0.7162192643577544899896537844517283192029_wp
      q(7,7) = zero
      q(8,7) = 0.8209721963136350137518635528607308559584_wp
      q(9,7) = -0.2166153552278720352907291696202456084323_wp
      q(10,7) = 0.04126006766245181624585317516576106827282_wp
      q(11,7) = -0.003868131343354857773048735171790100150577_wp
      q(12,7) = zero
  
      q(1,8) = 0.01090012273307913185972171872891927027272_wp
      q(2,8) = -0.03677379943661633440187210478144113487019_wp
      q(3,8) = 0.03963287609450569641558898819403594845401_wp
      q(4,8) = -0.02221139656912423686782339404035534965598_wp
      q(5,8) = 0.01981665906734246925389947980969432794418_wp
      q(6,8) = 0.1109698768905366566973791789571759073485_wp
      q(7,8) = -0.7509903604688148129224205834189298636300_wp
      q(8,8) = zero
      q(9,8) = 0.7926019635554773751160111698352821779736_wp
      q(10,8) = -0.1981504908888693437790027924588205444934_wp
      q(11,8) = 0.03774295064549892262457196046834677037969_wp
      q(12,8) = -0.003538401623015523996053621293907509723096_wp

    end subroutine coeffs_1_8_4_opt

    subroutine coeffs_up_8_4_opt ( a1, q1, a2, q2 )

      CCTK_REAL, dimension(-4:4), intent(OUT) :: a1, a2
      CCTK_REAL, dimension(12,8), intent(OUT) :: q1, q2

      a1(-4) = 0.007142857142857142857142857142857142857143_wp
      a1(-3) = -0.06666666666666666666666666666666666666667_wp
      a1(-2) = 0.3_wp
      a1(-1) = -1.0_wp
      a1(0) = 0.25_wp
      a1(1) = 0.6_wp
      a1(2) = -0.1_wp
      a1(3) = 0.009523809523809523809523809523809523809524_wp
      a1(4) = zero

      q1(1,1) = -1.683432578685956376544499542432311020539_wp
      q1(2,1) = 2.206194911574014373726323572790957192479_wp
      q1(3,1) = 0.01048240012328238642042347285056590509052_wp
      q1(4,1) = -0.7964177278055455600674641925035752355813_wp
      q1(5,1) = 0.04375559918421074726330123748705042737925_wp
      q1(6,1) = 0.3015760468426855559502835132688146245819_wp
      q1(7,1) = -0.05430139927462465206917597295687512552949_wp
      q1(8,1) = -0.02785725195806647467919208850462676788091_wp
      q1(9,1) = zero
      q1(10,1) = zero
      q1(11,1) = zero
      q1(12,1) = zero

      q1(1,2) = -0.4451389901846271204230470438663315357999_wp
      q1(2,2) = 0.03979384217746298000464189244717138065595_wp
      q1(3,2) = 0.04117629055092890361160914968898310705568_wp
      q1(4,2) = 0.5200531283081606472663885367960872642342_wp
      q1(5,2) = -0.04132093916257534380117029102933732094903_wp
      q1(6,2) = -0.1571076966553020594965523100929511511562_wp
      q1(7,2) = 0.02929412303865500620649355002594213960776_wp
      q1(8,2) = 0.01325024192729698663163651603043611635147_wp
      q1(9,2) = zero
      q1(10,2) = zero
      q1(11,2) = zero
      q1(12,2) = zero

      q1(1,3) = 0.1544592597396873173299548977341719367844_wp
      q1(2,3) = -1.020860666961436408391054234781341175626_wp
      q1(3,3) = 0.7352247001589518809305910119936419247628_wp
      q1(4,3) = -0.1443377614671085870438942853092495421459_wp
      q1(5,3) = 0.4029099592852542668089969628980513126845_wp
      q1(6,3) = -0.1256847607545176255327746274725689708913_wp
      q1(7,3) = -0.01218362768224537680217842793138395523727_wp
      q1(8,3) = 0.01047289768141453270035870286867846966855_wp
      q1(9,3) = zero
      q1(10,3) = zero
      q1(11,3) = zero
      q1(12,3) = zero

      q1(1,4) = 0.1147228529148158318034738218073197976023_wp
      q1(2,4) = -0.3538780254601914958115698874206531961307_wp
      q1(3,4) = -0.1858994785495101740456777336107169678288_wp
      q1(4,4) = 0.1370483808799853486477767415913728669524_wp
      q1(5,4) = 0.03972858582710281564215102156521588367881_wp
      q1(6,4) = 0.3165092134256211709442618752193670178495_wp
      q1(7,4) = -0.05069213274750316093259254157360772183295_wp
      q1(8,4) = -0.01753939629032033624782329757829768029059_wp
      q1(9,4) = zero
      q1(10,4) = zero
      q1(11,4) = zero
      q1(12,4) = zero

      q1(1,5) = -0.01395722952056041916292538026746562805620_wp
      q1(2,5) = 0.01429909551156929498323180030121632055909_wp
      q1(3,5) = 0.2332633686362482363893930312775209255800_wp
      q1(4,5) = -1.142300242635165935993053851544552454320_wp
      q1(5,5) = 0.6057550738330768177827040634066377989296_wp
      q1(6,5) = 0.2825179352075121763892365873707254541218_wp
      q1(7,5) = 0.03248906901967492013975174489420235738435_wp
      q1(8,5) = -0.01206707005235509052833799543828477419916_wp
      q1(9,5) = zero
      q1(10,5) = zero
      q1(11,5) = zero
      q1(12,5) = zero

      q1(1,6) = -0.06956044896551484587496684695381321009857_wp
      q1(2,6) = 0.1930764794969705526456194423505107491362_wp
      q1(3,6) = -0.01938619631973248083178041998353003089641_wp
      q1(4,6) = -0.2887164590533642022788579933885334942677_wp
      q1(5,6) = -0.4040701147451219074780049288053335780831_wp
      q1(6,6) = 0.1955440022503079624043376625896000853565_wp
      q1(7,6) = 0.4088327186474708718613662440475307694516_wp
      q1(8,6) = -0.02316927663483720615835459462174938908831_wp
      q1(9,6) = 0.007449295323821255710641434765318098489772_wp
      q1(10,6) = zero
      q1(11,6) = zero
      q1(12,6) = zero

      q1(1,7) = 0.01734328280101122760645584447264624051821_wp
      q1(2,7) = -0.04840773491893482968357096782164797435370_wp
      q1(3,7) = 0.01113356043999974094579726795398111697087_wp
      q1(4,7) = 0.03683258325116868830333840353178374785556_wp
      q1(5,7) = 0.2020929196538139300412937632456998385290_wp
      q1(6,7) = -0.9993401510385675084229179036326626049945_wp
      q1(7,7) = 0.2707691940348400441134114620253070105404_wp
      q1(8,7) = 0.6075690064749917706800934212435751620821_wp
      q1(9,7) = -0.1083076776139360176453645848101228042161_wp
      q1(10,7) = 0.01031501691561295406146329379144026706820_wp
      q1(11,7) = zero
      q1(12,7) = zero

      q1(1,8) = 0.008138876723484902470545040666116789759842_wp
      q1(2,8) = -0.02002921745130282774392822052620050919550_wp
      q1(3,8) = -0.002671343688148396036480262140129972226295_wp
      q1(4,8) = 0.03832298452325717814060194729385919629895_wp
      q1(5,8) = -0.05168030333817646678101782026786038350804_wp
      q1(6,8) = 0.2274981280312773565363597229555297756596_wp
      q1(7,8) = -0.9520792032950265164134838698189786312337_wp
      q1(8,8) = 0.2476881136110866797237534905735256806167_wp
      q1(9,8) = 0.5944514726666080313370083773764616334802_wp
      q1(10,8) = -0.09907524544443467188950139622941027224670_wp
      q1(11,8) = 0.009435737661374730656142990117086692594923_wp
      q1(12,8) = zero

      a2(-4) = zero
      a2(-3) = -0.009523809523809523809523809523809523809524_wp
      a2(-2) = 0.1000000000000000000000000000000000000000_wp
      a2(-1) = -0.6000000000000000000000000000000000000000_wp
      a2(0) = -0.2500000000000000000000000000000000000000_wp
      a2(1) = 1.000000000000000000000000000000000000000_wp
      a2(2) = -0.3000000000000000000000000000000000000000_wp
      a2(3) = 0.06666666666666666666666666666666666666667_wp
      a2(4) = -0.007142857142857142857142857142857142857143_wp

      q2(1,1) = -1.707654630177840640955211766064430603569_wp
      q2(2,1) = 2.303083117541551431369172467319435524597_wp
      q2(3,1) = -0.1348499088280232000438498689421515930856_wp
      q2(4,1) = -0.6995295218380085024246152979750969034639_wp
      q2(5,1) = 0.01953354769232648285258901385493084434991_wp
      q2(6,1) = 0.3015760468426855559502835132688146245819_wp
      q2(7,1) = -0.05430139927462465206917597295687512552949_wp
      q2(8,1) = -0.02785725195806647467919208850462676788091_wp
      q2(9,1) = zero
      q2(10,1) = zero
      q2(11,1) = zero
      q2(12,1) = zero

      q2(1,2) = -0.4264124762187621886561567415382508860795_wp
      q2(2,2) = -0.03979384217746298000464189244717138065595_wp
      q2(3,2) = 0.1722618883119834259798412659855476550988_wp
      q2(4,2) = 0.4170573014959035225484918739916436907718_wp
      q2(5,2) = -0.003867911230845480267389686373176021508129_wp
      q2(6,2) = -0.1617893251467682924382748856749713135863_wp
      q2(7,2) = 0.02929412303865500620649355002594213960776_wp
      q2(8,2) = 0.01325024192729698663163651603043611635147_wp
      q2(9,2) = zero
      q2(10,2) = zero
      q2(11,2) = zero
      q2(12,2) = zero

      q2(1,3) = -0.01200671010762254250338646347193566882227_wp
      q2(2,3) = -0.2440194743406570625021278824861723494614_wp
      q2(3,3) = -0.7352247001589518809305910119936419247628_wp
      q2(4,3) = 1.298367310542910198178397511810349706445_wp
      q2(5,3) = -0.3739312333355250790799293893971175134799_wp
      q2(6,3) = 0.09626986570856218757834718746890783658430_wp
      q2(7,3) = -0.03992795599013035344106865479906855617172_wp
      q2(8,3) = 0.01047289768141453270035870286867846966855_wp
      q2(9,3) = zero
      q2(10,3) = zero
      q2(11,3) = zero
      q2(12,3) = zero

      q2(1,4) = 0.1306125202632199301974189512671891155098_wp
      q2(2,4) = -0.4412711958764140369782680994499344446221_wp
      q2(3,4) = 0.02066619697974310507560894936758416496903_wp
      q2(4,4) = -0.1370483808799853486477767415913728669524_wp
      q2(5,4) = 0.2621839287047601931573828340033863343841_wp
      q2(6,4) = 0.2052815419867924821866459690002817924969_wp
      q2(7,4) = -0.01891279805069496414470228265386908601791_wp
      q2(8,4) = -0.02151181312742136084630957994326500976747_wp
      q2(9,4) = zero
      q2(10,4) = zero
      q2(11,4) = zero
      q2(12,4) = zero

      q2(1,5) = -0.03126451734436261395671692493622670802562_wp
      q2(2,5) = 0.1527573981019868533335641576513049603144_wp
      q2(3,5) = -0.2513406904302132178367702194477893135637_wp
      q2(4,5) = -0.1730921245022430275407273500939319760321_wp
      q2(5,5) = -0.6057550738330768177827040634066377989296_wp
      q2(6,5) = 1.251726053340435084841563088821345932409_wp
      q2(7,5) = -0.4521149900467865340864115058311078817594_wp
      q2(8,5) = 0.1263912325380624678219943619118038655562_wp
      q2(9,5) = -0.01730728782380219479379154466876107996942_wp
      q2(10,5) = zero
      q2(11,5) = zero
      q2(12,5) = zero

      q2(1,6) = -0.06956044896551484587496684695381321009857_wp
      q2(2,6) = 0.1874895080041046108626383662765221752689_wp
      q2(3,6) = 0.02530957562319505343206818860837856004222_wp
      q2(4,6) = -0.4451516608536105722023281234602135625529_wp
      q2(5,6) = -0.09119971114462916763106466866197344151272_wp
      q2(6,6) = -0.1955440022503079624043376625896000853565_wp
      q2(7,6) = 0.7217031222479636117083065041908909060220_wp
      q2(8,6) = -0.1796044784350835760818247246934294573735_wp
      q2(9,6) = 0.05214506726674878997449004335722668942840_wp
      q2(10,6) = -0.005586971492865941782981076073988573867329_wp
      q2(11,6) = zero
      q2(12,6) = zero

      q2(1,7) = 0.01734328280101122760645584447264624051821_wp
      q2(2,7) = -0.04840773491893482968357096782164797435370_wp
      q2(3,7) = 0.003397297753290025399699797610400916669720_wp
      q2(4,7) = 0.09872268474484641267211816628042535026479_wp
      q2(5,7) = -0.01452243557405810524943540637454576990329_wp
      q2(6,7) = -0.5661094405828234378414595643921713881299_wp
      q2(7,7) = -0.2707691940348400441134114620253070105404_wp
      q2(8,7) = 1.040799716930735841261551760484066378947_wp
      q2(9,7) = -0.3249230328418080529360937544303684126484_wp
      q2(10,7) = 0.07220511840929067843024305654008186947743_wp
      q2(11,7) = -0.007736262686709715546097470343580200301153_wp
      q2(12,7) = zero

      q2(1,8) = 0.008138876723484902470545040666116789759842_wp
      q2(2,8) = -0.02002921745130282774392822052620050919550_wp
      q2(3,8) = -0.002671343688148396036480262140129972226295_wp
      q2(4,8) = 0.03124618127722613014849470470604417685276_wp
      q2(5,8) = 0.004934122630071917155840120434659772061498_wp
      q2(6,8) = 0.02934763714240801275735693049670923116616_wp
      q2(7,8) = -0.5557782215172878288554782849013375422469_wp
      q2(8,8) = -0.2476881136110866797237534905735256806167_wp
      q2(9,8) = 0.9907524544443467188950139622941027224670_wp
      q2(10,8) = -0.2972257363333040156685041886882308167401_wp
      q2(11,8) = 0.06605016362962311459300093081960684816446_wp
      q2(12,8) = -0.007076803246031047992107242587815019446193_wp

    end subroutine coeffs_up_8_4_opt

    subroutine coeffs_1_4_3 ( a, q )

      CCTK_REAL, dimension(2), intent(OUT) :: a
      CCTK_REAL, dimension(7,5), intent(OUT) :: q
      CCTK_REAL :: f0, f1, f2, f3

      a(1) = 2.0_wp/3.0_wp
      a(2) = -1.0_wp/12.0_wp

      f0 = sqrt(26116897.0_wp)
      f1 = -56764003702447356523.0_wp + 8154993476273221.0_wp * f0
      f2 = -55804550303.0_wp + 9650225.0_wp * f0
      f3 = 3262210757.0_wp + 271861.0_wp * f0
  
      q(1,1) = -11._wp/6.0_wp
      q(2,1) = 3.0_wp
      q(3,1) = -3.0_wp/2.0_wp
      q(4,1) = 1.0_wp/3.0_wp
      q(5,1) = zero
      q(6,1) = zero
      q(7,1) = zero
  
      q(1,2) =  -24.0_wp * ( -779042810827742869.0_wp + &
                              104535124033147.0_wp * f0 ) / f1
      q(2,2) = -( -176530817412806109689.0_wp + &
                   29768274816875927.0_wp * f0 ) / ( 6.0_wp * f1 )
      q(3,2) = 343.0_wp * ( -171079116122226871.0_wp + &
                             27975630462649.0_wp * f0 ) / f1
      q(4,2) = -3.0_wp * ( -7475554291248533227.0_wp + &
                            1648464218793925.0_wp * f0 ) / ( 2.0_wp * f1 )
      q(5,2) = ( -2383792768180030915.0_wp + &
                  1179620587812973.0_wp * f0 ) / ( 3.0_wp * f1 )
      q(6,2) = -1232.0_wp * ( -115724529581315.0_wp + &
                                    37280576429.0_wp * f0 ) / f1
      q(7,2) = zero
  
      q(1,3) = -12.0_wp * ( -380966843.0_wp + 86315.0_wp * f0 ) / f2
      q(2,3) = ( 5024933015.0_wp + 2010631.0_wp * f0 ) / ( 3.0_wp * f2 )
      q(3,3) = -231.0_wp * ( -431968921.0_wp + &
                              86711.0_wp * f0 ) / ( 2.0_wp * f2 )
      q(4,3) = ( -65931742559.0_wp + 12256337.0_wp * f0 ) / f2
      q(5,3) = -( -50597298167.0_wp + 9716873.0_wp * f0 ) / ( 6.0_wp * f2 )
      q(6,3) = -88.0_wp * ( -15453061.0_wp + 2911.0_wp * f0 ) / f2
      q(7,3) = zero
  
      q(1,4) = 48.0_wp * ( -56020909845192541.0_wp + &
                            9790180507043.0_wp * f0 ) / f1
      q(2,4) = ( -9918249049237586011.0_wp + &
                  1463702013196501.0_wp * f0 ) / ( 6.0_wp * f1)
      q(3,4) = -13.0_wp * ( -4130451756851441723.0_wp + &
                             664278707201077.0_wp * f0 ) / f1
      q(4,4) = 3.0_wp * ( -26937108467782666617.0_wp + &
                           5169063172799767.0_wp * f0 ) / ( 2.0_wp * f1 )
      q(5,4) = -( 6548308508012371315.0_wp + &
                  3968886380989379.0_wp * f0 ) / ( 3.0_wp * f1 )
      q(6,4) = 88.0_wp * ( -91337851897923397.0_wp + &
                            19696768305507.0_wp * f0 ) / f1
      q(7,4) = 242.0_wp * ( -120683.0_wp + 15.0_wp * f0 ) / f3
  
      q(1,5) = 264.0_wp * ( -120683.0_wp + 15.0_wp * f0 ) / f3
      q(2,5) = ( -43118111.0_wp + 23357.0_wp * f0 ) / ( 3.0_wp * f3 )
      q(3,5) = -47.0_wp * ( -28770085.0_wp + 2259.0_wp * f0 ) / ( 2.0_wp * f3 )
      q(4,5) = -3.0_wp * ( 1003619433.0_wp + 11777.0_wp * f0 ) / f3
      q(5,5) = -11.0_wp * ( -384168269.0_wp + &
                              65747.0_wp * f0 ) / ( 6.0_wp * f3 )
      q(6,5) = 22.0_wp * ( 87290207.0_wp + 10221.0_wp * f0 ) / f3
      q(7,5) = -66.0_wp * ( 3692405.0_wp + 419.0_wp * f0 ) / f3

    end subroutine coeffs_1_4_3

    subroutine coeffs_1_4_3_opt ( a, q )

      CCTK_REAL, dimension(2), intent(OUT) :: a
      CCTK_REAL, dimension(7,5), intent(OUT) :: q

      a(1) = 2.0_wp/3.0_wp
      a(2) = -1.0_wp/12.0_wp

      q(1,1) = -2.09329763466349871588733_wp
      q(2,1) = 4.0398572053206615302160_wp
      q(3,1) = -3.0597858079809922953240_wp
      q(4,1) = 1.37319053865399486354933_wp
      q(5,1) = -0.25996430133016538255400_wp
      q(6,1) = zero
      q(7,1) = zero

      q(1,2) = -0.31641585285940445272297_wp
      q(2,2) = -0.53930788973980422327388_wp
      q(3,2) = 0.98517732028644343383297_wp
      q(4,2) = -0.05264665989297578146709_wp
      q(5,2) = -0.113807251750624235013258_wp
      q(6,2) = 0.039879767889849911803103_wp
      q(7,2) = -0.0028794339334846531588787_wp

      q(1,3) = 0.13026916185021164524452_wp
      q(2,3) = -0.87966858995059249256890_wp
      q(3,3) = 0.38609640961100070000134_wp
      q(4,3) = 0.31358369072435588745988_wp
      q(5,3) = 0.085318941913678384633511_wp
      q(6,3) = -0.039046615792734640274641_wp
      q(7,3) = 0.0034470016440805155042908_wp

      q(1,4) = -0.01724512193824647912172_wp
      q(2,4) = 0.16272288227127504381134_wp
      q(3,4) = -0.81349810248648813029217_wp
      q(4,4) = 0.13833269266479833215645_wp
      q(5,4) = 0.59743854328548053399616_wp
      q(6,4) = -0.066026434346299887619324_wp
      q(7,4) = -0.0017244594505194129307249_wp

      q(1,5) = -0.00883569468552192965061_wp
      q(2,5) = 0.03056074759203203857284_wp
      q(3,5) = 0.05021168274530854232278_wp
      q(4,5) = -0.66307364652444929534068_wp
      q(5,5) = 0.014878787464005191116088_wp
      q(6,5) = 0.65882706381707471953820_wp
      q(7,5) = -0.082568940408449266558615_wp

    end subroutine coeffs_1_4_3_opt

    subroutine coeffs_1_6_5 ( a, q )

      CCTK_REAL, dimension(3), intent(OUT) :: a
      CCTK_REAL, dimension(10,7), intent(OUT) :: q

      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) = -137.0_wp/60.0_wp
      q(2,1) = 5.0_wp
      q(3,1) = -5.0_wp
      q(4,1) = 10.0_wp/3.0_wp
      q(5,1) = -5.0_wp/4.0_wp
      q(6,1) = 1.0_wp/5.0_wp
      q(7,1) = zero
      q(8,1) = zero
      q(9,1) = zero
      q(10,1) = zero

      q(1,2) = -0.3209098489859668118347_wp
      q(2,2) = -0.359441670715467075894_wp
      q(3,2) = 0.195756852998105503889_wp
      q(4,2) = 1.394685510250317033169_wp
      q(5,2) = -1.448965775497476572821_wp
      q(6,2) = 0.6519476244467816674835_wp
      q(7,2) = -0.1115052611983591304249_wp
      q(8,2) = -0.00156743129793461356831_wp
      q(9,2) = zero
      q(10,2) = zero
  
      q(1,3) = 0.0980796259621024870281_wp
      q(2,3) = -0.784752101846956387607_wp
      q(3,3) = 0.364206221878666751300_wp
      q(4,3) = 0.102097753636344358379_wp
      q(5,3) = 0.377167650934576412896_wp
      q(6,3) = -0.173241400242683302020_wp
      q(7,3) = 0.0062120424243610783650_wp
      q(8,3) = 0.01153111791917461507898_wp
      q(9,3) = -0.00130091066558601341956_wp
      q(10,3) = zero
  
      q(1,4) = 0.0771947815148234674634_wp
      q(2,4) = -0.4124719432622107102831_wp
      q(3,4) = 0.654798717867622763906_wp
      q(4,4) = -1.873474657556531216053_wp
      q(5,4) = 2.161961221935152129759_wp
      q(6,4) = -0.729147814207596325277_wp
      q(7,4) = 0.1292509053479867634674_wp
      q(8,4) = -0.01092898521375837831427_wp
      q(9,4) = 0.00316983426828354261009_wp
      q(10,4) = -0.000352060693772037278516_wp
  
      q(1,5) = -0.01135836398682271847110_wp
      q(2,5) = 0.0114974985688122284189_wp
      q(3,5) = 0.228765565489867047381_wp
      q(4,5) = -1.179122428698534075609_wp
      q(5,5) = 0.774618306056840861215_wp
      q(6,5) = 0.016612403904600041087_wp
      q(7,5) = 0.2399304273448163041220_wp
      q(8,5) = -0.0959180672407791626759_wp
      q(9,5) = 0.01612460763754670343830_wp
      q(10,5) = -0.001149949076347228906108_wp
  
      q(1,6) = -0.01912046962023373131259_wp
      q(2,6) = 0.1569245474161812198493_wp
      q(3,6) = -0.567250839387340928718_wp
      q(4,6) = 1.230411310316952671353_wp
      q(5,6) = -2.048795717924254691221_wp
      q(6,6) = 0.982715833107347313532_wp
      q(7,6) = 0.2876584784911885215053_wp
      q(8,6) = -0.0207657038198929292487_wp
      q(9,6) = -0.00335290957957120811713_wp
      q(10,6) = 0.001575470999623762377395_wp
  
      q(1,7) = -0.0076072322555736921424_wp
      q(2,7) = 0.034503414345314189270_wp
      q(3,7) = -0.044377233426986476366_wp
      q(4,7) = -0.049408595803628143797_wp
      q(5,7) = 0.304693537653175477283_wp
      q(6,7) = -0.935863321236636595622_wp
      q(7,7) = 0.111375967229143365899_wp
      q(8,7) = 0.7149322477536284813995_wp
      q(9,7) = -0.1444768161656941519630_wp
      q(10,7) = 0.01622803190725754603783_wp

    end subroutine coeffs_1_6_5

    subroutine coeffs_1_6_5_opt ( a, q )
    
      CCTK_REAL, dimension(3), intent(OUT) :: a
      CCTK_REAL, dimension(10,7), intent(OUT) :: q

      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) = -2.465354921110524023660777656111276003457_wp
      q(2,1) = 6.092129526663144141964665936667656020742_wp
      q(3,1) = -7.730323816657860354911664841669140051855_wp
      q(4,1) = 6.973765088877147139882219788892186735807_wp
      q(5,1) = -3.980323816657860354911664841669140051855_wp
      q(6,1) = 1.292129526663144141964665936667656020742_wp
      q(7,1) = -0.1820215877771906903274443227779426701237_wp
      q(8,1) = zero
      q(9,1) = zero
      q(10,1) = zero

      q(1,2) = -0.2234725650784319828746535134412736890421_wp
      q(2,2) = -0.9329308121107134563129925525068570679651_wp
      q(3,2) = 1.586820596545839371759081303802027231274_wp
      q(4,2) = -0.3647002340377160216914505558624668821400_wp
      q(5,2) = -0.2666957784872806143914117440166232718819_wp
      q(6,2) = 0.3112949048634705032101261273629794071371_wp
      q(7,2) = -0.1404504214762266650000768489896480092493_wp
      q(8,2) = 0.03488568514730479833596013512958238764128_wp
      q(9,2) = -0.004964021886392518344179263072091597647654_wp
      q(10,2) = 0.0002126465201465853095969115943714918742904_wp

      q(1,3) = 0.1582216737061633151406179477554921935333_wp
      q(2,3) = -1.137049298003377811733609086574457439398_wp
      q(3,3) = 1.212364522932578587741649981040340946798_wp
      q(4,3) = -0.9562288729513894906148167047868730813830_wp
      q(5,3) = 1.066548057336766350478498057851678826640_wp
      q(6,3) = -0.3478788551267041838265477441805600110467_wp
      q(7,3) = -0.03133923293520187620333693909408071632123_wp
      q(8,3) = 0.04098845955755862691072597869183962277781_wp
      q(9,3) = -0.005963188634687155197078928402509551508436_wp
      q(10,3) = 0.0003367341182936373038974376991292099082999_wp

      q(1,4) = 0.02915734641890708196910927068736798144670_wp
      q(2,4) = -0.1169665089768926152768236581512624861308_wp
      q(3,4) = -0.1112219092451476301503253995474190870412_wp
      q(4,4) = -0.7924486261248032107393766820001361351677_wp
      q(5,4) = 1.266650704820613624987450232358951199911_wp
      q(6,4) = -0.2899273290506621673153239836530375587273_wp
      q(7,4) = 0.002515684257201926199329020583484434062150_wp
      q(8,4) = 0.01329713961871764653006682056620518602804_wp
      q(9,4) = -0.001124464399630667352932212208930962568134_wp
      q(10,4) = 0.00006796268169601114882659136477742818715059_wp

      q(1,5) = -0.04582150000326981674750984653096293434777_wp
      q(2,5) = 0.2240986548857151482718685516611524323427_wp
      q(3,5) = -0.3246718493011818141660859125588209338018_wp
      q(4,5) = -0.3929792921782506986152017485694441380503_wp
      q(5,5) = 0.1166355818729375628072830916953646214341_wp
      q(6,5) = 0.3449626905957060254933930895775644438105_wp
      q(7,5) = 0.1430419813354607083034935179267283951745_wp
      q(8,5) = -0.07764802499372607792980458731991885121073_wp
      q(9,5) = 0.01332439335504217034559288889042994978834_wp
      q(10,5) = -0.0009426355684332077630290447720929851395193_wp

      q(1,6) = 0.003172814452954821196677290327889903944225_wp
      q(2,6) = 0.00001061446045061551877105554145609103530766_wp
      q(3,6) = -0.08747763580209736614983637747947172321794_wp
      q(4,6) = 0.3975827322299876034907453299884380895682_wp
      q(5,6) = -1.148835072393422871630425744497391344782_wp
      q(6,6) = 0.3583006649535242306065761818925080902380_wp
      q(7,6) = 0.5647665154270147564019144982190032455071_wp
      q(8,6) = -0.09698196887272109736153117076061707705561_wp
      q(9,6) = 0.008843905091972988427261446924164441884143_wp
      q(10,6) = 0.0006174304523363194998474898440202828786385_wp

      q(1,7) = -0.008639107540858839028043929986084287776394_wp
      q(2,7) = 0.04722773954485212324714352753530343274219_wp
      q(3,7) = -0.1008747537650261142294540111407681552350_wp
      q(4,7) = 0.08043834953845218736895768965086958762389_wp
      q(5,7) = 0.1295138674713300902982857323205417604553_wp
      q(6,7) = -0.7909424166489541737614153656634872155367_wp
      q(7,7) = 0.03807866847647628589685997987877954466259_wp
      q(8,7) = 0.7367055699548196242687865288427927434250_wp
      q(9,7) = -0.1480235854665196220062411065981933720158_wp
      q(10,7) = 0.01651566843542843794512095516024596165494_wp

    end subroutine coeffs_1_6_5_opt

end module All_Coeffs_mod