aboutsummaryrefslogtreecommitdiff
path: root/ML_BSSN_UPW/src/ML_BSSN_UPW_convertFromADMBaseGamma.cc
blob: 332de3c9467e7262c1c9c4cf0bad5334988a9ac1 (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
/*  File produced by Kranc */

#define KRANC_C

#include <assert.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "GenericFD.h"
#include "Differencing.h"
#include "cctk_Loop.h"
#include "loopcontrol.h"
#include "vectors.h"

/* Define macros used in calculations */
#define INITVALUE (42)
#define QAD(x) (SQR(SQR(x)))
#define INV(x) (kdiv(ToReal(1.0),x))
#define SQR(x) (kmul(x,x))
#define CUB(x) (kmul(x,SQR(x)))

extern "C" void ML_BSSN_UPW_convertFromADMBaseGamma_SelectBCs(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_INT ierr = 0;
  ierr = Boundary_SelectGroupForBC(cctkGH, CCTK_ALL_FACES, GenericFD_GetBoundaryWidth(cctkGH), -1 /* no table */, "ML_BSSN_UPW::ML_dtlapse","flat");
  if (ierr < 0)
    CCTK_WARN(1, "Failed to register flat BC for ML_BSSN_UPW::ML_dtlapse.");
  ierr = Boundary_SelectGroupForBC(cctkGH, CCTK_ALL_FACES, GenericFD_GetBoundaryWidth(cctkGH), -1 /* no table */, "ML_BSSN_UPW::ML_dtshift","flat");
  if (ierr < 0)
    CCTK_WARN(1, "Failed to register flat BC for ML_BSSN_UPW::ML_dtshift.");
  ierr = Boundary_SelectGroupForBC(cctkGH, CCTK_ALL_FACES, GenericFD_GetBoundaryWidth(cctkGH), -1 /* no table */, "ML_BSSN_UPW::ML_Gamma","flat");
  if (ierr < 0)
    CCTK_WARN(1, "Failed to register flat BC for ML_BSSN_UPW::ML_Gamma.");
  return;
}

static void ML_BSSN_UPW_convertFromADMBaseGamma_Body(cGH const * restrict const cctkGH, int const dir, int const face, CCTK_REAL const normal[3], CCTK_REAL const tangentA[3], CCTK_REAL const tangentB[3], int const imin[3], int const imax[3], int const n_subblock_gfs, CCTK_REAL * restrict const subblock_gfs[])
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  
  /* Include user-supplied include files */
  
  /* Initialise finite differencing variables */
  ptrdiff_t const di = 1;
  ptrdiff_t const dj = CCTK_GFINDEX3D(cctkGH,0,1,0) - CCTK_GFINDEX3D(cctkGH,0,0,0);
  ptrdiff_t const dk = CCTK_GFINDEX3D(cctkGH,0,0,1) - CCTK_GFINDEX3D(cctkGH,0,0,0);
  ptrdiff_t const cdi = sizeof(CCTK_REAL) * di;
  ptrdiff_t const cdj = sizeof(CCTK_REAL) * dj;
  ptrdiff_t const cdk = sizeof(CCTK_REAL) * dk;
  CCTK_REAL_VEC const dx = ToReal(CCTK_DELTA_SPACE(0));
  CCTK_REAL_VEC const dy = ToReal(CCTK_DELTA_SPACE(1));
  CCTK_REAL_VEC const dz = ToReal(CCTK_DELTA_SPACE(2));
  CCTK_REAL_VEC const dt = ToReal(CCTK_DELTA_TIME);
  CCTK_REAL_VEC const t = ToReal(cctk_time);
  CCTK_REAL_VEC const dxi = INV(dx);
  CCTK_REAL_VEC const dyi = INV(dy);
  CCTK_REAL_VEC const dzi = INV(dz);
  CCTK_REAL_VEC const khalf = ToReal(0.5);
  CCTK_REAL_VEC const kthird = ToReal(1.0/3.0);
  CCTK_REAL_VEC const ktwothird = ToReal(2.0/3.0);
  CCTK_REAL_VEC const kfourthird = ToReal(4.0/3.0);
  CCTK_REAL_VEC const keightthird = ToReal(8.0/3.0);
  CCTK_REAL_VEC const hdxi = kmul(ToReal(0.5), dxi);
  CCTK_REAL_VEC const hdyi = kmul(ToReal(0.5), dyi);
  CCTK_REAL_VEC const hdzi = kmul(ToReal(0.5), dzi);
  
  /* Initialize predefined quantities */
  CCTK_REAL_VEC const p1o1024dx = kmul(INV(dx),ToReal(0.0009765625));
  CCTK_REAL_VEC const p1o1024dy = kmul(INV(dy),ToReal(0.0009765625));
  CCTK_REAL_VEC const p1o1024dz = kmul(INV(dz),ToReal(0.0009765625));
  CCTK_REAL_VEC const p1o120dx = kmul(INV(dx),ToReal(0.00833333333333333333333333333333));
  CCTK_REAL_VEC const p1o120dy = kmul(INV(dy),ToReal(0.00833333333333333333333333333333));
  CCTK_REAL_VEC const p1o120dz = kmul(INV(dz),ToReal(0.00833333333333333333333333333333));
  CCTK_REAL_VEC const p1o12dx = kmul(INV(dx),ToReal(0.0833333333333333333333333333333));
  CCTK_REAL_VEC const p1o12dy = kmul(INV(dy),ToReal(0.0833333333333333333333333333333));
  CCTK_REAL_VEC const p1o12dz = kmul(INV(dz),ToReal(0.0833333333333333333333333333333));
  CCTK_REAL_VEC const p1o144dxdy = kmul(INV(kmul(dx,dy)),ToReal(0.00694444444444444444444444444444));
  CCTK_REAL_VEC const p1o144dxdz = kmul(INV(kmul(dx,dz)),ToReal(0.00694444444444444444444444444444));
  CCTK_REAL_VEC const p1o144dydz = kmul(INV(kmul(dy,dz)),ToReal(0.00694444444444444444444444444444));
  CCTK_REAL_VEC const p1o1680dx = kmul(INV(dx),ToReal(0.000595238095238095238095238095238));
  CCTK_REAL_VEC const p1o1680dy = kmul(INV(dy),ToReal(0.000595238095238095238095238095238));
  CCTK_REAL_VEC const p1o1680dz = kmul(INV(dz),ToReal(0.000595238095238095238095238095238));
  CCTK_REAL_VEC const p1o16dx = kmul(INV(dx),ToReal(0.0625));
  CCTK_REAL_VEC const p1o16dy = kmul(INV(dy),ToReal(0.0625));
  CCTK_REAL_VEC const p1o16dz = kmul(INV(dz),ToReal(0.0625));
  CCTK_REAL_VEC const p1o180dx2 = kmul(INV(SQR(dx)),ToReal(0.00555555555555555555555555555556));
  CCTK_REAL_VEC const p1o180dy2 = kmul(INV(SQR(dy)),ToReal(0.00555555555555555555555555555556));
  CCTK_REAL_VEC const p1o180dz2 = kmul(INV(SQR(dz)),ToReal(0.00555555555555555555555555555556));
  CCTK_REAL_VEC const p1o24dx = kmul(INV(dx),ToReal(0.0416666666666666666666666666667));
  CCTK_REAL_VEC const p1o24dy = kmul(INV(dy),ToReal(0.0416666666666666666666666666667));
  CCTK_REAL_VEC const p1o24dz = kmul(INV(dz),ToReal(0.0416666666666666666666666666667));
  CCTK_REAL_VEC const p1o256dx = kmul(INV(dx),ToReal(0.00390625));
  CCTK_REAL_VEC const p1o256dy = kmul(INV(dy),ToReal(0.00390625));
  CCTK_REAL_VEC const p1o256dz = kmul(INV(dz),ToReal(0.00390625));
  CCTK_REAL_VEC const p1o2dx = kmul(INV(dx),ToReal(0.5));
  CCTK_REAL_VEC const p1o2dy = kmul(INV(dy),ToReal(0.5));
  CCTK_REAL_VEC const p1o2dz = kmul(INV(dz),ToReal(0.5));
  CCTK_REAL_VEC const p1o3600dxdy = kmul(INV(kmul(dx,dy)),ToReal(0.000277777777777777777777777777778));
  CCTK_REAL_VEC const p1o3600dxdz = kmul(INV(kmul(dx,dz)),ToReal(0.000277777777777777777777777777778));
  CCTK_REAL_VEC const p1o3600dydz = kmul(INV(kmul(dy,dz)),ToReal(0.000277777777777777777777777777778));
  CCTK_REAL_VEC const p1o4dx = kmul(INV(dx),ToReal(0.25));
  CCTK_REAL_VEC const p1o4dxdy = kmul(INV(kmul(dx,dy)),ToReal(0.25));
  CCTK_REAL_VEC const p1o4dxdz = kmul(INV(kmul(dx,dz)),ToReal(0.25));
  CCTK_REAL_VEC const p1o4dy = kmul(INV(dy),ToReal(0.25));
  CCTK_REAL_VEC const p1o4dydz = kmul(INV(kmul(dy,dz)),ToReal(0.25));
  CCTK_REAL_VEC const p1o4dz = kmul(INV(dz),ToReal(0.25));
  CCTK_REAL_VEC const p1o5040dx2 = kmul(INV(SQR(dx)),ToReal(0.000198412698412698412698412698413));
  CCTK_REAL_VEC const p1o5040dy2 = kmul(INV(SQR(dy)),ToReal(0.000198412698412698412698412698413));
  CCTK_REAL_VEC const p1o5040dz2 = kmul(INV(SQR(dz)),ToReal(0.000198412698412698412698412698413));
  CCTK_REAL_VEC const p1o560dx = kmul(INV(dx),ToReal(0.00178571428571428571428571428571));
  CCTK_REAL_VEC const p1o560dy = kmul(INV(dy),ToReal(0.00178571428571428571428571428571));
  CCTK_REAL_VEC const p1o560dz = kmul(INV(dz),ToReal(0.00178571428571428571428571428571));
  CCTK_REAL_VEC const p1o60dx = kmul(INV(dx),ToReal(0.0166666666666666666666666666667));
  CCTK_REAL_VEC const p1o60dy = kmul(INV(dy),ToReal(0.0166666666666666666666666666667));
  CCTK_REAL_VEC const p1o60dz = kmul(INV(dz),ToReal(0.0166666666666666666666666666667));
  CCTK_REAL_VEC const p1o64dx = kmul(INV(dx),ToReal(0.015625));
  CCTK_REAL_VEC const p1o64dy = kmul(INV(dy),ToReal(0.015625));
  CCTK_REAL_VEC const p1o64dz = kmul(INV(dz),ToReal(0.015625));
  CCTK_REAL_VEC const p1o705600dxdy = kmul(INV(kmul(dx,dy)),ToReal(1.41723356009070294784580498866e-6));
  CCTK_REAL_VEC const p1o705600dxdz = kmul(INV(kmul(dx,dz)),ToReal(1.41723356009070294784580498866e-6));
  CCTK_REAL_VEC const p1o705600dydz = kmul(INV(kmul(dy,dz)),ToReal(1.41723356009070294784580498866e-6));
  CCTK_REAL_VEC const p1o840dx = kmul(INV(dx),ToReal(0.00119047619047619047619047619048));
  CCTK_REAL_VEC const p1o840dy = kmul(INV(dy),ToReal(0.00119047619047619047619047619048));
  CCTK_REAL_VEC const p1o840dz = kmul(INV(dz),ToReal(0.00119047619047619047619047619048));
  CCTK_REAL_VEC const p1odx = INV(dx);
  CCTK_REAL_VEC const p1odx2 = INV(SQR(dx));
  CCTK_REAL_VEC const p1ody = INV(dy);
  CCTK_REAL_VEC const p1ody2 = INV(SQR(dy));
  CCTK_REAL_VEC const p1odz = INV(dz);
  CCTK_REAL_VEC const p1odz2 = INV(SQR(dz));
  CCTK_REAL_VEC const pm1o120dx = kmul(INV(dx),ToReal(-0.00833333333333333333333333333333));
  CCTK_REAL_VEC const pm1o120dy = kmul(INV(dy),ToReal(-0.00833333333333333333333333333333));
  CCTK_REAL_VEC const pm1o120dz = kmul(INV(dz),ToReal(-0.00833333333333333333333333333333));
  CCTK_REAL_VEC const pm1o12dx2 = kmul(INV(SQR(dx)),ToReal(-0.0833333333333333333333333333333));
  CCTK_REAL_VEC const pm1o12dy2 = kmul(INV(SQR(dy)),ToReal(-0.0833333333333333333333333333333));
  CCTK_REAL_VEC const pm1o12dz2 = kmul(INV(SQR(dz)),ToReal(-0.0833333333333333333333333333333));
  CCTK_REAL_VEC const pm1o2dx = kmul(INV(dx),ToReal(-0.5));
  CCTK_REAL_VEC const pm1o2dy = kmul(INV(dy),ToReal(-0.5));
  CCTK_REAL_VEC const pm1o2dz = kmul(INV(dz),ToReal(-0.5));
  CCTK_REAL_VEC const pm1o4dx = kmul(INV(dx),ToReal(-0.25));
  CCTK_REAL_VEC const pm1o4dy = kmul(INV(dy),ToReal(-0.25));
  CCTK_REAL_VEC const pm1o4dz = kmul(INV(dz),ToReal(-0.25));
  CCTK_REAL_VEC const pm1o60dx = kmul(INV(dx),ToReal(-0.0166666666666666666666666666667));
  CCTK_REAL_VEC const pm1o60dy = kmul(INV(dy),ToReal(-0.0166666666666666666666666666667));
  CCTK_REAL_VEC const pm1o60dz = kmul(INV(dz),ToReal(-0.0166666666666666666666666666667));
  
  /* Jacobian variable pointers */
  bool const use_jacobian = (!CCTK_IsFunctionAliased("MultiPatch_GetMap") || MultiPatch_GetMap(cctkGH) != jacobian_identity_map)
                       && strlen(jacobian_group) > 0;
  if (use_jacobian && strlen(jacobian_derivative_group) == 0)
  {
    CCTK_WARN (1, "GenericFD::jacobian_group and GenericFD::jacobian_derivative_group must both be set to valid group names");
  }
  
  CCTK_REAL const *restrict jacobian_ptrs[9];
  if (use_jacobian) GenericFD_GroupDataPointers(cctkGH, jacobian_group,
                                                9, jacobian_ptrs);
  
  CCTK_REAL const *restrict const J11 = use_jacobian ? jacobian_ptrs[0] : 0;
  CCTK_REAL const *restrict const J12 = use_jacobian ? jacobian_ptrs[1] : 0;
  CCTK_REAL const *restrict const J13 = use_jacobian ? jacobian_ptrs[2] : 0;
  CCTK_REAL const *restrict const J21 = use_jacobian ? jacobian_ptrs[3] : 0;
  CCTK_REAL const *restrict const J22 = use_jacobian ? jacobian_ptrs[4] : 0;
  CCTK_REAL const *restrict const J23 = use_jacobian ? jacobian_ptrs[5] : 0;
  CCTK_REAL const *restrict const J31 = use_jacobian ? jacobian_ptrs[6] : 0;
  CCTK_REAL const *restrict const J32 = use_jacobian ? jacobian_ptrs[7] : 0;
  CCTK_REAL const *restrict const J33 = use_jacobian ? jacobian_ptrs[8] : 0;
  
  CCTK_REAL const *restrict jacobian_derivative_ptrs[18];
  if (use_jacobian) GenericFD_GroupDataPointers(cctkGH, jacobian_derivative_group,
                                                18, jacobian_derivative_ptrs);
  
  CCTK_REAL const *restrict const dJ111 = use_jacobian ? jacobian_derivative_ptrs[0] : 0;
  CCTK_REAL const *restrict const dJ112 = use_jacobian ? jacobian_derivative_ptrs[1] : 0;
  CCTK_REAL const *restrict const dJ113 = use_jacobian ? jacobian_derivative_ptrs[2] : 0;
  CCTK_REAL const *restrict const dJ122 = use_jacobian ? jacobian_derivative_ptrs[3] : 0;
  CCTK_REAL const *restrict const dJ123 = use_jacobian ? jacobian_derivative_ptrs[4] : 0;
  CCTK_REAL const *restrict const dJ133 = use_jacobian ? jacobian_derivative_ptrs[5] : 0;
  CCTK_REAL const *restrict const dJ211 = use_jacobian ? jacobian_derivative_ptrs[6] : 0;
  CCTK_REAL const *restrict const dJ212 = use_jacobian ? jacobian_derivative_ptrs[7] : 0;
  CCTK_REAL const *restrict const dJ213 = use_jacobian ? jacobian_derivative_ptrs[8] : 0;
  CCTK_REAL const *restrict const dJ222 = use_jacobian ? jacobian_derivative_ptrs[9] : 0;
  CCTK_REAL const *restrict const dJ223 = use_jacobian ? jacobian_derivative_ptrs[10] : 0;
  CCTK_REAL const *restrict const dJ233 = use_jacobian ? jacobian_derivative_ptrs[11] : 0;
  CCTK_REAL const *restrict const dJ311 = use_jacobian ? jacobian_derivative_ptrs[12] : 0;
  CCTK_REAL const *restrict const dJ312 = use_jacobian ? jacobian_derivative_ptrs[13] : 0;
  CCTK_REAL const *restrict const dJ313 = use_jacobian ? jacobian_derivative_ptrs[14] : 0;
  CCTK_REAL const *restrict const dJ322 = use_jacobian ? jacobian_derivative_ptrs[15] : 0;
  CCTK_REAL const *restrict const dJ323 = use_jacobian ? jacobian_derivative_ptrs[16] : 0;
  CCTK_REAL const *restrict const dJ333 = use_jacobian ? jacobian_derivative_ptrs[17] : 0;
  
  /* Assign local copies of arrays functions */
  
  
  
  /* Calculate temporaries and arrays functions */
  
  /* Copy local copies back to grid functions */
  
  /* Loop over the grid points */
  #pragma omp parallel
  LC_LOOP3VEC(ML_BSSN_UPW_convertFromADMBaseGamma,
    i,j,k, imin[0],imin[1],imin[2], imax[0],imax[1],imax[2],
    cctk_lsh[0],cctk_lsh[1],cctk_lsh[2],
    CCTK_REAL_VEC_SIZE)
  {
    ptrdiff_t const index = di*i + dj*j + dk*k;
    
    /* Assign local copies of grid functions */
    
    CCTK_REAL_VEC alphaL = vec_load(alpha[index]);
    CCTK_REAL_VEC beta1L = vec_load(beta1[index]);
    CCTK_REAL_VEC beta2L = vec_load(beta2[index]);
    CCTK_REAL_VEC beta3L = vec_load(beta3[index]);
    CCTK_REAL_VEC dtalpL = vec_load(dtalp[index]);
    CCTK_REAL_VEC dtbetaxL = vec_load(dtbetax[index]);
    CCTK_REAL_VEC dtbetayL = vec_load(dtbetay[index]);
    CCTK_REAL_VEC dtbetazL = vec_load(dtbetaz[index]);
    CCTK_REAL_VEC gt11L = vec_load(gt11[index]);
    CCTK_REAL_VEC gt12L = vec_load(gt12[index]);
    CCTK_REAL_VEC gt13L = vec_load(gt13[index]);
    CCTK_REAL_VEC gt22L = vec_load(gt22[index]);
    CCTK_REAL_VEC gt23L = vec_load(gt23[index]);
    CCTK_REAL_VEC gt33L = vec_load(gt33[index]);
    CCTK_REAL_VEC rL = vec_load(r[index]);
    
    
    CCTK_REAL_VEC J11L, J12L, J13L, J21L, J22L, J23L, J31L, J32L, J33L;
    
    if (use_jacobian)
    {
      J11L = vec_load(J11[index]);
      J12L = vec_load(J12[index]);
      J13L = vec_load(J13[index]);
      J21L = vec_load(J21[index]);
      J22L = vec_load(J22[index]);
      J23L = vec_load(J23[index]);
      J31L = vec_load(J31[index]);
      J32L = vec_load(J32[index]);
      J33L = vec_load(J33[index]);
    }
    
    /* Include user supplied include files */
    
    /* Precompute derivatives */
    CCTK_REAL_VEC PDupwindNth1alpha;
    CCTK_REAL_VEC PDupwindNth2alpha;
    CCTK_REAL_VEC PDupwindNth3alpha;
    CCTK_REAL_VEC PDupwindNth1beta1;
    CCTK_REAL_VEC PDupwindNth2beta1;
    CCTK_REAL_VEC PDupwindNth3beta1;
    CCTK_REAL_VEC PDupwindNth1beta2;
    CCTK_REAL_VEC PDupwindNth2beta2;
    CCTK_REAL_VEC PDupwindNth3beta2;
    CCTK_REAL_VEC PDupwindNth1beta3;
    CCTK_REAL_VEC PDupwindNth2beta3;
    CCTK_REAL_VEC PDupwindNth3beta3;
    CCTK_REAL_VEC PDstandardNth1gt11;
    CCTK_REAL_VEC PDstandardNth2gt11;
    CCTK_REAL_VEC PDstandardNth3gt11;
    CCTK_REAL_VEC PDstandardNth1gt12;
    CCTK_REAL_VEC PDstandardNth2gt12;
    CCTK_REAL_VEC PDstandardNth3gt12;
    CCTK_REAL_VEC PDstandardNth1gt13;
    CCTK_REAL_VEC PDstandardNth2gt13;
    CCTK_REAL_VEC PDstandardNth3gt13;
    CCTK_REAL_VEC PDstandardNth1gt22;
    CCTK_REAL_VEC PDstandardNth2gt22;
    CCTK_REAL_VEC PDstandardNth3gt22;
    CCTK_REAL_VEC PDstandardNth1gt23;
    CCTK_REAL_VEC PDstandardNth2gt23;
    CCTK_REAL_VEC PDstandardNth3gt23;
    CCTK_REAL_VEC PDstandardNth1gt33;
    CCTK_REAL_VEC PDstandardNth2gt33;
    CCTK_REAL_VEC PDstandardNth3gt33;
    
    switch(fdOrder)
    {
      case 2:
        PDupwindNth1alpha = PDupwindNthfdOrder21(&alpha[index]);
        PDupwindNth2alpha = PDupwindNthfdOrder22(&alpha[index]);
        PDupwindNth3alpha = PDupwindNthfdOrder23(&alpha[index]);
        PDupwindNth1beta1 = PDupwindNthfdOrder21(&beta1[index]);
        PDupwindNth2beta1 = PDupwindNthfdOrder22(&beta1[index]);
        PDupwindNth3beta1 = PDupwindNthfdOrder23(&beta1[index]);
        PDupwindNth1beta2 = PDupwindNthfdOrder21(&beta2[index]);
        PDupwindNth2beta2 = PDupwindNthfdOrder22(&beta2[index]);
        PDupwindNth3beta2 = PDupwindNthfdOrder23(&beta2[index]);
        PDupwindNth1beta3 = PDupwindNthfdOrder21(&beta3[index]);
        PDupwindNth2beta3 = PDupwindNthfdOrder22(&beta3[index]);
        PDupwindNth3beta3 = PDupwindNthfdOrder23(&beta3[index]);
        PDstandardNth1gt11 = PDstandardNthfdOrder21(&gt11[index]);
        PDstandardNth2gt11 = PDstandardNthfdOrder22(&gt11[index]);
        PDstandardNth3gt11 = PDstandardNthfdOrder23(&gt11[index]);
        PDstandardNth1gt12 = PDstandardNthfdOrder21(&gt12[index]);
        PDstandardNth2gt12 = PDstandardNthfdOrder22(&gt12[index]);
        PDstandardNth3gt12 = PDstandardNthfdOrder23(&gt12[index]);
        PDstandardNth1gt13 = PDstandardNthfdOrder21(&gt13[index]);
        PDstandardNth2gt13 = PDstandardNthfdOrder22(&gt13[index]);
        PDstandardNth3gt13 = PDstandardNthfdOrder23(&gt13[index]);
        PDstandardNth1gt22 = PDstandardNthfdOrder21(&gt22[index]);
        PDstandardNth2gt22 = PDstandardNthfdOrder22(&gt22[index]);
        PDstandardNth3gt22 = PDstandardNthfdOrder23(&gt22[index]);
        PDstandardNth1gt23 = PDstandardNthfdOrder21(&gt23[index]);
        PDstandardNth2gt23 = PDstandardNthfdOrder22(&gt23[index]);
        PDstandardNth3gt23 = PDstandardNthfdOrder23(&gt23[index]);
        PDstandardNth1gt33 = PDstandardNthfdOrder21(&gt33[index]);
        PDstandardNth2gt33 = PDstandardNthfdOrder22(&gt33[index]);
        PDstandardNth3gt33 = PDstandardNthfdOrder23(&gt33[index]);
        break;
      
      case 4:
        PDupwindNth1alpha = PDupwindNthfdOrder41(&alpha[index]);
        PDupwindNth2alpha = PDupwindNthfdOrder42(&alpha[index]);
        PDupwindNth3alpha = PDupwindNthfdOrder43(&alpha[index]);
        PDupwindNth1beta1 = PDupwindNthfdOrder41(&beta1[index]);
        PDupwindNth2beta1 = PDupwindNthfdOrder42(&beta1[index]);
        PDupwindNth3beta1 = PDupwindNthfdOrder43(&beta1[index]);
        PDupwindNth1beta2 = PDupwindNthfdOrder41(&beta2[index]);
        PDupwindNth2beta2 = PDupwindNthfdOrder42(&beta2[index]);
        PDupwindNth3beta2 = PDupwindNthfdOrder43(&beta2[index]);
        PDupwindNth1beta3 = PDupwindNthfdOrder41(&beta3[index]);
        PDupwindNth2beta3 = PDupwindNthfdOrder42(&beta3[index]);
        PDupwindNth3beta3 = PDupwindNthfdOrder43(&beta3[index]);
        PDstandardNth1gt11 = PDstandardNthfdOrder41(&gt11[index]);
        PDstandardNth2gt11 = PDstandardNthfdOrder42(&gt11[index]);
        PDstandardNth3gt11 = PDstandardNthfdOrder43(&gt11[index]);
        PDstandardNth1gt12 = PDstandardNthfdOrder41(&gt12[index]);
        PDstandardNth2gt12 = PDstandardNthfdOrder42(&gt12[index]);
        PDstandardNth3gt12 = PDstandardNthfdOrder43(&gt12[index]);
        PDstandardNth1gt13 = PDstandardNthfdOrder41(&gt13[index]);
        PDstandardNth2gt13 = PDstandardNthfdOrder42(&gt13[index]);
        PDstandardNth3gt13 = PDstandardNthfdOrder43(&gt13[index]);
        PDstandardNth1gt22 = PDstandardNthfdOrder41(&gt22[index]);
        PDstandardNth2gt22 = PDstandardNthfdOrder42(&gt22[index]);
        PDstandardNth3gt22 = PDstandardNthfdOrder43(&gt22[index]);
        PDstandardNth1gt23 = PDstandardNthfdOrder41(&gt23[index]);
        PDstandardNth2gt23 = PDstandardNthfdOrder42(&gt23[index]);
        PDstandardNth3gt23 = PDstandardNthfdOrder43(&gt23[index]);
        PDstandardNth1gt33 = PDstandardNthfdOrder41(&gt33[index]);
        PDstandardNth2gt33 = PDstandardNthfdOrder42(&gt33[index]);
        PDstandardNth3gt33 = PDstandardNthfdOrder43(&gt33[index]);
        break;
      
      case 6:
        PDupwindNth1alpha = PDupwindNthfdOrder61(&alpha[index]);
        PDupwindNth2alpha = PDupwindNthfdOrder62(&alpha[index]);
        PDupwindNth3alpha = PDupwindNthfdOrder63(&alpha[index]);
        PDupwindNth1beta1 = PDupwindNthfdOrder61(&beta1[index]);
        PDupwindNth2beta1 = PDupwindNthfdOrder62(&beta1[index]);
        PDupwindNth3beta1 = PDupwindNthfdOrder63(&beta1[index]);
        PDupwindNth1beta2 = PDupwindNthfdOrder61(&beta2[index]);
        PDupwindNth2beta2 = PDupwindNthfdOrder62(&beta2[index]);
        PDupwindNth3beta2 = PDupwindNthfdOrder63(&beta2[index]);
        PDupwindNth1beta3 = PDupwindNthfdOrder61(&beta3[index]);
        PDupwindNth2beta3 = PDupwindNthfdOrder62(&beta3[index]);
        PDupwindNth3beta3 = PDupwindNthfdOrder63(&beta3[index]);
        PDstandardNth1gt11 = PDstandardNthfdOrder61(&gt11[index]);
        PDstandardNth2gt11 = PDstandardNthfdOrder62(&gt11[index]);
        PDstandardNth3gt11 = PDstandardNthfdOrder63(&gt11[index]);
        PDstandardNth1gt12 = PDstandardNthfdOrder61(&gt12[index]);
        PDstandardNth2gt12 = PDstandardNthfdOrder62(&gt12[index]);
        PDstandardNth3gt12 = PDstandardNthfdOrder63(&gt12[index]);
        PDstandardNth1gt13 = PDstandardNthfdOrder61(&gt13[index]);
        PDstandardNth2gt13 = PDstandardNthfdOrder62(&gt13[index]);
        PDstandardNth3gt13 = PDstandardNthfdOrder63(&gt13[index]);
        PDstandardNth1gt22 = PDstandardNthfdOrder61(&gt22[index]);
        PDstandardNth2gt22 = PDstandardNthfdOrder62(&gt22[index]);
        PDstandardNth3gt22 = PDstandardNthfdOrder63(&gt22[index]);
        PDstandardNth1gt23 = PDstandardNthfdOrder61(&gt23[index]);
        PDstandardNth2gt23 = PDstandardNthfdOrder62(&gt23[index]);
        PDstandardNth3gt23 = PDstandardNthfdOrder63(&gt23[index]);
        PDstandardNth1gt33 = PDstandardNthfdOrder61(&gt33[index]);
        PDstandardNth2gt33 = PDstandardNthfdOrder62(&gt33[index]);
        PDstandardNth3gt33 = PDstandardNthfdOrder63(&gt33[index]);
        break;
      
      case 8:
        PDupwindNth1alpha = PDupwindNthfdOrder81(&alpha[index]);
        PDupwindNth2alpha = PDupwindNthfdOrder82(&alpha[index]);
        PDupwindNth3alpha = PDupwindNthfdOrder83(&alpha[index]);
        PDupwindNth1beta1 = PDupwindNthfdOrder81(&beta1[index]);
        PDupwindNth2beta1 = PDupwindNthfdOrder82(&beta1[index]);
        PDupwindNth3beta1 = PDupwindNthfdOrder83(&beta1[index]);
        PDupwindNth1beta2 = PDupwindNthfdOrder81(&beta2[index]);
        PDupwindNth2beta2 = PDupwindNthfdOrder82(&beta2[index]);
        PDupwindNth3beta2 = PDupwindNthfdOrder83(&beta2[index]);
        PDupwindNth1beta3 = PDupwindNthfdOrder81(&beta3[index]);
        PDupwindNth2beta3 = PDupwindNthfdOrder82(&beta3[index]);
        PDupwindNth3beta3 = PDupwindNthfdOrder83(&beta3[index]);
        PDstandardNth1gt11 = PDstandardNthfdOrder81(&gt11[index]);
        PDstandardNth2gt11 = PDstandardNthfdOrder82(&gt11[index]);
        PDstandardNth3gt11 = PDstandardNthfdOrder83(&gt11[index]);
        PDstandardNth1gt12 = PDstandardNthfdOrder81(&gt12[index]);
        PDstandardNth2gt12 = PDstandardNthfdOrder82(&gt12[index]);
        PDstandardNth3gt12 = PDstandardNthfdOrder83(&gt12[index]);
        PDstandardNth1gt13 = PDstandardNthfdOrder81(&gt13[index]);
        PDstandardNth2gt13 = PDstandardNthfdOrder82(&gt13[index]);
        PDstandardNth3gt13 = PDstandardNthfdOrder83(&gt13[index]);
        PDstandardNth1gt22 = PDstandardNthfdOrder81(&gt22[index]);
        PDstandardNth2gt22 = PDstandardNthfdOrder82(&gt22[index]);
        PDstandardNth3gt22 = PDstandardNthfdOrder83(&gt22[index]);
        PDstandardNth1gt23 = PDstandardNthfdOrder81(&gt23[index]);
        PDstandardNth2gt23 = PDstandardNthfdOrder82(&gt23[index]);
        PDstandardNth3gt23 = PDstandardNthfdOrder83(&gt23[index]);
        PDstandardNth1gt33 = PDstandardNthfdOrder81(&gt33[index]);
        PDstandardNth2gt33 = PDstandardNthfdOrder82(&gt33[index]);
        PDstandardNth3gt33 = PDstandardNthfdOrder83(&gt33[index]);
        break;
    }
    
    /* Calculate temporaries and grid functions */
    ptrdiff_t dir1 = Sign(beta1L);
    
    ptrdiff_t dir2 = Sign(beta2L);
    
    ptrdiff_t dir3 = Sign(beta3L);
    
    CCTK_REAL_VEC JacPDstandardNth1gt11;
    CCTK_REAL_VEC JacPDstandardNth1gt12;
    CCTK_REAL_VEC JacPDstandardNth1gt13;
    CCTK_REAL_VEC JacPDstandardNth1gt22;
    CCTK_REAL_VEC JacPDstandardNth1gt23;
    CCTK_REAL_VEC JacPDstandardNth1gt33;
    CCTK_REAL_VEC JacPDstandardNth2gt11;
    CCTK_REAL_VEC JacPDstandardNth2gt12;
    CCTK_REAL_VEC JacPDstandardNth2gt13;
    CCTK_REAL_VEC JacPDstandardNth2gt22;
    CCTK_REAL_VEC JacPDstandardNth2gt23;
    CCTK_REAL_VEC JacPDstandardNth2gt33;
    CCTK_REAL_VEC JacPDstandardNth3gt11;
    CCTK_REAL_VEC JacPDstandardNth3gt12;
    CCTK_REAL_VEC JacPDstandardNth3gt13;
    CCTK_REAL_VEC JacPDstandardNth3gt22;
    CCTK_REAL_VEC JacPDstandardNth3gt23;
    CCTK_REAL_VEC JacPDstandardNth3gt33;
    CCTK_REAL_VEC JacPDupwindNth1alpha;
    CCTK_REAL_VEC JacPDupwindNth1beta1;
    CCTK_REAL_VEC JacPDupwindNth1beta2;
    CCTK_REAL_VEC JacPDupwindNth1beta3;
    CCTK_REAL_VEC JacPDupwindNth2alpha;
    CCTK_REAL_VEC JacPDupwindNth2beta1;
    CCTK_REAL_VEC JacPDupwindNth2beta2;
    CCTK_REAL_VEC JacPDupwindNth2beta3;
    CCTK_REAL_VEC JacPDupwindNth3alpha;
    CCTK_REAL_VEC JacPDupwindNth3beta1;
    CCTK_REAL_VEC JacPDupwindNth3beta2;
    CCTK_REAL_VEC JacPDupwindNth3beta3;
    
    if (use_jacobian)
    {
      JacPDstandardNth1gt11 = 
        kmadd(J11L,PDstandardNth1gt11,kmadd(J21L,PDstandardNth2gt11,kmul(J31L,PDstandardNth3gt11)));
      
      JacPDstandardNth1gt12 = 
        kmadd(J11L,PDstandardNth1gt12,kmadd(J21L,PDstandardNth2gt12,kmul(J31L,PDstandardNth3gt12)));
      
      JacPDstandardNth1gt13 = 
        kmadd(J11L,PDstandardNth1gt13,kmadd(J21L,PDstandardNth2gt13,kmul(J31L,PDstandardNth3gt13)));
      
      JacPDstandardNth1gt22 = 
        kmadd(J11L,PDstandardNth1gt22,kmadd(J21L,PDstandardNth2gt22,kmul(J31L,PDstandardNth3gt22)));
      
      JacPDstandardNth1gt23 = 
        kmadd(J11L,PDstandardNth1gt23,kmadd(J21L,PDstandardNth2gt23,kmul(J31L,PDstandardNth3gt23)));
      
      JacPDstandardNth1gt33 = 
        kmadd(J11L,PDstandardNth1gt33,kmadd(J21L,PDstandardNth2gt33,kmul(J31L,PDstandardNth3gt33)));
      
      JacPDstandardNth2gt11 = 
        kmadd(J12L,PDstandardNth1gt11,kmadd(J22L,PDstandardNth2gt11,kmul(J32L,PDstandardNth3gt11)));
      
      JacPDstandardNth2gt12 = 
        kmadd(J12L,PDstandardNth1gt12,kmadd(J22L,PDstandardNth2gt12,kmul(J32L,PDstandardNth3gt12)));
      
      JacPDstandardNth2gt13 = 
        kmadd(J12L,PDstandardNth1gt13,kmadd(J22L,PDstandardNth2gt13,kmul(J32L,PDstandardNth3gt13)));
      
      JacPDstandardNth2gt22 = 
        kmadd(J12L,PDstandardNth1gt22,kmadd(J22L,PDstandardNth2gt22,kmul(J32L,PDstandardNth3gt22)));
      
      JacPDstandardNth2gt23 = 
        kmadd(J12L,PDstandardNth1gt23,kmadd(J22L,PDstandardNth2gt23,kmul(J32L,PDstandardNth3gt23)));
      
      JacPDstandardNth2gt33 = 
        kmadd(J12L,PDstandardNth1gt33,kmadd(J22L,PDstandardNth2gt33,kmul(J32L,PDstandardNth3gt33)));
      
      JacPDstandardNth3gt11 = 
        kmadd(J13L,PDstandardNth1gt11,kmadd(J23L,PDstandardNth2gt11,kmul(J33L,PDstandardNth3gt11)));
      
      JacPDstandardNth3gt12 = 
        kmadd(J13L,PDstandardNth1gt12,kmadd(J23L,PDstandardNth2gt12,kmul(J33L,PDstandardNth3gt12)));
      
      JacPDstandardNth3gt13 = 
        kmadd(J13L,PDstandardNth1gt13,kmadd(J23L,PDstandardNth2gt13,kmul(J33L,PDstandardNth3gt13)));
      
      JacPDstandardNth3gt22 = 
        kmadd(J13L,PDstandardNth1gt22,kmadd(J23L,PDstandardNth2gt22,kmul(J33L,PDstandardNth3gt22)));
      
      JacPDstandardNth3gt23 = 
        kmadd(J13L,PDstandardNth1gt23,kmadd(J23L,PDstandardNth2gt23,kmul(J33L,PDstandardNth3gt23)));
      
      JacPDstandardNth3gt33 = 
        kmadd(J13L,PDstandardNth1gt33,kmadd(J23L,PDstandardNth2gt33,kmul(J33L,PDstandardNth3gt33)));
      
      JacPDupwindNth1alpha = 
        kmadd(J11L,PDupwindNth1alpha,kmadd(J21L,PDupwindNth2alpha,kmul(J31L,PDupwindNth3alpha)));
      
      JacPDupwindNth1beta1 = 
        kmadd(J11L,PDupwindNth1beta1,kmadd(J21L,PDupwindNth2beta1,kmul(J31L,PDupwindNth3beta1)));
      
      JacPDupwindNth1beta2 = 
        kmadd(J11L,PDupwindNth1beta2,kmadd(J21L,PDupwindNth2beta2,kmul(J31L,PDupwindNth3beta2)));
      
      JacPDupwindNth1beta3 = 
        kmadd(J11L,PDupwindNth1beta3,kmadd(J21L,PDupwindNth2beta3,kmul(J31L,PDupwindNth3beta3)));
      
      JacPDupwindNth2alpha = 
        kmadd(J12L,PDupwindNth1alpha,kmadd(J22L,PDupwindNth2alpha,kmul(J32L,PDupwindNth3alpha)));
      
      JacPDupwindNth2beta1 = 
        kmadd(J12L,PDupwindNth1beta1,kmadd(J22L,PDupwindNth2beta1,kmul(J32L,PDupwindNth3beta1)));
      
      JacPDupwindNth2beta2 = 
        kmadd(J12L,PDupwindNth1beta2,kmadd(J22L,PDupwindNth2beta2,kmul(J32L,PDupwindNth3beta2)));
      
      JacPDupwindNth2beta3 = 
        kmadd(J12L,PDupwindNth1beta3,kmadd(J22L,PDupwindNth2beta3,kmul(J32L,PDupwindNth3beta3)));
      
      JacPDupwindNth3alpha = 
        kmadd(J13L,PDupwindNth1alpha,kmadd(J23L,PDupwindNth2alpha,kmul(J33L,PDupwindNth3alpha)));
      
      JacPDupwindNth3beta1 = 
        kmadd(J13L,PDupwindNth1beta1,kmadd(J23L,PDupwindNth2beta1,kmul(J33L,PDupwindNth3beta1)));
      
      JacPDupwindNth3beta2 = 
        kmadd(J13L,PDupwindNth1beta2,kmadd(J23L,PDupwindNth2beta2,kmul(J33L,PDupwindNth3beta2)));
      
      JacPDupwindNth3beta3 = 
        kmadd(J13L,PDupwindNth1beta3,kmadd(J23L,PDupwindNth2beta3,kmul(J33L,PDupwindNth3beta3)));
    }
    else
    {
      JacPDstandardNth1gt11 = PDstandardNth1gt11;
      
      JacPDstandardNth1gt12 = PDstandardNth1gt12;
      
      JacPDstandardNth1gt13 = PDstandardNth1gt13;
      
      JacPDstandardNth1gt22 = PDstandardNth1gt22;
      
      JacPDstandardNth1gt23 = PDstandardNth1gt23;
      
      JacPDstandardNth1gt33 = PDstandardNth1gt33;
      
      JacPDstandardNth2gt11 = PDstandardNth2gt11;
      
      JacPDstandardNth2gt12 = PDstandardNth2gt12;
      
      JacPDstandardNth2gt13 = PDstandardNth2gt13;
      
      JacPDstandardNth2gt22 = PDstandardNth2gt22;
      
      JacPDstandardNth2gt23 = PDstandardNth2gt23;
      
      JacPDstandardNth2gt33 = PDstandardNth2gt33;
      
      JacPDstandardNth3gt11 = PDstandardNth3gt11;
      
      JacPDstandardNth3gt12 = PDstandardNth3gt12;
      
      JacPDstandardNth3gt13 = PDstandardNth3gt13;
      
      JacPDstandardNth3gt22 = PDstandardNth3gt22;
      
      JacPDstandardNth3gt23 = PDstandardNth3gt23;
      
      JacPDstandardNth3gt33 = PDstandardNth3gt33;
      
      JacPDupwindNth1alpha = PDupwindNth1alpha;
      
      JacPDupwindNth1beta1 = PDupwindNth1beta1;
      
      JacPDupwindNth1beta2 = PDupwindNth1beta2;
      
      JacPDupwindNth1beta3 = PDupwindNth1beta3;
      
      JacPDupwindNth2alpha = PDupwindNth2alpha;
      
      JacPDupwindNth2beta1 = PDupwindNth2beta1;
      
      JacPDupwindNth2beta2 = PDupwindNth2beta2;
      
      JacPDupwindNth2beta3 = PDupwindNth2beta3;
      
      JacPDupwindNth3alpha = PDupwindNth3alpha;
      
      JacPDupwindNth3beta1 = PDupwindNth3beta1;
      
      JacPDupwindNth3beta2 = PDupwindNth3beta2;
      
      JacPDupwindNth3beta3 = PDupwindNth3beta3;
    }
    
    CCTK_REAL_VEC detgt = ToReal(1);
    
    CCTK_REAL_VEC gtu11 = 
      kmul(INV(detgt),kmsub(gt22L,gt33L,SQR(gt23L)));
    
    CCTK_REAL_VEC gtu12 = 
      kmul(INV(detgt),kmsub(gt13L,gt23L,kmul(gt12L,gt33L)));
    
    CCTK_REAL_VEC gtu13 = 
      kmul(INV(detgt),kmsub(gt12L,gt23L,kmul(gt13L,gt22L)));
    
    CCTK_REAL_VEC gtu22 = 
      kmul(INV(detgt),kmsub(gt11L,gt33L,SQR(gt13L)));
    
    CCTK_REAL_VEC gtu23 = 
      kmul(INV(detgt),kmsub(gt12L,gt13L,kmul(gt11L,gt23L)));
    
    CCTK_REAL_VEC gtu33 = 
      kmul(INV(detgt),kmsub(gt11L,gt22L,SQR(gt12L)));
    
    CCTK_REAL_VEC Gt111 = 
      kmul(ToReal(0.5),kmadd(gtu11,JacPDstandardNth1gt11,knmsub(gtu12,JacPDstandardNth2gt11,kmsub(kmadd(gtu12,JacPDstandardNth1gt12,kmul(gtu13,JacPDstandardNth1gt13)),ToReal(2),kmul(gtu13,JacPDstandardNth3gt11)))));
    
    CCTK_REAL_VEC Gt211 = 
      kmul(ToReal(0.5),kmadd(gtu12,JacPDstandardNth1gt11,knmsub(gtu22,JacPDstandardNth2gt11,kmsub(kmadd(gtu22,JacPDstandardNth1gt12,kmul(gtu23,JacPDstandardNth1gt13)),ToReal(2),kmul(gtu23,JacPDstandardNth3gt11)))));
    
    CCTK_REAL_VEC Gt311 = 
      kmul(ToReal(0.5),kmadd(gtu13,JacPDstandardNth1gt11,knmsub(gtu23,JacPDstandardNth2gt11,kmsub(kmadd(gtu23,JacPDstandardNth1gt12,kmul(gtu33,JacPDstandardNth1gt13)),ToReal(2),kmul(gtu33,JacPDstandardNth3gt11)))));
    
    CCTK_REAL_VEC Gt112 = 
      kmul(kmadd(gtu12,JacPDstandardNth1gt22,kmadd(gtu11,JacPDstandardNth2gt11,kmul(gtu13,kadd(JacPDstandardNth1gt23,ksub(JacPDstandardNth2gt13,JacPDstandardNth3gt12))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt212 = 
      kmul(kmadd(gtu22,JacPDstandardNth1gt22,kmadd(gtu12,JacPDstandardNth2gt11,kmul(gtu23,kadd(JacPDstandardNth1gt23,ksub(JacPDstandardNth2gt13,JacPDstandardNth3gt12))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt312 = 
      kmul(kmadd(gtu23,JacPDstandardNth1gt22,kmadd(gtu13,JacPDstandardNth2gt11,kmul(gtu33,kadd(JacPDstandardNth1gt23,ksub(JacPDstandardNth2gt13,JacPDstandardNth3gt12))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt113 = 
      kmul(kmadd(gtu13,JacPDstandardNth1gt33,kmadd(gtu11,JacPDstandardNth3gt11,kmul(gtu12,kadd(JacPDstandardNth1gt23,ksub(JacPDstandardNth3gt12,JacPDstandardNth2gt13))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt213 = 
      kmul(kmadd(gtu23,JacPDstandardNth1gt33,kmadd(gtu12,JacPDstandardNth3gt11,kmul(gtu22,kadd(JacPDstandardNth1gt23,ksub(JacPDstandardNth3gt12,JacPDstandardNth2gt13))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt313 = 
      kmul(kmadd(gtu33,JacPDstandardNth1gt33,kmadd(gtu13,JacPDstandardNth3gt11,kmul(gtu23,kadd(JacPDstandardNth1gt23,ksub(JacPDstandardNth3gt12,JacPDstandardNth2gt13))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt122 = 
      kmul(ToReal(0.5),kmadd(gtu12,JacPDstandardNth2gt22,kmadd(gtu11,kmsub(JacPDstandardNth2gt12,ToReal(2),JacPDstandardNth1gt22),kmul(gtu13,kmsub(JacPDstandardNth2gt23,ToReal(2),JacPDstandardNth3gt22)))));
    
    CCTK_REAL_VEC Gt222 = 
      kmul(ToReal(0.5),kmadd(gtu22,JacPDstandardNth2gt22,kmadd(gtu12,kmsub(JacPDstandardNth2gt12,ToReal(2),JacPDstandardNth1gt22),kmul(gtu23,kmsub(JacPDstandardNth2gt23,ToReal(2),JacPDstandardNth3gt22)))));
    
    CCTK_REAL_VEC Gt322 = 
      kmul(ToReal(0.5),kmadd(gtu23,JacPDstandardNth2gt22,kmadd(gtu13,kmsub(JacPDstandardNth2gt12,ToReal(2),JacPDstandardNth1gt22),kmul(gtu33,kmsub(JacPDstandardNth2gt23,ToReal(2),JacPDstandardNth3gt22)))));
    
    CCTK_REAL_VEC Gt123 = 
      kmul(kmadd(gtu13,JacPDstandardNth2gt33,kmadd(gtu12,JacPDstandardNth3gt22,kmul(gtu11,kadd(JacPDstandardNth2gt13,ksub(JacPDstandardNth3gt12,JacPDstandardNth1gt23))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt223 = 
      kmul(kmadd(gtu23,JacPDstandardNth2gt33,kmadd(gtu22,JacPDstandardNth3gt22,kmul(gtu12,kadd(JacPDstandardNth2gt13,ksub(JacPDstandardNth3gt12,JacPDstandardNth1gt23))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt323 = 
      kmul(kmadd(gtu33,JacPDstandardNth2gt33,kmadd(gtu23,JacPDstandardNth3gt22,kmul(gtu13,kadd(JacPDstandardNth2gt13,ksub(JacPDstandardNth3gt12,JacPDstandardNth1gt23))))),ToReal(0.5));
    
    CCTK_REAL_VEC Gt133 = 
      kmul(ToReal(0.5),kmadd(gtu13,JacPDstandardNth3gt33,kmadd(gtu11,kmsub(JacPDstandardNth3gt13,ToReal(2),JacPDstandardNth1gt33),kmul(gtu12,kmsub(JacPDstandardNth3gt23,ToReal(2),JacPDstandardNth2gt33)))));
    
    CCTK_REAL_VEC Gt233 = 
      kmul(ToReal(0.5),kmadd(gtu23,JacPDstandardNth3gt33,kmadd(gtu12,kmsub(JacPDstandardNth3gt13,ToReal(2),JacPDstandardNth1gt33),kmul(gtu22,kmsub(JacPDstandardNth3gt23,ToReal(2),JacPDstandardNth2gt33)))));
    
    CCTK_REAL_VEC Gt333 = 
      kmul(ToReal(0.5),kmadd(gtu33,JacPDstandardNth3gt33,kmadd(gtu13,kmsub(JacPDstandardNth3gt13,ToReal(2),JacPDstandardNth1gt33),kmul(gtu23,kmsub(JacPDstandardNth3gt23,ToReal(2),JacPDstandardNth2gt33)))));
    
    CCTK_REAL_VEC Xt1L = 
      kmadd(Gt111,gtu11,kmadd(Gt122,gtu22,kmadd(Gt133,gtu33,kmul(kmadd(Gt112,gtu12,kmadd(Gt113,gtu13,kmul(Gt123,gtu23))),ToReal(2)))));
    
    CCTK_REAL_VEC Xt2L = 
      kmadd(Gt211,gtu11,kmadd(Gt222,gtu22,kmadd(Gt233,gtu33,kmul(kmadd(Gt212,gtu12,kmadd(Gt213,gtu13,kmul(Gt223,gtu23))),ToReal(2)))));
    
    CCTK_REAL_VEC Xt3L = 
      kmadd(Gt311,gtu11,kmadd(Gt322,gtu22,kmadd(Gt333,gtu33,kmul(kmadd(Gt312,gtu12,kmadd(Gt313,gtu13,kmul(Gt323,gtu23))),ToReal(2)))));
    
    CCTK_REAL_VEC AL = IfThen(LapseACoeff != 
      0,kmul(INV(ToReal(harmonicF)),kmul(kpow(alphaL,-harmonicN),kmsub(kmadd(beta1L,JacPDupwindNth1alpha,kmadd(beta2L,JacPDupwindNth2alpha,kmul(beta3L,JacPDupwindNth3alpha))),ToReal(LapseAdvectionCoeff),dtalpL))),ToReal(0));
    
    CCTK_REAL_VEC theta = 
      kfmin(ToReal(1),kexp(knmsub(rL,INV(ToReal(SpatialShiftGammaCoeffRadius)),ToReal(1))));
    
    CCTK_REAL_VEC B1L;
    CCTK_REAL_VEC B2L;
    CCTK_REAL_VEC B3L;
    
    if (ShiftBCoeff*ShiftGammaCoeff != 0)
    {
      B1L = 
        kmul(INV(kmul(theta,ToReal(ShiftGammaCoeff))),knmsub(kmadd(beta1L,JacPDupwindNth1beta1,kmadd(beta2L,JacPDupwindNth2beta1,kmul(beta3L,JacPDupwindNth3beta1))),ToReal(ShiftAdvectionCoeff),dtbetaxL));
      
      B2L = 
        kmul(INV(kmul(theta,ToReal(ShiftGammaCoeff))),knmsub(kmadd(beta1L,JacPDupwindNth1beta2,kmadd(beta2L,JacPDupwindNth2beta2,kmul(beta3L,JacPDupwindNth3beta2))),ToReal(ShiftAdvectionCoeff),dtbetayL));
      
      B3L = 
        kmul(INV(kmul(theta,ToReal(ShiftGammaCoeff))),knmsub(kmadd(beta1L,JacPDupwindNth1beta3,kmadd(beta2L,JacPDupwindNth2beta3,kmul(beta3L,JacPDupwindNth3beta3))),ToReal(ShiftAdvectionCoeff),dtbetazL));
    }
    else
    {
      B1L = ToReal(0);
      
      B2L = ToReal(0);
      
      B3L = ToReal(0);
    }
    
    /* Copy local copies back to grid functions */
    vec_store_partial_prepare(i,lc_imin,lc_imax);
    vec_store_nta_partial(A[index],AL);
    vec_store_nta_partial(B1[index],B1L);
    vec_store_nta_partial(B2[index],B2L);
    vec_store_nta_partial(B3[index],B3L);
    vec_store_nta_partial(Xt1[index],Xt1L);
    vec_store_nta_partial(Xt2[index],Xt2L);
    vec_store_nta_partial(Xt3[index],Xt3L);
  }
  LC_ENDLOOP3VEC(ML_BSSN_UPW_convertFromADMBaseGamma);
}

extern "C" void ML_BSSN_UPW_convertFromADMBaseGamma(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  
  if (verbose > 1)
  {
    CCTK_VInfo(CCTK_THORNSTRING,"Entering ML_BSSN_UPW_convertFromADMBaseGamma_Body");
  }
  
  if (cctk_iteration % ML_BSSN_UPW_convertFromADMBaseGamma_calc_every != ML_BSSN_UPW_convertFromADMBaseGamma_calc_offset)
  {
    return;
  }
  
  const char *const groups[] = {
    "ADMBase::dtlapse",
    "ADMBase::dtshift",
    "grid::coordinates",
    "Grid::coordinates",
    "ML_BSSN_UPW::ML_dtlapse",
    "ML_BSSN_UPW::ML_dtshift",
    "ML_BSSN_UPW::ML_Gamma",
    "ML_BSSN_UPW::ML_lapse",
    "ML_BSSN_UPW::ML_metric",
    "ML_BSSN_UPW::ML_shift"};
  GenericFD_AssertGroupStorage(cctkGH, "ML_BSSN_UPW_convertFromADMBaseGamma", 10, groups);
  
  switch(fdOrder)
  {
    case 2:
      GenericFD_EnsureStencilFits(cctkGH, "ML_BSSN_UPW_convertFromADMBaseGamma", 2, 2, 2);
      break;
    
    case 4:
      GenericFD_EnsureStencilFits(cctkGH, "ML_BSSN_UPW_convertFromADMBaseGamma", 3, 3, 3);
      break;
    
    case 6:
      GenericFD_EnsureStencilFits(cctkGH, "ML_BSSN_UPW_convertFromADMBaseGamma", 4, 4, 4);
      break;
    
    case 8:
      GenericFD_EnsureStencilFits(cctkGH, "ML_BSSN_UPW_convertFromADMBaseGamma", 5, 5, 5);
      break;
  }
  
  GenericFD_LoopOverInterior(cctkGH, ML_BSSN_UPW_convertFromADMBaseGamma_Body);
  
  if (verbose > 1)
  {
    CCTK_VInfo(CCTK_THORNSTRING,"Leaving ML_BSSN_UPW_convertFromADMBaseGamma_Body");
  }
}