/* File produced by Kranc */ #define KRANC_C #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 INV(x) ((1.0) / (x)) #define SQR(x) ((x) * (x)) #define CUB(x) ((x) * (x) * (x)) #define QAD(x) ((x) * (x) * (x) * (x)) void ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary_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 */, "ADMBase::dtlapse","flat"); if (ierr < 0) CCTK_WARN(1, "Failed to register flat BC for ADMBase::dtlapse."); ierr = Boundary_SelectGroupForBC(cctkGH, CCTK_ALL_FACES, GenericFD_GetBoundaryWidth(cctkGH), -1 /* no table */, "ADMBase::dtshift","flat"); if (ierr < 0) CCTK_WARN(1, "Failed to register flat BC for ADMBase::dtshift."); return; } void ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary_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_UPW_convertToADMBaseDtLapseShiftBoundary_Body"); } if (cctk_iteration % ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary_calc_every != ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary_calc_offset) { return; } const char *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_shift","ML_BSSN_UPW::ML_trace_curv"}; GenericFD_AssertGroupStorage(cctkGH, "ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary", 10, groups); /* Include user-supplied include files */ /* Initialise finite differencing variables */ CCTK_REAL const dx = CCTK_DELTA_SPACE(0); CCTK_REAL const dy = CCTK_DELTA_SPACE(1); CCTK_REAL const dz = CCTK_DELTA_SPACE(2); int const di = 1; int const dj = CCTK_GFINDEX3D(cctkGH,0,1,0) - CCTK_GFINDEX3D(cctkGH,0,0,0); int const dk = CCTK_GFINDEX3D(cctkGH,0,0,1) - CCTK_GFINDEX3D(cctkGH,0,0,0); CCTK_REAL const dxi = 1.0 / dx; CCTK_REAL const dyi = 1.0 / dy; CCTK_REAL const dzi = 1.0 / 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 p1o12dx = INV(dx)/12.; CCTK_REAL const p1o12dy = INV(dy)/12.; CCTK_REAL const p1o12dz = INV(dz)/12.; CCTK_REAL const p1o144dxdy = (INV(dx)*INV(dy))/144.; CCTK_REAL const p1o144dxdz = (INV(dx)*INV(dz))/144.; CCTK_REAL const p1o144dydz = (INV(dy)*INV(dz))/144.; CCTK_REAL const p1o24dx = INV(dx)/24.; CCTK_REAL const p1o24dy = INV(dy)/24.; CCTK_REAL const p1o24dz = INV(dz)/24.; CCTK_REAL const p1o64dx = INV(dx)/64.; CCTK_REAL const p1o64dy = INV(dy)/64.; CCTK_REAL const p1o64dz = INV(dz)/64.; CCTK_REAL const p1odx = INV(dx); CCTK_REAL const p1ody = INV(dy); CCTK_REAL const p1odz = INV(dz); CCTK_REAL const pm1o12dx2 = -pow(dx,-2)/12.; CCTK_REAL const pm1o12dy2 = -pow(dy,-2)/12.; CCTK_REAL const pm1o12dz2 = -pow(dz,-2)/12.; /* Loop over the grid points */ #pragma omp parallel LC_LOOP3 (ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary, i,j,k, min[0],min[1],min[2], max[0],max[1],max[2], cctk_lsh[0],cctk_lsh[1],cctk_lsh[2]) { // int index = INITVALUE; int const index = CCTK_GFINDEX3D(cctkGH,i,j,k); /* Declare derivatives */ /* Assign local copies of grid functions */ CCTK_REAL AL = A[index]; CCTK_REAL alphaL = alpha[index]; CCTK_REAL B1L = B1[index]; CCTK_REAL B2L = B2[index]; CCTK_REAL B3L = B3[index]; CCTK_REAL beta1L = beta1[index]; CCTK_REAL beta2L = beta2[index]; CCTK_REAL beta3L = beta3[index]; CCTK_REAL rL = r[index]; CCTK_REAL trKL = trK[index]; CCTK_REAL Xt1L = Xt1[index]; CCTK_REAL Xt2L = Xt2[index]; CCTK_REAL Xt3L = Xt3[index]; /* Include user supplied include files */ /* Precompute derivatives */ /* Calculate temporaries and grid functions */ CCTK_REAL eta = fmin(1,SpatialBetaDriverRadius*INV(rL)); CCTK_REAL theta = fmin(1,exp(1 - rL*INV(SpatialShiftGammaCoeffRadius))); CCTK_REAL dtalpL = -(harmonicF*(LapseACoeff*(AL - trKL) + trKL)*pow(alphaL,harmonicN)); CCTK_REAL dtbetaxL = ShiftGammaCoeff*theta*(beta1L*BetaDriver*eta*(-1 + ShiftBCoeff) + ShiftBCoeff*(B1L - Xt1L) + Xt1L); CCTK_REAL dtbetayL = ShiftGammaCoeff*theta*(beta2L*BetaDriver*eta*(-1 + ShiftBCoeff) + ShiftBCoeff*(B2L - Xt2L) + Xt2L); CCTK_REAL dtbetazL = ShiftGammaCoeff*theta*(beta3L*BetaDriver*eta*(-1 + ShiftBCoeff) + ShiftBCoeff*(B3L - Xt3L) + Xt3L); /* Copy local copies back to grid functions */ dtalp[index] = dtalpL; dtbetax[index] = dtbetaxL; dtbetay[index] = dtbetayL; dtbetaz[index] = dtbetazL; } LC_ENDLOOP3 (ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary); } void ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary(CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; GenericFD_LoopOverBoundaryWithGhosts(cctkGH, &ML_BSSN_UPW_convertToADMBaseDtLapseShiftBoundary_Body); }