aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2006-08-24 11:45:17 +0000
committerhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2006-08-24 11:45:17 +0000
commite76b9ee4215de4ca9258261646421ef2adfd2c19 (patch)
tree99e73891daa44a25f3848d816620e6fd40b4d05e
parent8b15c1ee0334cd9008376280b91c8e40723f3c7f (diff)
Theta ICN. Patch from Frank Loeffler.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@120 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
-rw-r--r--param.ccl9
-rw-r--r--src/ICN.c22
2 files changed, 27 insertions, 4 deletions
diff --git a/param.ccl b/param.ccl
index bd8e86f..7e29e56 100644
--- a/param.ccl
+++ b/param.ccl
@@ -102,6 +102,15 @@ KEYWORD Generic_Type "If using the generic method, which sort"
"Classic RK3" :: "Efficient RK3 - classical version"
} "RK"
+CCTK_REAL ICN_avg_theta "theta of averaged ICN method, usually 0.5"
+{
+ 0:1 :: "0 <= theta <= 1"
+} 0.5
+
+BOOLEAN ICN_avg_swapped "Use swapped averages in ICN method?"
+{
+} "no"
+
CCTK_INT MoL_Intermediate_Steps "Number of intermediate steps taken by the ODE method"
{
1:* :: "Anything greater than 1"
diff --git a/src/ICN.c b/src/ICN.c
index 7e0acde..7f75a35 100644
--- a/src/ICN.c
+++ b/src/ICN.c
@@ -221,6 +221,8 @@ void MoL_ICNAverage(CCTK_ARGUMENTS)
/* FIXME */
+ CCTK_REAL theta;
+
#ifdef MOLDOESCOMPLEX
CCTK_COMPLEX *OldComplexVar;
@@ -230,15 +232,23 @@ void MoL_ICNAverage(CCTK_ARGUMENTS)
#endif
+
+ theta = ICN_avg_theta;
+ if (ICN_avg_swapped && (*MoL_Intermediate_Step%2))
+ {
+ theta = 1.0 - theta;
+ }
+
#ifdef MOLDEBUG
printf("Inside ICN.\nProcessor %d.\nStep %d.\nRefinement %d.\n"
- "Timestep %g.\nSpacestep %g.\nTime %g\n",
+ "Timestep %g.\nSpacestep %g.\nTime %g Theta %g\n",
CCTK_MyProc(cctkGH),
MoL_Intermediate_Steps - *MoL_Intermediate_Step + 1,
*cctk_levfac,
CCTK_DELTA_TIME,
CCTK_DELTA_SPACE(0),
- cctk_time);
+ cctk_time,
+ theta);
#endif
totalsize = 1;
@@ -261,7 +271,9 @@ void MoL_ICNAverage(CCTK_ARGUMENTS)
for (index = 0; index < totalsize; index++)
{
- UpdateVar[index] = 0.5 * (UpdateVar[index] + OldVar[index]);
+/* UpdateVar[index] = 0.5 * (UpdateVar[index] + OldVar[index]); */
+ UpdateVar[index] = (1.0 - theta) * UpdateVar[index] +
+ theta * OldVar[index];
}
}
@@ -292,7 +304,9 @@ void MoL_ICNAverage(CCTK_ARGUMENTS)
for (index = 0; index < arraytotalsize; index++)
{
- UpdateVar[index] = 0.5 * (UpdateVar[index] + OldVar[index]);
+/* UpdateVar[index] = 0.5 * (UpdateVar[index] + OldVar[index]); */
+ UpdateVar[index] = (1.0 - theta) * UpdateVar[index] +
+ theta * OldVar[index];
}
}