aboutsummaryrefslogtreecommitdiff
path: root/ML_BSSN
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@cct.lsu.edu>2011-06-12 22:22:13 -0400
committerErik Schnetter <schnetter@cct.lsu.edu>2011-06-12 22:22:13 -0400
commita8b69dc06e565427bc96d16c443bb4b3392c4365 (patch)
tree0bbd7ff18f442eacee9454384bd5548315900d35 /ML_BSSN
parent05347a08d0c9bd2a87846ab4ad8990fe26274a4a (diff)
Rename tensor indices
Tensor indices with numbers are not supported any more in Kranc; rename them to use letters instead.
Diffstat (limited to 'ML_BSSN')
-rw-r--r--ML_BSSN/src/ML_BSSN_Minkowski.cc2
-rw-r--r--ML_BSSN/src/ML_BSSN_RHS1.cc141
-rw-r--r--ML_BSSN/src/ML_BSSN_RHS2.cc62
-rw-r--r--ML_BSSN/src/ML_BSSN_boundary.cc2
-rw-r--r--ML_BSSN/src/ML_BSSN_constraints1.cc62
-rw-r--r--ML_BSSN/src/ML_BSSN_constraints2.cc2
-rw-r--r--ML_BSSN/src/ML_BSSN_convertFromADMBaseGamma.cc59
-rw-r--r--ML_BSSN/src/ML_BSSN_convertToADMBase.cc3
8 files changed, 178 insertions, 155 deletions
diff --git a/ML_BSSN/src/ML_BSSN_Minkowski.cc b/ML_BSSN/src/ML_BSSN_Minkowski.cc
index 10d5506..77b7e2d 100644
--- a/ML_BSSN/src/ML_BSSN_Minkowski.cc
+++ b/ML_BSSN/src/ML_BSSN_Minkowski.cc
@@ -105,7 +105,7 @@ static void ML_BSSN_Minkowski_Body(cGH const * restrict const cctkGH, int const
/* Precompute derivatives */
/* Calculate temporaries and grid functions */
- CCTK_REAL phiL = IfThen(ToReal(conformalMethod),1,0);
+ CCTK_REAL phiL = IfThen(conformalMethod,1,0);
CCTK_REAL gt11L = 1;
diff --git a/ML_BSSN/src/ML_BSSN_RHS1.cc b/ML_BSSN/src/ML_BSSN_RHS1.cc
index f5456e3..0d11a97 100644
--- a/ML_BSSN/src/ML_BSSN_RHS1.cc
+++ b/ML_BSSN/src/ML_BSSN_RHS1.cc
@@ -534,7 +534,7 @@ static void ML_BSSN_RHS1_Body(cGH const * restrict const cctkGH, int const dir,
CCTK_REAL Xtn3 = Gt311*gtu11 + Gt322*gtu22 + 2*(Gt312*gtu12 +
Gt313*gtu13 + Gt323*gtu23) + Gt333*gtu33;
- CCTK_REAL fac1 = IfThen(ToReal(conformalMethod),-0.5*INV(phiL),1);
+ CCTK_REAL fac1 = IfThen(conformalMethod,-0.5*INV(phiL),1);
CCTK_REAL cdphi1 = fac1*PDstandardNth1phi;
@@ -572,8 +572,7 @@ static void ML_BSSN_RHS1_Body(cGH const * restrict const cctkGH, int const dir,
CCTK_REAL Atu33 = Atm31*gtu13 + Atm32*gtu23 + Atm33*gtu33;
- CCTK_REAL e4phi =
- IfThen(ToReal(conformalMethod),INV(SQR(phiL)),exp(4*phiL));
+ CCTK_REAL e4phi = IfThen(conformalMethod,INV(SQR(phiL)),exp(4*phiL));
CCTK_REAL em4phi = INV(e4phi);
@@ -599,11 +598,14 @@ static void ML_BSSN_RHS1_Body(cGH const * restrict const cctkGH, int const dir,
beta1L*PDupwindNthAnti1phi + beta2L*PDupwindNthAnti2phi +
beta3L*PDupwindNthAnti3phi + PDupwindNthSymm1phi*Abs(beta1L) +
PDupwindNthSymm2phi*Abs(beta2L) + PDupwindNthSymm3phi*Abs(beta3L) +
- IfThen(ToReal(conformalMethod),phiL*(-0.333333333333333333333333333333*(PDstandardNth1beta1
- + PDstandardNth2beta2 + PDstandardNth3beta3) +
- 0.333333333333333333333333333333*alphaL*trKL),0.166666666666666666666666666667*(PDstandardNth1beta1
- + PDstandardNth2beta2 + PDstandardNth3beta3) -
- 0.166666666666666666666666666667*alphaL*trKL);
+ IfThen(conformalMethod,-0.333333333333333333333333333333*(PDstandardNth1beta1
+ +
+ PDstandardNth2beta2)*phiL,0.166666666666666666666666666667*(PDstandardNth1beta1
+ + PDstandardNth2beta2)) +
+ IfThen(conformalMethod,phiL*(-0.333333333333333333333333333333*PDstandardNth3beta3
+ +
+ 0.333333333333333333333333333333*alphaL*trKL),0.166666666666666666666666666667*PDstandardNth3beta3
+ - 0.166666666666666666666666666667*alphaL*trKL);
CCTK_REAL gt11rhsL = -2*alphaL*At11L + epsdiss1*PDdissipationNth1gt11
+ epsdiss2*PDdissipationNth2gt11 + epsdiss3*PDdissipationNth3gt11 +
@@ -670,72 +672,77 @@ static void ML_BSSN_RHS1_Body(cGH const * restrict const cctkGH, int const dir,
CCTK_REAL dotXt1 =
0.333333333333333333333333333333*(7*(gtu12*PDstandardNth12beta1 +
- gtu13*PDstandardNth13beta1) + gtu11*(4*PDstandardNth11beta1 +
- PDstandardNth12beta2 + PDstandardNth13beta3) +
- gtu12*(PDstandardNth22beta2 + PDstandardNth23beta3) +
- gtu13*(PDstandardNth23beta2 + PDstandardNth33beta3) -
+ gtu13*PDstandardNth13beta1) + 6*gtu23*PDstandardNth23beta1 -
6*(Atu11*PDstandardNth1alpha + Atu12*PDstandardNth2alpha +
- Atu13*PDstandardNth3alpha) + 6*(gtu23*PDstandardNth23beta1 +
- alphaL*(6*(Atu11*cdphi1 + Atu12*cdphi2 + Atu13*cdphi3) + Atu11*Gt111 +
- Atu22*Gt122 + 2*(Atu12*Gt112 + Atu13*Gt113 + Atu23*Gt123) + Atu33*Gt133
- - 0.666666666666666666666666666667*(gtu11*PDstandardNth1trK +
- gtu12*PDstandardNth2trK + gtu13*PDstandardNth3trK))) -
- 150.7964473723100754462068823974161384415*alphaL*(gtu11*S1 + gtu12*S2 +
- gtu13*S3) + (-3*PDstandardNth1beta1 + 2*(PDstandardNth1beta1 +
- PDstandardNth2beta2 + PDstandardNth3beta3))*Xtn1 -
- 3*(PDstandardNth2beta1*Xtn2 + PDstandardNth3beta1*Xtn3) +
- 3*(epsdiss1*PDdissipationNth1Xt1 + epsdiss2*PDdissipationNth2Xt1 +
- epsdiss3*PDdissipationNth3Xt1 + gtu22*PDstandardNth22beta1 +
- gtu33*PDstandardNth33beta1 + beta1L*PDupwindNthAnti1Xt1 +
- beta2L*PDupwindNthAnti2Xt1 + beta3L*PDupwindNthAnti3Xt1 +
- PDupwindNthSymm1Xt1*Abs(beta1L) + PDupwindNthSymm2Xt1*Abs(beta2L) +
- PDupwindNthSymm3Xt1*Abs(beta3L)));
+ Atu13*PDstandardNth3alpha) + gtu11*(4*PDstandardNth11beta1 +
+ PDstandardNth12beta2 + PDstandardNth13beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S1) +
+ gtu12*(PDstandardNth22beta2 + PDstandardNth23beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S2) +
+ gtu13*(PDstandardNth23beta2 + PDstandardNth33beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S3) +
+ (-PDstandardNth1beta1 + 2*PDstandardNth3beta3)*Xtn1 +
+ 2*(alphaL*(18*(Atu11*cdphi1 + Atu12*cdphi2 + Atu13*cdphi3) +
+ 6*(Atu12*Gt112 + Atu13*Gt113 + Atu23*Gt123) + 3*(Atu11*Gt111 +
+ Atu22*Gt122 + Atu33*Gt133) - 2*(gtu11*PDstandardNth1trK +
+ gtu12*PDstandardNth2trK + gtu13*PDstandardNth3trK)) +
+ PDstandardNth2beta2*Xtn1) - 3*(PDstandardNth2beta1*Xtn2 +
+ PDstandardNth3beta1*Xtn3) + 3*(epsdiss1*PDdissipationNth1Xt1 +
+ epsdiss2*PDdissipationNth2Xt1 + epsdiss3*PDdissipationNth3Xt1 +
+ gtu22*PDstandardNth22beta1 + gtu33*PDstandardNth33beta1 +
+ beta1L*PDupwindNthAnti1Xt1 + beta2L*PDupwindNthAnti2Xt1 +
+ beta3L*PDupwindNthAnti3Xt1 + PDupwindNthSymm1Xt1*Abs(beta1L) +
+ PDupwindNthSymm2Xt1*Abs(beta2L) + PDupwindNthSymm3Xt1*Abs(beta3L)));
CCTK_REAL dotXt2 =
- 0.333333333333333333333333333333*(gtu12*(PDstandardNth11beta1 +
- 7*PDstandardNth12beta2 + PDstandardNth13beta3) +
- gtu22*(PDstandardNth12beta1 + 4*PDstandardNth22beta2 +
- PDstandardNth23beta3) + gtu23*(PDstandardNth13beta1 +
- 7*PDstandardNth23beta2 + PDstandardNth33beta3) -
+ 0.333333333333333333333333333333*(6*gtu13*PDstandardNth13beta2 -
6*(Atu12*PDstandardNth1alpha + Atu22*PDstandardNth2alpha +
- Atu23*PDstandardNth3alpha) + 6*(gtu13*PDstandardNth13beta2 +
- alphaL*(6*(Atu12*cdphi1 + Atu22*cdphi2 + Atu23*cdphi3) + Atu11*Gt211 +
- Atu22*Gt222 + 2*(Atu12*Gt212 + Atu13*Gt213 + Atu23*Gt223) + Atu33*Gt233
- - 0.666666666666666666666666666667*(gtu12*PDstandardNth1trK +
- gtu22*PDstandardNth2trK + gtu23*PDstandardNth3trK))) -
- 150.7964473723100754462068823974161384415*alphaL*(gtu12*S1 + gtu22*S2 +
- gtu23*S3) + 2*(PDstandardNth1beta1 + PDstandardNth2beta2 +
- PDstandardNth3beta3)*Xtn2 - 3*(PDstandardNth1beta2*Xtn1 +
- PDstandardNth2beta2*Xtn2 + PDstandardNth3beta2*Xtn3) +
- 3*(epsdiss1*PDdissipationNth1Xt2 + epsdiss2*PDdissipationNth2Xt2 +
- epsdiss3*PDdissipationNth3Xt2 + gtu11*PDstandardNth11beta2 +
- gtu33*PDstandardNth33beta2 + beta1L*PDupwindNthAnti1Xt2 +
- beta2L*PDupwindNthAnti2Xt2 + beta3L*PDupwindNthAnti3Xt2 +
- PDupwindNthSymm1Xt2*Abs(beta1L) + PDupwindNthSymm2Xt2*Abs(beta2L) +
- PDupwindNthSymm3Xt2*Abs(beta3L)));
+ Atu23*PDstandardNth3alpha) + gtu12*(PDstandardNth11beta1 +
+ 7*PDstandardNth12beta2 + PDstandardNth13beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S1) +
+ gtu22*(PDstandardNth12beta1 + 4*PDstandardNth22beta2 +
+ PDstandardNth23beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S2) +
+ gtu23*(PDstandardNth13beta1 + 7*PDstandardNth23beta2 +
+ PDstandardNth33beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S3) +
+ (-PDstandardNth2beta2 + 2*PDstandardNth3beta3)*Xtn2 +
+ 2*(alphaL*(18*(Atu12*cdphi1 + Atu22*cdphi2 + Atu23*cdphi3) +
+ 6*(Atu12*Gt212 + Atu13*Gt213 + Atu23*Gt223) + 3*(Atu11*Gt211 +
+ Atu22*Gt222 + Atu33*Gt233) - 2*(gtu12*PDstandardNth1trK +
+ gtu22*PDstandardNth2trK + gtu23*PDstandardNth3trK)) +
+ PDstandardNth1beta1*Xtn2) - 3*(PDstandardNth1beta2*Xtn1 +
+ PDstandardNth3beta2*Xtn3) + 3*(epsdiss1*PDdissipationNth1Xt2 +
+ epsdiss2*PDdissipationNth2Xt2 + epsdiss3*PDdissipationNth3Xt2 +
+ gtu11*PDstandardNth11beta2 + gtu33*PDstandardNth33beta2 +
+ beta1L*PDupwindNthAnti1Xt2 + beta2L*PDupwindNthAnti2Xt2 +
+ beta3L*PDupwindNthAnti3Xt2 + PDupwindNthSymm1Xt2*Abs(beta1L) +
+ PDupwindNthSymm2Xt2*Abs(beta2L) + PDupwindNthSymm3Xt2*Abs(beta3L)));
CCTK_REAL dotXt3 =
- 0.333333333333333333333333333333*(gtu13*(PDstandardNth11beta1 +
- PDstandardNth12beta2 + 7*PDstandardNth13beta3) +
- gtu23*(PDstandardNth12beta1 + PDstandardNth22beta2 +
- 7*PDstandardNth23beta3) + gtu33*(PDstandardNth13beta1 +
- PDstandardNth23beta2 + 4*PDstandardNth33beta3) -
+ 0.333333333333333333333333333333*(6*gtu12*PDstandardNth12beta3 -
6*(Atu13*PDstandardNth1alpha + Atu23*PDstandardNth2alpha +
- Atu33*PDstandardNth3alpha) + 6*(gtu12*PDstandardNth12beta3 +
- alphaL*(6*(Atu13*cdphi1 + Atu23*cdphi2 + Atu33*cdphi3) + Atu11*Gt311 +
- Atu22*Gt322 + 2*(Atu12*Gt312 + Atu13*Gt313 + Atu23*Gt323) + Atu33*Gt333
- - 0.666666666666666666666666666667*(gtu13*PDstandardNth1trK +
- gtu23*PDstandardNth2trK + gtu33*PDstandardNth3trK))) -
- 150.7964473723100754462068823974161384415*alphaL*(gtu13*S1 + gtu23*S2 +
- gtu33*S3) + 2*(PDstandardNth1beta1 + PDstandardNth2beta2 +
- PDstandardNth3beta3)*Xtn3 - 3*(PDstandardNth1beta3*Xtn1 +
- PDstandardNth2beta3*Xtn2 + PDstandardNth3beta3*Xtn3) +
- 3*(epsdiss1*PDdissipationNth1Xt3 + epsdiss2*PDdissipationNth2Xt3 +
- epsdiss3*PDdissipationNth3Xt3 + gtu11*PDstandardNth11beta3 +
- gtu22*PDstandardNth22beta3 + beta1L*PDupwindNthAnti1Xt3 +
- beta2L*PDupwindNthAnti2Xt3 + beta3L*PDupwindNthAnti3Xt3 +
- PDupwindNthSymm1Xt3*Abs(beta1L) + PDupwindNthSymm2Xt3*Abs(beta2L) +
- PDupwindNthSymm3Xt3*Abs(beta3L)));
+ Atu33*PDstandardNth3alpha) + gtu13*(PDstandardNth11beta1 +
+ PDstandardNth12beta2 + 7*PDstandardNth13beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S1) +
+ gtu23*(PDstandardNth12beta1 + PDstandardNth22beta2 +
+ 7*PDstandardNth23beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S2) +
+ gtu33*(PDstandardNth13beta1 + PDstandardNth23beta2 +
+ 4*PDstandardNth33beta3 -
+ 150.7964473723100754462068823974161384415*alphaL*S3) -
+ 3*(PDstandardNth1beta3*Xtn1 + PDstandardNth2beta3*Xtn2) +
+ (2*PDstandardNth2beta2 - PDstandardNth3beta3)*Xtn3 +
+ 2*(alphaL*(18*(Atu13*cdphi1 + Atu23*cdphi2 + Atu33*cdphi3) +
+ 6*(Atu12*Gt312 + Atu13*Gt313 + Atu23*Gt323) + 3*(Atu11*Gt311 +
+ Atu22*Gt322 + Atu33*Gt333) - 2*(gtu13*PDstandardNth1trK +
+ gtu23*PDstandardNth2trK + gtu33*PDstandardNth3trK)) +
+ PDstandardNth1beta1*Xtn3) + 3*(epsdiss1*PDdissipationNth1Xt3 +
+ epsdiss2*PDdissipationNth2Xt3 + epsdiss3*PDdissipationNth3Xt3 +
+ gtu11*PDstandardNth11beta3 + gtu22*PDstandardNth22beta3 +
+ beta1L*PDupwindNthAnti1Xt3 + beta2L*PDupwindNthAnti2Xt3 +
+ beta3L*PDupwindNthAnti3Xt3 + PDupwindNthSymm1Xt3*Abs(beta1L) +
+ PDupwindNthSymm2Xt3*Abs(beta2L) + PDupwindNthSymm3Xt3*Abs(beta3L)));
CCTK_REAL Xt1rhsL = dotXt1;
diff --git a/ML_BSSN/src/ML_BSSN_RHS2.cc b/ML_BSSN/src/ML_BSSN_RHS2.cc
index 14e28ba..11b0ec3 100644
--- a/ML_BSSN/src/ML_BSSN_RHS2.cc
+++ b/ML_BSSN/src/ML_BSSN_RHS2.cc
@@ -467,15 +467,16 @@ static void ML_BSSN_RHS2_Body(cGH const * restrict const cctkGH, int const dir,
CCTK_REAL Xtn3 = Gt311*gtu11 + Gt322*gtu22 + 2*(Gt312*gtu12 +
Gt313*gtu13 + Gt323*gtu23) + Gt333*gtu33;
- CCTK_REAL Rt11 = 3*(Gt111*Gtlu111 + Gt112*Gtlu112 + Gt113*Gtlu113) +
- 2*(Gt211*Gtlu121 + Gt212*Gtlu122 + Gt213*Gtlu123 + Gt311*Gtlu131 +
- Gt312*Gtlu132 + Gt313*Gtlu133) + Gt211*Gtlu211 + Gt212*Gtlu212 +
+ CCTK_REAL Rt11 = 0.5*(6*(Gt111*Gtlu111 + Gt112*Gtlu112 +
+ Gt113*Gtlu113) + 4*(Gt211*Gtlu121 + Gt212*Gtlu122 + Gt213*Gtlu123 +
+ Gt311*Gtlu131 + Gt312*Gtlu132 + Gt313*Gtlu133) -
+ gtu11*PDstandardNth11gt11 - 2*gtu12*PDstandardNth12gt11 -
+ 2*gtu13*PDstandardNth13gt11 + 2*(Gt211*Gtlu211 + Gt212*Gtlu212 +
Gt213*Gtlu213 + Gt311*Gtlu311 + Gt312*Gtlu312 + Gt313*Gtlu313 +
- gt11L*PDstandardNth1Xt1 + gt12L*PDstandardNth1Xt2 +
- gt13L*PDstandardNth1Xt3 + 0.5*(-(gtu11*PDstandardNth11gt11) -
- 2*gtu12*PDstandardNth12gt11 - 2*gtu13*PDstandardNth13gt11 -
- gtu22*PDstandardNth22gt11 - 2*gtu23*PDstandardNth23gt11 -
- gtu33*PDstandardNth33gt11) + Gtl111*Xtn1 + Gtl112*Xtn2 + Gtl113*Xtn3;
+ gt11L*PDstandardNth1Xt1) + 2*gt12L*PDstandardNth1Xt2 +
+ 2*gt13L*PDstandardNth1Xt3 - gtu22*PDstandardNth22gt11 -
+ 2*gtu23*PDstandardNth23gt11 - gtu33*PDstandardNth33gt11 + 2*Gtl111*Xtn1
+ + 2*Gtl112*Xtn2 + 2*Gtl113*Xtn3);
CCTK_REAL Rt12 = 0.5*(4*(Gt211*Gtlu221 + Gt212*Gtlu222 +
Gt213*Gtlu223) + 2*(Gt112*Gtlu111 + Gt122*Gtlu112 + Gt123*Gtlu113 +
@@ -507,15 +508,16 @@ static void ML_BSSN_RHS2_Body(cGH const * restrict const cctkGH, int const dir,
gt13L*PDstandardNth3Xt3 + Gtl113*Xtn1 + Gtl311*Xtn1 + Gtl123*Xtn2 +
Gtl312*Xtn2 + Gtl133*Xtn3 + Gtl313*Xtn3);
- CCTK_REAL Rt22 = Gt112*(Gtlu121 + 2*Gtlu211) + Gt122*(Gtlu122 +
- 2*Gtlu212) + Gt123*(Gtlu123 + 2*Gtlu213) + 3*(Gt212*Gtlu221 +
- Gt222*Gtlu222 + Gt223*Gtlu223) + 2*(Gt312*Gtlu231 + Gt322*Gtlu232 +
- Gt323*Gtlu233) + Gt312*Gtlu321 + Gt322*Gtlu322 + Gt323*Gtlu323 +
- gt12L*PDstandardNth2Xt1 + gt22L*PDstandardNth2Xt2 +
- gt23L*PDstandardNth2Xt3 + 0.5*(-(gtu11*PDstandardNth11gt22) -
- 2*gtu12*PDstandardNth12gt22 - 2*gtu13*PDstandardNth13gt22 -
- gtu22*PDstandardNth22gt22 - 2*gtu23*PDstandardNth23gt22 -
- gtu33*PDstandardNth33gt22) + Gtl212*Xtn1 + Gtl222*Xtn2 + Gtl223*Xtn3;
+ CCTK_REAL Rt22 = 0.5*(6*(Gt212*Gtlu221 + Gt222*Gtlu222 +
+ Gt223*Gtlu223) + 4*(Gt112*Gtlu211 + Gt122*Gtlu212 + Gt123*Gtlu213 +
+ Gt312*Gtlu231 + Gt322*Gtlu232 + Gt323*Gtlu233) -
+ gtu11*PDstandardNth11gt22 - 2*gtu12*PDstandardNth12gt22 -
+ 2*gtu13*PDstandardNth13gt22 - gtu22*PDstandardNth22gt22 -
+ 2*gtu23*PDstandardNth23gt22 + 2*(Gt112*Gtlu121 + Gt122*Gtlu122 +
+ Gt123*Gtlu123 + Gt312*Gtlu321 + Gt322*Gtlu322 + Gt323*Gtlu323 +
+ gt12L*PDstandardNth2Xt1) + 2*gt22L*PDstandardNth2Xt2 +
+ 2*gt23L*PDstandardNth2Xt3 - gtu33*PDstandardNth33gt22 + 2*Gtl212*Xtn1 +
+ 2*Gtl222*Xtn2 + 2*Gtl223*Xtn3);
CCTK_REAL Rt23 = 0.5*(2*(Gt112*Gtlu131 + Gt122*Gtlu132 + Gt123*Gtlu133
+ Gt113*Gtlu211 + Gt123*Gtlu212 + Gt133*Gtlu213 + Gt213*Gtlu221 +
@@ -532,17 +534,18 @@ static void ML_BSSN_RHS2_Body(cGH const * restrict const cctkGH, int const dir,
gt23L*PDstandardNth3Xt3 + Gtl213*Xtn1 + Gtl312*Xtn1 + Gtl223*Xtn2 +
Gtl322*Xtn2 + Gtl233*Xtn3 + Gtl323*Xtn3);
- CCTK_REAL Rt33 = Gt113*(Gtlu131 + 2*Gtlu311) + Gt123*(Gtlu132 +
- 2*Gtlu312) + Gt133*(Gtlu133 + 2*Gtlu313) + Gt213*(Gtlu231 + 2*Gtlu321)
- + Gt223*(Gtlu232 + 2*Gtlu322) + Gt233*(Gtlu233 + 2*Gtlu323) +
- 3*(Gt313*Gtlu331 + Gt323*Gtlu332 + Gt333*Gtlu333) +
- 0.5*(-(gtu11*PDstandardNth11gt33) - 2*gtu12*PDstandardNth12gt33 -
- 2*gtu13*PDstandardNth13gt33 - gtu22*PDstandardNth22gt33 -
- 2*gtu23*PDstandardNth23gt33 - gtu33*PDstandardNth33gt33) +
- gt13L*PDstandardNth3Xt1 + gt23L*PDstandardNth3Xt2 +
- gt33L*PDstandardNth3Xt3 + Gtl313*Xtn1 + Gtl323*Xtn2 + Gtl333*Xtn3;
+ CCTK_REAL Rt33 = 0.5*(4*(Gt113*Gtlu311 + Gt123*Gtlu312 + Gt133*Gtlu313
+ + Gt213*Gtlu321 + Gt223*Gtlu322 + Gt233*Gtlu323) + 6*(Gt313*Gtlu331 +
+ Gt323*Gtlu332 + Gt333*Gtlu333) - gtu11*PDstandardNth11gt33 -
+ 2*gtu12*PDstandardNth12gt33 - 2*gtu13*PDstandardNth13gt33 -
+ gtu22*PDstandardNth22gt33 - 2*gtu23*PDstandardNth23gt33 -
+ gtu33*PDstandardNth33gt33 + 2*(Gt113*Gtlu131 + Gt123*Gtlu132 +
+ Gt133*Gtlu133 + Gt213*Gtlu231 + Gt223*Gtlu232 + Gt233*Gtlu233 +
+ gt13L*PDstandardNth3Xt1) + 2*gt23L*PDstandardNth3Xt2 +
+ 2*gt33L*PDstandardNth3Xt3 + 2*Gtl313*Xtn1 + 2*Gtl323*Xtn2 +
+ 2*Gtl333*Xtn3);
- CCTK_REAL fac1 = IfThen(ToReal(conformalMethod),-0.5*INV(phiL),1);
+ CCTK_REAL fac1 = IfThen(conformalMethod,-0.5*INV(phiL),1);
CCTK_REAL cdphi1 = fac1*PDstandardNth1phi;
@@ -550,7 +553,7 @@ static void ML_BSSN_RHS2_Body(cGH const * restrict const cctkGH, int const dir,
CCTK_REAL cdphi3 = fac1*PDstandardNth3phi;
- CCTK_REAL fac2 = IfThen(ToReal(conformalMethod),0.5*INV(SQR(phiL)),0);
+ CCTK_REAL fac2 = IfThen(conformalMethod,0.5*INV(SQR(phiL)),0);
CCTK_REAL cdphi211 = -(fac1*(-PDstandardNth11phi +
Gt111*PDstandardNth1phi + Gt211*PDstandardNth2phi +
@@ -630,8 +633,7 @@ static void ML_BSSN_RHS2_Body(cGH const * restrict const cctkGH, int const dir,
CCTK_REAL Atm33 = At13L*gtu13 + At23L*gtu23 + At33L*gtu33;
- CCTK_REAL e4phi =
- IfThen(ToReal(conformalMethod),INV(SQR(phiL)),exp(4*phiL));
+ CCTK_REAL e4phi = IfThen(conformalMethod,INV(SQR(phiL)),exp(4*phiL));
CCTK_REAL em4phi = INV(e4phi);
diff --git a/ML_BSSN/src/ML_BSSN_boundary.cc b/ML_BSSN/src/ML_BSSN_boundary.cc
index abfef01..ac1e90f 100644
--- a/ML_BSSN/src/ML_BSSN_boundary.cc
+++ b/ML_BSSN/src/ML_BSSN_boundary.cc
@@ -141,7 +141,7 @@ static void ML_BSSN_boundary_Body(cGH const * restrict const cctkGH, int const d
/* Precompute derivatives */
/* Calculate temporaries and grid functions */
- CCTK_REAL phiL = IfThen(ToReal(conformalMethod),1,0);
+ CCTK_REAL phiL = IfThen(conformalMethod,1,0);
CCTK_REAL gt11L = 1;
diff --git a/ML_BSSN/src/ML_BSSN_constraints1.cc b/ML_BSSN/src/ML_BSSN_constraints1.cc
index de541c8..77156fc 100644
--- a/ML_BSSN/src/ML_BSSN_constraints1.cc
+++ b/ML_BSSN/src/ML_BSSN_constraints1.cc
@@ -391,15 +391,16 @@ static void ML_BSSN_constraints1_Body(cGH const * restrict const cctkGH, int con
CCTK_REAL Xtn3 = Gt311*gtu11 + Gt322*gtu22 + 2*(Gt312*gtu12 +
Gt313*gtu13 + Gt323*gtu23) + Gt333*gtu33;
- CCTK_REAL Rt11 = 3*(Gt111*Gtlu111 + Gt112*Gtlu112 + Gt113*Gtlu113) +
- 2*(Gt211*Gtlu121 + Gt212*Gtlu122 + Gt213*Gtlu123 + Gt311*Gtlu131 +
- Gt312*Gtlu132 + Gt313*Gtlu133) + Gt211*Gtlu211 + Gt212*Gtlu212 +
+ CCTK_REAL Rt11 = 0.5*(6*(Gt111*Gtlu111 + Gt112*Gtlu112 +
+ Gt113*Gtlu113) + 4*(Gt211*Gtlu121 + Gt212*Gtlu122 + Gt213*Gtlu123 +
+ Gt311*Gtlu131 + Gt312*Gtlu132 + Gt313*Gtlu133) -
+ gtu11*PDstandardNth11gt11 - 2*gtu12*PDstandardNth12gt11 -
+ 2*gtu13*PDstandardNth13gt11 + 2*(Gt211*Gtlu211 + Gt212*Gtlu212 +
Gt213*Gtlu213 + Gt311*Gtlu311 + Gt312*Gtlu312 + Gt313*Gtlu313 +
- gt11L*PDstandardNth1Xt1 + gt12L*PDstandardNth1Xt2 +
- gt13L*PDstandardNth1Xt3 + 0.5*(-(gtu11*PDstandardNth11gt11) -
- 2*gtu12*PDstandardNth12gt11 - 2*gtu13*PDstandardNth13gt11 -
- gtu22*PDstandardNth22gt11 - 2*gtu23*PDstandardNth23gt11 -
- gtu33*PDstandardNth33gt11) + Gtl111*Xtn1 + Gtl112*Xtn2 + Gtl113*Xtn3;
+ gt11L*PDstandardNth1Xt1) + 2*gt12L*PDstandardNth1Xt2 +
+ 2*gt13L*PDstandardNth1Xt3 - gtu22*PDstandardNth22gt11 -
+ 2*gtu23*PDstandardNth23gt11 - gtu33*PDstandardNth33gt11 + 2*Gtl111*Xtn1
+ + 2*Gtl112*Xtn2 + 2*Gtl113*Xtn3);
CCTK_REAL Rt12 = 0.5*(4*(Gt211*Gtlu221 + Gt212*Gtlu222 +
Gt213*Gtlu223) + 2*(Gt112*Gtlu111 + Gt122*Gtlu112 + Gt123*Gtlu113 +
@@ -431,15 +432,16 @@ static void ML_BSSN_constraints1_Body(cGH const * restrict const cctkGH, int con
gt13L*PDstandardNth3Xt3 + Gtl113*Xtn1 + Gtl311*Xtn1 + Gtl123*Xtn2 +
Gtl312*Xtn2 + Gtl133*Xtn3 + Gtl313*Xtn3);
- CCTK_REAL Rt22 = Gt112*(Gtlu121 + 2*Gtlu211) + Gt122*(Gtlu122 +
- 2*Gtlu212) + Gt123*(Gtlu123 + 2*Gtlu213) + 3*(Gt212*Gtlu221 +
- Gt222*Gtlu222 + Gt223*Gtlu223) + 2*(Gt312*Gtlu231 + Gt322*Gtlu232 +
- Gt323*Gtlu233) + Gt312*Gtlu321 + Gt322*Gtlu322 + Gt323*Gtlu323 +
- gt12L*PDstandardNth2Xt1 + gt22L*PDstandardNth2Xt2 +
- gt23L*PDstandardNth2Xt3 + 0.5*(-(gtu11*PDstandardNth11gt22) -
- 2*gtu12*PDstandardNth12gt22 - 2*gtu13*PDstandardNth13gt22 -
- gtu22*PDstandardNth22gt22 - 2*gtu23*PDstandardNth23gt22 -
- gtu33*PDstandardNth33gt22) + Gtl212*Xtn1 + Gtl222*Xtn2 + Gtl223*Xtn3;
+ CCTK_REAL Rt22 = 0.5*(6*(Gt212*Gtlu221 + Gt222*Gtlu222 +
+ Gt223*Gtlu223) + 4*(Gt112*Gtlu211 + Gt122*Gtlu212 + Gt123*Gtlu213 +
+ Gt312*Gtlu231 + Gt322*Gtlu232 + Gt323*Gtlu233) -
+ gtu11*PDstandardNth11gt22 - 2*gtu12*PDstandardNth12gt22 -
+ 2*gtu13*PDstandardNth13gt22 - gtu22*PDstandardNth22gt22 -
+ 2*gtu23*PDstandardNth23gt22 + 2*(Gt112*Gtlu121 + Gt122*Gtlu122 +
+ Gt123*Gtlu123 + Gt312*Gtlu321 + Gt322*Gtlu322 + Gt323*Gtlu323 +
+ gt12L*PDstandardNth2Xt1) + 2*gt22L*PDstandardNth2Xt2 +
+ 2*gt23L*PDstandardNth2Xt3 - gtu33*PDstandardNth33gt22 + 2*Gtl212*Xtn1 +
+ 2*Gtl222*Xtn2 + 2*Gtl223*Xtn3);
CCTK_REAL Rt23 = 0.5*(2*(Gt112*Gtlu131 + Gt122*Gtlu132 + Gt123*Gtlu133
+ Gt113*Gtlu211 + Gt123*Gtlu212 + Gt133*Gtlu213 + Gt213*Gtlu221 +
@@ -456,17 +458,18 @@ static void ML_BSSN_constraints1_Body(cGH const * restrict const cctkGH, int con
gt23L*PDstandardNth3Xt3 + Gtl213*Xtn1 + Gtl312*Xtn1 + Gtl223*Xtn2 +
Gtl322*Xtn2 + Gtl233*Xtn3 + Gtl323*Xtn3);
- CCTK_REAL Rt33 = Gt113*(Gtlu131 + 2*Gtlu311) + Gt123*(Gtlu132 +
- 2*Gtlu312) + Gt133*(Gtlu133 + 2*Gtlu313) + Gt213*(Gtlu231 + 2*Gtlu321)
- + Gt223*(Gtlu232 + 2*Gtlu322) + Gt233*(Gtlu233 + 2*Gtlu323) +
- 3*(Gt313*Gtlu331 + Gt323*Gtlu332 + Gt333*Gtlu333) +
- 0.5*(-(gtu11*PDstandardNth11gt33) - 2*gtu12*PDstandardNth12gt33 -
- 2*gtu13*PDstandardNth13gt33 - gtu22*PDstandardNth22gt33 -
- 2*gtu23*PDstandardNth23gt33 - gtu33*PDstandardNth33gt33) +
- gt13L*PDstandardNth3Xt1 + gt23L*PDstandardNth3Xt2 +
- gt33L*PDstandardNth3Xt3 + Gtl313*Xtn1 + Gtl323*Xtn2 + Gtl333*Xtn3;
+ CCTK_REAL Rt33 = 0.5*(4*(Gt113*Gtlu311 + Gt123*Gtlu312 + Gt133*Gtlu313
+ + Gt213*Gtlu321 + Gt223*Gtlu322 + Gt233*Gtlu323) + 6*(Gt313*Gtlu331 +
+ Gt323*Gtlu332 + Gt333*Gtlu333) - gtu11*PDstandardNth11gt33 -
+ 2*gtu12*PDstandardNth12gt33 - 2*gtu13*PDstandardNth13gt33 -
+ gtu22*PDstandardNth22gt33 - 2*gtu23*PDstandardNth23gt33 -
+ gtu33*PDstandardNth33gt33 + 2*(Gt113*Gtlu131 + Gt123*Gtlu132 +
+ Gt133*Gtlu133 + Gt213*Gtlu231 + Gt223*Gtlu232 + Gt233*Gtlu233 +
+ gt13L*PDstandardNth3Xt1) + 2*gt23L*PDstandardNth3Xt2 +
+ 2*gt33L*PDstandardNth3Xt3 + 2*Gtl313*Xtn1 + 2*Gtl323*Xtn2 +
+ 2*Gtl333*Xtn3);
- CCTK_REAL fac1 = IfThen(ToReal(conformalMethod),-0.5*INV(phiL),1);
+ CCTK_REAL fac1 = IfThen(conformalMethod,-0.5*INV(phiL),1);
CCTK_REAL cdphi1 = fac1*PDstandardNth1phi;
@@ -474,7 +477,7 @@ static void ML_BSSN_constraints1_Body(cGH const * restrict const cctkGH, int con
CCTK_REAL cdphi3 = fac1*PDstandardNth3phi;
- CCTK_REAL fac2 = IfThen(ToReal(conformalMethod),0.5*INV(SQR(phiL)),0);
+ CCTK_REAL fac2 = IfThen(conformalMethod,0.5*INV(SQR(phiL)),0);
CCTK_REAL cdphi211 = -(fac1*(-PDstandardNth11phi +
Gt111*PDstandardNth1phi + Gt211*PDstandardNth2phi +
@@ -536,8 +539,7 @@ static void ML_BSSN_constraints1_Body(cGH const * restrict const cctkGH, int con
+ 2*SQR(cdphi1)) + gtu22*(cdphi222 + 2*SQR(cdphi2))) + 2*(-1 +
gt33L*gtu33)*SQR(cdphi3));
- CCTK_REAL e4phi =
- IfThen(ToReal(conformalMethod),INV(SQR(phiL)),exp(4*phiL));
+ CCTK_REAL e4phi = IfThen(conformalMethod,INV(SQR(phiL)),exp(4*phiL));
CCTK_REAL em4phi = INV(e4phi);
diff --git a/ML_BSSN/src/ML_BSSN_constraints2.cc b/ML_BSSN/src/ML_BSSN_constraints2.cc
index f9cd342..bc123f0 100644
--- a/ML_BSSN/src/ML_BSSN_constraints2.cc
+++ b/ML_BSSN/src/ML_BSSN_constraints2.cc
@@ -305,7 +305,7 @@ static void ML_BSSN_constraints2_Body(cGH const * restrict const cctkGH, int con
CCTK_REAL Gt333 = Gtl133*gtu13 + Gtl233*gtu23 + Gtl333*gtu33;
- CCTK_REAL fac1 = IfThen(ToReal(conformalMethod),-0.5*INV(phiL),1);
+ CCTK_REAL fac1 = IfThen(conformalMethod,-0.5*INV(phiL),1);
CCTK_REAL cdphi1 = fac1*PDstandardNth1phi;
diff --git a/ML_BSSN/src/ML_BSSN_convertFromADMBaseGamma.cc b/ML_BSSN/src/ML_BSSN_convertFromADMBaseGamma.cc
index f89f8b1..917a7e9 100644
--- a/ML_BSSN/src/ML_BSSN_convertFromADMBaseGamma.cc
+++ b/ML_BSSN/src/ML_BSSN_convertFromADMBaseGamma.cc
@@ -282,36 +282,49 @@ static void ML_BSSN_convertFromADMBaseGamma_Body(cGH const * restrict const cctk
CCTK_REAL Xt3L = Gt311*gtu11 + Gt322*gtu22 + 2*(Gt312*gtu12 +
Gt313*gtu13 + Gt323*gtu23) + Gt333*gtu33;
- CCTK_REAL AL = IfThen(ToReal(LapseACoeff) !=
- 0,-3*INV(ToReal(harmonicF))*pow(alphaL,-ToReal(harmonicN))*(9*dtalpL -
+ CCTK_REAL AL = IfThen(LapseACoeff !=
+ 0,-(INV(ToReal(harmonicF))*pow(alphaL,-ToReal(harmonicN))*(dtalpL -
(beta1L*PDupwindNthAnti1alpha + beta2L*PDupwindNthAnti2alpha +
beta3L*PDupwindNthAnti3alpha + PDupwindNthSymm1alpha*Abs(beta1L) +
PDupwindNthSymm2alpha*Abs(beta2L) +
- PDupwindNthSymm3alpha*Abs(beta3L))*ToReal(LapseAdvectionCoeff)),0);
+ PDupwindNthSymm3alpha*Abs(beta3L))*ToReal(LapseAdvectionCoeff))),0);
CCTK_REAL theta = fmin(1,exp(1 -
rL*INV(ToReal(SpatialShiftGammaCoeffRadius))));
- CCTK_REAL B1L = IfThen(ToReal(ShiftBCoeff)*ToReal(ShiftGammaCoeff) !=
- 0,INV(theta)*INV(ToReal(ShiftGammaCoeff))*(27*dtbetaxL -
- 3*(beta1L*PDupwindNthAnti1beta1 + beta2L*PDupwindNthAnti2beta1 +
- beta3L*PDupwindNthAnti3beta1 + PDupwindNthSymm1beta1*Abs(beta1L) +
- PDupwindNthSymm2beta1*Abs(beta2L) +
- PDupwindNthSymm3beta1*Abs(beta3L))*ToReal(ShiftAdvectionCoeff)),0);
-
- CCTK_REAL B2L = IfThen(ToReal(ShiftBCoeff)*ToReal(ShiftGammaCoeff) !=
- 0,INV(theta)*INV(ToReal(ShiftGammaCoeff))*(27*dtbetayL -
- 3*(beta1L*PDupwindNthAnti1beta2 + beta2L*PDupwindNthAnti2beta2 +
- beta3L*PDupwindNthAnti3beta2 + PDupwindNthSymm1beta2*Abs(beta1L) +
- PDupwindNthSymm2beta2*Abs(beta2L) +
- PDupwindNthSymm3beta2*Abs(beta3L))*ToReal(ShiftAdvectionCoeff)),0);
-
- CCTK_REAL B3L = IfThen(ToReal(ShiftBCoeff)*ToReal(ShiftGammaCoeff) !=
- 0,INV(theta)*INV(ToReal(ShiftGammaCoeff))*(27*dtbetazL -
- 3*(beta1L*PDupwindNthAnti1beta3 + beta2L*PDupwindNthAnti2beta3 +
- beta3L*PDupwindNthAnti3beta3 + PDupwindNthSymm1beta3*Abs(beta1L) +
- PDupwindNthSymm2beta3*Abs(beta2L) +
- PDupwindNthSymm3beta3*Abs(beta3L))*ToReal(ShiftAdvectionCoeff)),0);
+ CCTK_REAL B1L;
+ CCTK_REAL B2L;
+ CCTK_REAL B3L;
+
+ if (ShiftBCoeff*ShiftGammaCoeff != 0)
+ {
+ B1L = INV(theta)*INV(ToReal(ShiftGammaCoeff))*(dtbetaxL -
+ (beta1L*PDupwindNthAnti1beta1 + beta2L*PDupwindNthAnti2beta1 +
+ beta3L*PDupwindNthAnti3beta1 + PDupwindNthSymm1beta1*Abs(beta1L) +
+ PDupwindNthSymm2beta1*Abs(beta2L) +
+ PDupwindNthSymm3beta1*Abs(beta3L))*ToReal(ShiftAdvectionCoeff));
+
+ B2L = INV(theta)*INV(ToReal(ShiftGammaCoeff))*(dtbetayL -
+ (beta1L*PDupwindNthAnti1beta2 + beta2L*PDupwindNthAnti2beta2 +
+ beta3L*PDupwindNthAnti3beta2 + PDupwindNthSymm1beta2*Abs(beta1L) +
+ PDupwindNthSymm2beta2*Abs(beta2L) +
+ PDupwindNthSymm3beta2*Abs(beta3L))*ToReal(ShiftAdvectionCoeff));
+
+ B3L = INV(theta)*INV(ToReal(ShiftGammaCoeff))*(dtbetazL -
+ (beta1L*PDupwindNthAnti1beta3 + beta2L*PDupwindNthAnti2beta3 +
+ beta3L*PDupwindNthAnti3beta3 + PDupwindNthSymm1beta3*Abs(beta1L) +
+ PDupwindNthSymm2beta3*Abs(beta2L) +
+ PDupwindNthSymm3beta3*Abs(beta3L))*ToReal(ShiftAdvectionCoeff));
+ }
+ else
+ {
+ B1L = 0;
+
+ B2L = 0;
+
+ B3L = 0;
+ }
+
/* Copy local copies back to grid functions */
A[index] = AL;
diff --git a/ML_BSSN/src/ML_BSSN_convertToADMBase.cc b/ML_BSSN/src/ML_BSSN_convertToADMBase.cc
index 804af11..589ccc2 100644
--- a/ML_BSSN/src/ML_BSSN_convertToADMBase.cc
+++ b/ML_BSSN/src/ML_BSSN_convertToADMBase.cc
@@ -129,8 +129,7 @@ static void ML_BSSN_convertToADMBase_Body(cGH const * restrict const cctkGH, int
/* Precompute derivatives */
/* Calculate temporaries and grid functions */
- CCTK_REAL e4phi =
- IfThen(ToReal(conformalMethod),INV(SQR(phiL)),exp(4*phiL));
+ CCTK_REAL e4phi = IfThen(conformalMethod,INV(SQR(phiL)),exp(4*phiL));
gxxL = e4phi*gt11L;