aboutsummaryrefslogtreecommitdiff
path: root/ML_BSSN_MP
diff options
context:
space:
mode:
authorPeter Diener <diener@linux-hn3d.site>2010-01-25 17:46:15 -0600
committerPeter Diener <diener@linux-hn3d.site>2010-01-25 17:46:15 -0600
commit77f234967fbcc6344888c8ef5f84d42adf7d1078 (patch)
treec36ec06e01546a9762fa3a9818489bd39bc7de15 /ML_BSSN_MP
parent2ea07083c273b392d07f19fdac918576ca7bc73e (diff)
Make BetaDriver spatially varying.
Make BetaDriver spatially varying. For now additional storage is unconditional, but that can be changed when we decide on exactly how to do it. Signed-off-by: Peter Diener <diener@linux-hn3d.site>
Diffstat (limited to 'ML_BSSN_MP')
-rw-r--r--ML_BSSN_MP/interface.ccl6
-rw-r--r--ML_BSSN_MP/param.ccl25
-rw-r--r--ML_BSSN_MP/schedule.ccl11
-rw-r--r--ML_BSSN_MP/src/ML_BSSN_MP_Minkowski.c4
-rw-r--r--ML_BSSN_MP/src/ML_BSSN_MP_RHS.c8
-rw-r--r--ML_BSSN_MP/src/ML_BSSN_MP_convertFromADMBase.c4
-rw-r--r--ML_BSSN_MP/src/ML_BSSN_MP_setBetaDriver.c142
-rw-r--r--ML_BSSN_MP/src/RegisterSymmetries.c5
-rw-r--r--ML_BSSN_MP/src/make.code.defn2
9 files changed, 203 insertions, 4 deletions
diff --git a/ML_BSSN_MP/interface.ccl b/ML_BSSN_MP/interface.ccl
index 8573da7..ee13843 100644
--- a/ML_BSSN_MP/interface.ccl
+++ b/ML_BSSN_MP/interface.ccl
@@ -25,6 +25,12 @@ CCTK_INT FUNCTION Boundary_SelectVarForBC(CCTK_POINTER_TO_CONST IN GH, CCTK_INT
USES FUNCTION Boundary_SelectVarForBC
public:
+CCTK_REAL ML_BetaDriver type=GF timelevels=1 tags='tensortypealias="Scalar" tensorweight=1.0000000000000000000'
+{
+ eta
+} "ML_BetaDriver"
+
+public:
CCTK_REAL ML_cons_detg type=GF timelevels=1 tags='tensortypealias="Scalar" tensorweight=2.0000000000000000000'
{
cS
diff --git a/ML_BSSN_MP/param.ccl b/ML_BSSN_MP/param.ccl
index e006904..23ac896 100644
--- a/ML_BSSN_MP/param.ccl
+++ b/ML_BSSN_MP/param.ccl
@@ -97,6 +97,12 @@ CCTK_REAL MinimumLapse "Minimum value of the lapse function"
} -1
restricted:
+CCTK_REAL SpatialBetaDriverRadius "Radius at which the BetaDriver starts to be reduced"
+{
+ "*:*" :: ""
+} 1000000000000
+
+restricted:
CCTK_INT harmonicN "d/dt alpha = - f alpha^n K (harmonic=2, 1+log=1)"
{
*:* :: ""
@@ -156,6 +162,13 @@ KEYWORD calculate_ADMBase_variables_at "calculate_ADMBase_variables_at"
"CCTK_ANALYSIS" :: "CCTK_ANALYSIS"
} "MoL_PostStep"
+restricted:
+KEYWORD UseSpatialBetaDriver "UseSpatialBetaDriver"
+{
+ "no" :: "no"
+ "yes" :: "yes"
+} "no"
+
private:
KEYWORD dt_lapse_shift_method "Treatment of ADMBase dtlapse and dtshift"
{
@@ -200,6 +213,12 @@ CCTK_INT ML_BSSN_MP_convertFromADMBaseGamma_calc_every "ML_BSSN_MP_convertFromAD
} 1
restricted:
+CCTK_INT ML_BSSN_MP_setBetaDriver_calc_every "ML_BSSN_MP_setBetaDriver_calc_every"
+{
+ *:* :: ""
+} 1
+
+restricted:
CCTK_INT ML_BSSN_MP_RHS_calc_every "ML_BSSN_MP_RHS_calc_every"
{
*:* :: ""
@@ -284,6 +303,12 @@ CCTK_INT ML_BSSN_MP_convertFromADMBaseGamma_calc_offset "ML_BSSN_MP_convertFromA
} 0
restricted:
+CCTK_INT ML_BSSN_MP_setBetaDriver_calc_offset "ML_BSSN_MP_setBetaDriver_calc_offset"
+{
+ *:* :: ""
+} 0
+
+restricted:
CCTK_INT ML_BSSN_MP_RHS_calc_offset "ML_BSSN_MP_RHS_calc_offset"
{
*:* :: ""
diff --git a/ML_BSSN_MP/schedule.ccl b/ML_BSSN_MP/schedule.ccl
index 0624678..737f61f 100644
--- a/ML_BSSN_MP/schedule.ccl
+++ b/ML_BSSN_MP/schedule.ccl
@@ -1,6 +1,8 @@
# File produced by Kranc
+STORAGE: ML_BetaDriver[1]
+
STORAGE: ML_cons_detg[1]
STORAGE: ML_cons_Gamma[1]
@@ -365,6 +367,15 @@ if (CCTK_EQUALS(my_initial_data, "ADMBase"))
} "ML_BSSN_MP_convertFromADMBaseGamma"
}
+
+if (CCTK_EQUALS(UseSpatialBetaDriver, "yes"))
+{
+ schedule ML_BSSN_MP_setBetaDriver AT initial AFTER ADMBase_PostInitial AFTER ML_BSSN_MP_convertFromADMBase
+ {
+ LANG: C
+ } "ML_BSSN_MP_setBetaDriver"
+}
+
schedule ML_BSSN_MP_RHS IN ML_BSSN_MP_evolCalcGroup
{
LANG: C
diff --git a/ML_BSSN_MP/src/ML_BSSN_MP_Minkowski.c b/ML_BSSN_MP/src/ML_BSSN_MP_Minkowski.c
index 277ef96..8874902 100644
--- a/ML_BSSN_MP/src/ML_BSSN_MP_Minkowski.c
+++ b/ML_BSSN_MP/src/ML_BSSN_MP_Minkowski.c
@@ -108,6 +108,7 @@ void ML_BSSN_MP_Minkowski_Body(cGH const * const cctkGH, CCTK_INT const dir, CCT
CCTK_REAL At11L = INITVALUE, At12L = INITVALUE, At13L = INITVALUE, At22L = INITVALUE, At23L = INITVALUE, At33L = INITVALUE;
CCTK_REAL B1L = INITVALUE, B2L = INITVALUE, B3L = INITVALUE;
CCTK_REAL beta1L = INITVALUE, beta2L = INITVALUE, beta3L = INITVALUE;
+ CCTK_REAL etaL = INITVALUE;
CCTK_REAL gt11L = INITVALUE, gt12L = INITVALUE, gt13L = INITVALUE, gt22L = INITVALUE, gt23L = INITVALUE, gt33L = INITVALUE;
CCTK_REAL phiL = INITVALUE;
CCTK_REAL trKL = INITVALUE;
@@ -177,6 +178,8 @@ void ML_BSSN_MP_Minkowski_Body(cGH const * const cctkGH, CCTK_INT const dir, CCT
B3L = 0;
+ etaL = BetaDriver;
+
/* Copy local copies back to grid functions */
A[index] = AL;
@@ -193,6 +196,7 @@ void ML_BSSN_MP_Minkowski_Body(cGH const * const cctkGH, CCTK_INT const dir, CCT
beta1[index] = beta1L;
beta2[index] = beta2L;
beta3[index] = beta3L;
+ eta[index] = etaL;
gt11[index] = gt11L;
gt12[index] = gt12L;
gt13[index] = gt13L;
diff --git a/ML_BSSN_MP/src/ML_BSSN_MP_RHS.c b/ML_BSSN_MP/src/ML_BSSN_MP_RHS.c
index 9abbddd..0142ad7 100644
--- a/ML_BSSN_MP/src/ML_BSSN_MP_RHS.c
+++ b/ML_BSSN_MP/src/ML_BSSN_MP_RHS.c
@@ -146,6 +146,7 @@ void ML_BSSN_MP_RHS_Body(cGH const * const cctkGH, CCTK_INT const dir, CCTK_INT
CCTK_REAL dJ111L = INITVALUE, dJ112L = INITVALUE, dJ113L = INITVALUE, dJ122L = INITVALUE, dJ123L = INITVALUE, dJ133L = INITVALUE;
CCTK_REAL dJ211L = INITVALUE, dJ212L = INITVALUE, dJ213L = INITVALUE, dJ222L = INITVALUE, dJ223L = INITVALUE, dJ233L = INITVALUE;
CCTK_REAL dJ311L = INITVALUE, dJ312L = INITVALUE, dJ313L = INITVALUE, dJ322L = INITVALUE, dJ323L = INITVALUE, dJ333L = INITVALUE;
+ CCTK_REAL etaL = INITVALUE;
CCTK_REAL gt11L = INITVALUE, gt11rhsL = INITVALUE, gt12L = INITVALUE, gt12rhsL = INITVALUE, gt13L = INITVALUE, gt13rhsL = INITVALUE;
CCTK_REAL gt22L = INITVALUE, gt22rhsL = INITVALUE, gt23L = INITVALUE, gt23rhsL = INITVALUE, gt33L = INITVALUE, gt33rhsL = INITVALUE;
CCTK_REAL J11L = INITVALUE, J12L = INITVALUE, J13L = INITVALUE, J21L = INITVALUE, J22L = INITVALUE, J23L = INITVALUE;
@@ -301,6 +302,7 @@ void ML_BSSN_MP_RHS_Body(cGH const * const cctkGH, CCTK_INT const dir, CCTK_INT
dJ322L = dJ322[index];
dJ323L = dJ323[index];
dJ333L = dJ333[index];
+ etaL = eta[index];
gt11L = gt11[index];
gt12L = gt12[index];
gt13L = gt13[index];
@@ -1584,7 +1586,7 @@ void ML_BSSN_MP_RHS_Body(cGH const * const cctkGH, CCTK_INT const dir, CCTK_INT
PDupwindNth3(beta3, i, j, k)*(beta1L*J31L + beta2L*J32L + beta3L*J33L))*ShiftAdvectionCoeff +
B3L*ShiftGammaCoeff;
- B1rhsL = -(B1L*BetaDriver) + (beta1L*((PDupwindNth1(B1, i, j, k) - PDupwindNth1(Xt1, i, j, k))*J11L +
+ B1rhsL = -(B1L*etaL) + (beta1L*((PDupwindNth1(B1, i, j, k) - PDupwindNth1(Xt1, i, j, k))*J11L +
(PDupwindNth2(B1, i, j, k) - PDupwindNth2(Xt1, i, j, k))*J21L +
(PDupwindNth3(B1, i, j, k) - PDupwindNth3(Xt1, i, j, k))*J31L) +
beta2L*((PDupwindNth1(B1, i, j, k) - PDupwindNth1(Xt1, i, j, k))*J12L +
@@ -1594,7 +1596,7 @@ void ML_BSSN_MP_RHS_Body(cGH const * const cctkGH, CCTK_INT const dir, CCTK_INT
(PDupwindNth2(B1, i, j, k) - PDupwindNth2(Xt1, i, j, k))*J23L +
(PDupwindNth3(B1, i, j, k) - PDupwindNth3(Xt1, i, j, k))*J33L))*ShiftAdvectionCoeff + Xt1rhsL;
- B2rhsL = -(B2L*BetaDriver) + (beta1L*((PDupwindNth1(B2, i, j, k) - PDupwindNth1(Xt2, i, j, k))*J11L +
+ B2rhsL = -(B2L*etaL) + (beta1L*((PDupwindNth1(B2, i, j, k) - PDupwindNth1(Xt2, i, j, k))*J11L +
(PDupwindNth2(B2, i, j, k) - PDupwindNth2(Xt2, i, j, k))*J21L +
(PDupwindNth3(B2, i, j, k) - PDupwindNth3(Xt2, i, j, k))*J31L) +
beta2L*((PDupwindNth1(B2, i, j, k) - PDupwindNth1(Xt2, i, j, k))*J12L +
@@ -1604,7 +1606,7 @@ void ML_BSSN_MP_RHS_Body(cGH const * const cctkGH, CCTK_INT const dir, CCTK_INT
(PDupwindNth2(B2, i, j, k) - PDupwindNth2(Xt2, i, j, k))*J23L +
(PDupwindNth3(B2, i, j, k) - PDupwindNth3(Xt2, i, j, k))*J33L))*ShiftAdvectionCoeff + Xt2rhsL;
- B3rhsL = -(B3L*BetaDriver) + (beta1L*((PDupwindNth1(B3, i, j, k) - PDupwindNth1(Xt3, i, j, k))*J11L +
+ B3rhsL = -(B3L*etaL) + (beta1L*((PDupwindNth1(B3, i, j, k) - PDupwindNth1(Xt3, i, j, k))*J11L +
(PDupwindNth2(B3, i, j, k) - PDupwindNth2(Xt3, i, j, k))*J21L +
(PDupwindNth3(B3, i, j, k) - PDupwindNth3(Xt3, i, j, k))*J31L) +
beta2L*((PDupwindNth1(B3, i, j, k) - PDupwindNth1(Xt3, i, j, k))*J12L +
diff --git a/ML_BSSN_MP/src/ML_BSSN_MP_convertFromADMBase.c b/ML_BSSN_MP/src/ML_BSSN_MP_convertFromADMBase.c
index 0ed9bfd..d64668c 100644
--- a/ML_BSSN_MP/src/ML_BSSN_MP_convertFromADMBase.c
+++ b/ML_BSSN_MP/src/ML_BSSN_MP_convertFromADMBase.c
@@ -115,6 +115,7 @@ void ML_BSSN_MP_convertFromADMBase_Body(cGH const * const cctkGH, CCTK_INT const
CCTK_REAL betaxL = INITVALUE;
CCTK_REAL betayL = INITVALUE;
CCTK_REAL betazL = INITVALUE;
+ CCTK_REAL etaL = INITVALUE;
CCTK_REAL gt11L = INITVALUE, gt12L = INITVALUE, gt13L = INITVALUE, gt22L = INITVALUE, gt23L = INITVALUE, gt33L = INITVALUE;
CCTK_REAL gxxL = INITVALUE;
CCTK_REAL gxyL = INITVALUE;
@@ -239,6 +240,8 @@ void ML_BSSN_MP_convertFromADMBase_Body(cGH const * const cctkGH, CCTK_INT const
beta3L = betazL;
+ etaL = BetaDriver;
+
/* Copy local copies back to grid functions */
alpha[index] = alphaL;
@@ -251,6 +254,7 @@ void ML_BSSN_MP_convertFromADMBase_Body(cGH const * const cctkGH, CCTK_INT const
beta1[index] = beta1L;
beta2[index] = beta2L;
beta3[index] = beta3L;
+ eta[index] = etaL;
gt11[index] = gt11L;
gt12[index] = gt12L;
gt13[index] = gt13L;
diff --git a/ML_BSSN_MP/src/ML_BSSN_MP_setBetaDriver.c b/ML_BSSN_MP/src/ML_BSSN_MP_setBetaDriver.c
new file mode 100644
index 0000000..bc24ee7
--- /dev/null
+++ b/ML_BSSN_MP/src/ML_BSSN_MP_setBetaDriver.c
@@ -0,0 +1,142 @@
+/* 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_MP_setBetaDriver_Body(cGH const * const cctkGH, CCTK_INT const dir, CCTK_INT const face, CCTK_REAL const normal[3], CCTK_REAL const tangentA[3], CCTK_REAL const tangentB[3], CCTK_INT const min[3], CCTK_INT const max[3], CCTK_INT const n_subblock_gfs, CCTK_REAL * const subblock_gfs[])
+{
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+
+ /* Declare finite differencing variables */
+ CCTK_REAL dx = INITVALUE, dy = INITVALUE, dz = INITVALUE;
+ CCTK_REAL dxi = INITVALUE, dyi = INITVALUE, dzi = INITVALUE;
+ CCTK_REAL khalf = INITVALUE, kthird = INITVALUE, ktwothird = INITVALUE, kfourthird = INITVALUE, keightthird = INITVALUE;
+ CCTK_REAL hdxi = INITVALUE, hdyi = INITVALUE, hdzi = INITVALUE;
+
+
+ /* Declare predefined quantities */
+ CCTK_REAL p1o12dx = INITVALUE;
+ CCTK_REAL p1o12dy = INITVALUE;
+ CCTK_REAL p1o12dz = INITVALUE;
+ CCTK_REAL p1o144dxdy = INITVALUE;
+ CCTK_REAL p1o144dxdz = INITVALUE;
+ CCTK_REAL p1o144dydz = INITVALUE;
+ CCTK_REAL p1odx = INITVALUE;
+ CCTK_REAL p1ody = INITVALUE;
+ CCTK_REAL p1odz = INITVALUE;
+ CCTK_REAL pm1o12dx2 = INITVALUE;
+ CCTK_REAL pm1o12dy2 = INITVALUE;
+ CCTK_REAL pm1o12dz2 = INITVALUE;
+
+ if (verbose > 1)
+ {
+ CCTK_VInfo(CCTK_THORNSTRING,"Entering ML_BSSN_MP_setBetaDriver_Body");
+ }
+
+ if (cctk_iteration % ML_BSSN_MP_setBetaDriver_calc_every != ML_BSSN_MP_setBetaDriver_calc_offset)
+ {
+ return;
+ }
+
+ /* Include user-supplied include files */
+
+ /* Initialise finite differencing variables */
+ dx = CCTK_DELTA_SPACE(0);
+ dy = CCTK_DELTA_SPACE(1);
+ dz = CCTK_DELTA_SPACE(2);
+ dxi = 1.0 / dx;
+ dyi = 1.0 / dy;
+ dzi = 1.0 / dz;
+ khalf = 0.5;
+ kthird = 1/3.0;
+ ktwothird = 2.0/3.0;
+ kfourthird = 4.0/3.0;
+ keightthird = 8.0/3.0;
+ hdxi = 0.5 * dxi;
+ hdyi = 0.5 * dyi;
+ hdzi = 0.5 * dzi;
+
+ /* Initialize predefined quantities */
+ p1o12dx = INV(dx)/12.;
+ p1o12dy = INV(dy)/12.;
+ p1o12dz = INV(dz)/12.;
+ p1o144dxdy = (INV(dx)*INV(dy))/144.;
+ p1o144dxdz = (INV(dx)*INV(dz))/144.;
+ p1o144dydz = (INV(dy)*INV(dz))/144.;
+ p1odx = INV(dx);
+ p1ody = INV(dy);
+ p1odz = INV(dz);
+ pm1o12dx2 = -pow(dx,-2)/12.;
+ pm1o12dy2 = -pow(dy,-2)/12.;
+ pm1o12dz2 = -pow(dz,-2)/12.;
+
+ /* Loop over the grid points */
+ #pragma omp parallel
+ LC_LOOP3 (ML_BSSN_MP_setBetaDriver,
+ 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 subblock_index = INITVALUE;
+ index = CCTK_GFINDEX3D(cctkGH,i,j,k);
+ subblock_index = i - min[0] + (max[0] - min[0]) * (j - min[1] + (max[1]-min[1]) * (k - min[2]));
+
+ /* Declare shorthands */
+
+ /* Declare local copies of grid functions */
+ CCTK_REAL etaL = INITVALUE;
+ CCTK_REAL rL = INITVALUE;
+ /* Declare precomputed derivatives*/
+
+ /* Declare derivatives */
+
+ /* Assign local copies of grid functions */
+ etaL = eta[index];
+ rL = r[index];
+
+ /* Assign local copies of subblock grid functions */
+
+ /* Include user supplied include files */
+
+ /* Precompute derivatives (new style) */
+
+ /* Precompute derivatives (old style) */
+
+ /* Calculate temporaries and grid functions */
+ etaL = etaL*IfThen(rL > SpatialBetaDriverRadius,SpatialBetaDriverRadius*INV(rL),1);
+
+
+ /* Copy local copies back to grid functions */
+ eta[index] = etaL;
+
+ /* Copy local copies back to subblock grid functions */
+ }
+ LC_ENDLOOP3 (ML_BSSN_MP_setBetaDriver);
+}
+
+void ML_BSSN_MP_setBetaDriver(CCTK_ARGUMENTS)
+{
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+ GenericFD_LoopOverEverything(cctkGH, &ML_BSSN_MP_setBetaDriver_Body);
+}
diff --git a/ML_BSSN_MP/src/RegisterSymmetries.c b/ML_BSSN_MP/src/RegisterSymmetries.c
index e2081a6..d6d178b 100644
--- a/ML_BSSN_MP/src/RegisterSymmetries.c
+++ b/ML_BSSN_MP/src/RegisterSymmetries.c
@@ -144,6 +144,11 @@ void ML_BSSN_MP_RegisterSymmetries(CCTK_ARGUMENTS)
sym[0] = 1;
sym[1] = 1;
sym[2] = 1;
+ SetCartSymVN(cctkGH, sym, "ML_BSSN_MP::eta");
+
+ sym[0] = 1;
+ sym[1] = 1;
+ sym[2] = 1;
SetCartSymVN(cctkGH, sym, "ML_BSSN_MP::cS");
sym[0] = -1;
diff --git a/ML_BSSN_MP/src/make.code.defn b/ML_BSSN_MP/src/make.code.defn
index 8d63814..ed62c1a 100644
--- a/ML_BSSN_MP/src/make.code.defn
+++ b/ML_BSSN_MP/src/make.code.defn
@@ -1,3 +1,3 @@
# File produced by Kranc
-SRCS = Startup.c RegisterMoL.c RegisterSymmetries.c ML_BSSN_MP_Minkowski.c ML_BSSN_MP_convertFromADMBase.c ML_BSSN_MP_convertFromADMBaseGamma.c ML_BSSN_MP_RHS.c ML_BSSN_MP_RHSStaticBoundary.c ML_BSSN_MP_RHSRadiativeBoundary.c ML_BSSN_MP_enforce.c ML_BSSN_MP_boundary.c ML_BSSN_MP_convertToADMBase.c ML_BSSN_MP_convertToADMBaseDtLapseShift.c ML_BSSN_MP_convertToADMBaseDtLapseShiftBoundary.c ML_BSSN_MP_convertToADMBaseFakeDtLapseShift.c ML_BSSN_MP_constraints.c ML_BSSN_MP_constraints_boundary.c Boundaries.c
+SRCS = Startup.c RegisterMoL.c RegisterSymmetries.c ML_BSSN_MP_Minkowski.c ML_BSSN_MP_convertFromADMBase.c ML_BSSN_MP_convertFromADMBaseGamma.c ML_BSSN_MP_setBetaDriver.c ML_BSSN_MP_RHS.c ML_BSSN_MP_RHSStaticBoundary.c ML_BSSN_MP_RHSRadiativeBoundary.c ML_BSSN_MP_enforce.c ML_BSSN_MP_boundary.c ML_BSSN_MP_convertToADMBase.c ML_BSSN_MP_convertToADMBaseDtLapseShift.c ML_BSSN_MP_convertToADMBaseDtLapseShiftBoundary.c ML_BSSN_MP_convertToADMBaseFakeDtLapseShift.c ML_BSSN_MP_constraints.c ML_BSSN_MP_constraints_boundary.c Boundaries.c