aboutsummaryrefslogtreecommitdiff
path: root/ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@cct.lsu.edu>2010-05-03 12:55:21 -0500
committerErik Schnetter <schnetter@cct.lsu.edu>2010-05-03 12:55:21 -0500
commit376eca95c980bdead0f6c12d8c47a73312cdbfae (patch)
tree98d5cd055250c7499501b7d8aa16fd61eac2591b /ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c
parentcd01945fcf1737cf09e03cfce0583803ed137960 (diff)
Regenerate code
Diffstat (limited to 'ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c')
-rw-r--r--ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c151
1 files changed, 151 insertions, 0 deletions
diff --git a/ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c b/ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c
new file mode 100644
index 0000000..3b57c9d
--- /dev/null
+++ b/ML_BSSN_O2/src/ML_BSSN_O2_convertToADMBaseFakeDtLapseShift.c
@@ -0,0 +1,151 @@
+/* File produced by Kranc */
+
+#define KRANC_C
+
+#include <assert.h>
+#include <math.h>
+#include <stdio.h>
+#include <stdlib.h>
+#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_O2_convertToADMBaseFakeDtLapseShift_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_O2_convertToADMBaseFakeDtLapseShift_Body");
+ }
+
+ if (cctk_iteration % ML_BSSN_O2_convertToADMBaseFakeDtLapseShift_calc_every != ML_BSSN_O2_convertToADMBaseFakeDtLapseShift_calc_offset)
+ {
+ return;
+ }
+
+ /* 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 p1o16dx = INV(dx)/16.;
+ CCTK_REAL const p1o16dy = INV(dy)/16.;
+ CCTK_REAL const p1o16dz = INV(dz)/16.;
+ CCTK_REAL const p1o2dx = khalf*INV(dx);
+ CCTK_REAL const p1o2dy = khalf*INV(dy);
+ CCTK_REAL const p1o2dz = khalf*INV(dz);
+ CCTK_REAL const p1o4dx = INV(dx)/4.;
+ CCTK_REAL const p1o4dxdy = (INV(dx)*INV(dy))/4.;
+ CCTK_REAL const p1o4dxdz = (INV(dx)*INV(dz))/4.;
+ CCTK_REAL const p1o4dy = INV(dy)/4.;
+ CCTK_REAL const p1o4dydz = (INV(dy)*INV(dz))/4.;
+ CCTK_REAL const p1o4dz = INV(dz)/4.;
+ CCTK_REAL const p1odx = INV(dx);
+ CCTK_REAL const p1odx2 = pow(dx,-2);
+ CCTK_REAL const p1ody = INV(dy);
+ CCTK_REAL const p1ody2 = pow(dy,-2);
+ CCTK_REAL const p1odz = INV(dz);
+ CCTK_REAL const p1odz2 = pow(dz,-2);
+ CCTK_REAL const pm1o2dx = -(khalf*INV(dx));
+ CCTK_REAL const pm1o2dy = -(khalf*INV(dy));
+ CCTK_REAL const pm1o2dz = -(khalf*INV(dz));
+ CCTK_REAL const pm1o4dx = -INV(dx)/4.;
+ CCTK_REAL const pm1o4dy = -INV(dy)/4.;
+ CCTK_REAL const pm1o4dz = -INV(dz)/4.;
+
+ /* Loop over the grid points */
+ #pragma omp parallel
+ LC_LOOP3 (ML_BSSN_O2_convertToADMBaseFakeDtLapseShift,
+ 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_O2_convertToADMBaseFakeDtLapseShift);
+}
+
+void ML_BSSN_O2_convertToADMBaseFakeDtLapseShift(CCTK_ARGUMENTS)
+{
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+ GenericFD_LoopOverEverything(cctkGH, &ML_BSSN_O2_convertToADMBaseFakeDtLapseShift_Body);
+}