/* File produced by Kranc */ #define KRANC_C #include #include #include #include #include #include "cctk.h" #include "cctk_Arguments.h" #include "cctk_Parameters.h" #include "GenericFD.h" #include "Differencing.h" #include "loopcontrol.h" /* Define macros used in calculations */ #define INITVALUE (42) #define QAD(x) (SQR(SQR(x))) #define INV(x) ((1.0) / (x)) #define SQR(x) ((x) * (x)) #define CUB(x) ((x) * (x) * (x)) static void ML_BSSN_O8_enforce_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 min[3], int const max[3], int const n_subblock_gfs, CCTK_REAL * restrict const subblock_gfs[]) { DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; /* Declare finite differencing variables */ if (verbose > 1) { CCTK_VInfo(CCTK_THORNSTRING,"Entering ML_BSSN_O8_enforce_Body"); } if (cctk_iteration % ML_BSSN_O8_enforce_calc_every != ML_BSSN_O8_enforce_calc_offset) { return; } const char *groups[] = {"ML_BSSN_O8::ML_curv","ML_BSSN_O8::ML_lapse","ML_BSSN_O8::ML_metric"}; GenericFD_AssertGroupStorage(cctkGH, "ML_BSSN_O8_enforce", 3, groups); /* 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 const dx = ToReal(CCTK_DELTA_SPACE(0)); CCTK_REAL const dy = ToReal(CCTK_DELTA_SPACE(1)); CCTK_REAL const dz = ToReal(CCTK_DELTA_SPACE(2)); CCTK_REAL const dt = ToReal(CCTK_DELTA_TIME); CCTK_REAL const dxi = INV(dx); CCTK_REAL const dyi = INV(dy); CCTK_REAL const dzi = INV(dz); CCTK_REAL const khalf = 0.5; CCTK_REAL const kthird = 1/3.0; CCTK_REAL const ktwothird = 2.0/3.0; CCTK_REAL const kfourthird = 4.0/3.0; CCTK_REAL const keightthird = 8.0/3.0; CCTK_REAL const hdxi = 0.5 * dxi; CCTK_REAL const hdyi = 0.5 * dyi; CCTK_REAL const hdzi = 0.5 * dzi; /* Initialize predefined quantities */ CCTK_REAL const p1o1024dx = 0.0009765625*INV(dx); CCTK_REAL const p1o1024dy = 0.0009765625*INV(dy); CCTK_REAL const p1o1024dz = 0.0009765625*INV(dz); CCTK_REAL const p1o1680dx = 0.000595238095238095238095238095238*INV(dx); CCTK_REAL const p1o1680dy = 0.000595238095238095238095238095238*INV(dy); CCTK_REAL const p1o1680dz = 0.000595238095238095238095238095238*INV(dz); CCTK_REAL const p1o5040dx2 = 0.000198412698412698412698412698413*INV(SQR(dx)); CCTK_REAL const p1o5040dy2 = 0.000198412698412698412698412698413*INV(SQR(dy)); CCTK_REAL const p1o5040dz2 = 0.000198412698412698412698412698413*INV(SQR(dz)); CCTK_REAL const p1o560dx = 0.00178571428571428571428571428571*INV(dx); CCTK_REAL const p1o560dy = 0.00178571428571428571428571428571*INV(dy); CCTK_REAL const p1o560dz = 0.00178571428571428571428571428571*INV(dz); CCTK_REAL const p1o705600dxdy = 1.41723356009070294784580498866e-6*INV(dx)*INV(dy); CCTK_REAL const p1o705600dxdz = 1.41723356009070294784580498866e-6*INV(dx)*INV(dz); CCTK_REAL const p1o705600dydz = 1.41723356009070294784580498866e-6*INV(dy)*INV(dz); CCTK_REAL const p1o840dx = 0.00119047619047619047619047619048*INV(dx); CCTK_REAL const p1o840dy = 0.00119047619047619047619047619048*INV(dy); CCTK_REAL const p1o840dz = 0.00119047619047619047619047619048*INV(dz); CCTK_REAL const p1odx = INV(dx); CCTK_REAL const p1ody = INV(dy); CCTK_REAL const p1odz = INV(dz); CCTK_REAL const pm1o840dx = -0.00119047619047619047619047619048*INV(dx); CCTK_REAL const pm1o840dy = -0.00119047619047619047619047619048*INV(dy); CCTK_REAL const pm1o840dz = -0.00119047619047619047619047619048*INV(dz); /* Loop over the grid points */ #pragma omp parallel LC_LOOP3 (ML_BSSN_O8_enforce, i,j,k, min[0],min[1],min[2], max[0],max[1],max[2], cctk_lsh[0],cctk_lsh[1],cctk_lsh[2]) { ptrdiff_t const index = di*i + dj*j + dk*k; /* Assign local copies of grid functions */ CCTK_REAL alphaL = alpha[index]; CCTK_REAL At11L = At11[index]; CCTK_REAL At12L = At12[index]; CCTK_REAL At13L = At13[index]; CCTK_REAL At22L = At22[index]; CCTK_REAL At23L = At23[index]; CCTK_REAL At33L = At33[index]; CCTK_REAL gt11L = gt11[index]; CCTK_REAL gt12L = gt12[index]; CCTK_REAL gt13L = gt13[index]; CCTK_REAL gt22L = gt22[index]; CCTK_REAL gt23L = gt23[index]; CCTK_REAL gt33L = gt33[index]; /* Include user supplied include files */ /* Precompute derivatives */ /* Calculate temporaries and grid functions */ CCTK_REAL detgt = 1; CCTK_REAL gtu11 = INV(detgt)*(gt22L*gt33L - SQR(gt23L)); CCTK_REAL gtu12 = (gt13L*gt23L - gt12L*gt33L)*INV(detgt); CCTK_REAL gtu13 = (-(gt13L*gt22L) + gt12L*gt23L)*INV(detgt); CCTK_REAL gtu22 = INV(detgt)*(gt11L*gt33L - SQR(gt13L)); CCTK_REAL gtu23 = (gt12L*gt13L - gt11L*gt23L)*INV(detgt); CCTK_REAL gtu33 = INV(detgt)*(gt11L*gt22L - SQR(gt12L)); CCTK_REAL trAt = At11L*gtu11 + At22L*gtu22 + 2*(At12L*gtu12 + At13L*gtu13 + At23L*gtu23) + At33L*gtu33; At11L = At11L - 0.333333333333333333333333333333*gt11L*trAt; At12L = At12L - 0.333333333333333333333333333333*gt12L*trAt; At13L = At13L - 0.333333333333333333333333333333*gt13L*trAt; At22L = At22L - 0.333333333333333333333333333333*gt22L*trAt; At23L = At23L - 0.333333333333333333333333333333*gt23L*trAt; At33L = At33L - 0.333333333333333333333333333333*gt33L*trAt; alphaL = fmax(alphaL,ToReal(MinimumLapse)); /* Copy local copies back to grid functions */ alpha[index] = alphaL; At11[index] = At11L; At12[index] = At12L; At13[index] = At13L; At22[index] = At22L; At23[index] = At23L; At33[index] = At33L; } LC_ENDLOOP3 (ML_BSSN_O8_enforce); } extern "C" void ML_BSSN_O8_enforce(CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; GenericFD_LoopOverEverything(cctkGH, &ML_BSSN_O8_enforce_Body); }