/* 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" /* 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)) extern "C" void burgers_reconstruct_1_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 */, "Burgers::uLeft_group","flat"); if (ierr < 0) CCTK_WARN(1, "Failed to register flat BC for Burgers::uLeft_group."); ierr = Boundary_SelectGroupForBC(cctkGH, CCTK_ALL_FACES, GenericFD_GetBoundaryWidth(cctkGH), -1 /* no table */, "Burgers::uR_group","flat"); if (ierr < 0) CCTK_WARN(1, "Failed to register flat BC for Burgers::uR_group."); return; } static void burgers_reconstruct_1_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 the variables used for looping over grid points */ CCTK_INT i, j, k; // CCTK_INT index = INITVALUE; /* Declare finite differencing variables */ if (verbose > 1) { CCTK_VInfo(CCTK_THORNSTRING,"Entering burgers_reconstruct_1_Body"); } if (cctk_iteration % burgers_reconstruct_1_calc_every != burgers_reconstruct_1_calc_offset) { return; } const char *groups[] = {"Burgers::u_group","Burgers::uLeft_group","Burgers::uR_group"}; GenericFD_AssertGroupStorage(cctkGH, "burgers_reconstruct_1", 3, groups); GenericFD_EnsureStencilFits(cctkGH, "burgers_reconstruct_1", 1, 1, 1); /* 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); 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 p1o1 = 1; CCTK_REAL const p1odx = INV(dx); CCTK_REAL const p1ody = INV(dy); CCTK_REAL const p1odz = INV(dz); /* Loop over the grid points */ for (k = min[2]; k < max[2]; k++) { for (j = min[1]; j < max[1]; j++) { for (i = min[0]; i < max[0]; i++) { int const index = CCTK_GFINDEX3D(cctkGH,i,j,k) ; /* Assign local copies of grid functions */ CCTK_REAL uL = u[index]; /* Include user supplied include files */ /* Precompute derivatives */ CCTK_REAL const DiffPlus1u = DiffPlus1(&u[index]); CCTK_REAL const DiffMinus1u = DiffMinus1(&u[index]); /* Calculate temporaries and grid functions */ CCTK_REAL slopeL = DiffMinus1u; CCTK_REAL slopeR = DiffPlus1u; CCTK_REAL slope = IfThen(slopeL*slopeR < 0,0,IfThen(Abs(slopeL) < Abs(slopeR),slopeL,slopeR)); CCTK_REAL uLeftL = -0.5*slope + uL; CCTK_REAL uRL = 0.5*slope + uL; /* Copy local copies back to grid functions */ uLeft[index] = uLeftL; uR[index] = uRL; } } } } extern "C" void burgers_reconstruct_1(CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; GenericFD_LoopOverInterior(cctkGH, &burgers_reconstruct_1_Body); }