aboutsummaryrefslogtreecommitdiff
path: root/src/SetTime.c
diff options
context:
space:
mode:
authorhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2006-01-23 10:39:58 +0000
committerhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2006-01-23 10:39:58 +0000
commitc05ef7225ab51b7e66e884cf4d53805e8d261982 (patch)
tree5885d1761161ada1381942a09dfff646f872fc9c /src/SetTime.c
parent0863f0decdf6024f548f6178ce6221bc1a6fa722 (diff)
Peter Diener's RK65 and RK87 adaptive timestep integrators.
As yet not altered to do grid arrays. As with RK45, adaptive timestepping does not work with mesh refinement. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@106 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/SetTime.c')
-rw-r--r--src/SetTime.c59
1 files changed, 59 insertions, 0 deletions
diff --git a/src/SetTime.c b/src/SetTime.c
index af044f4..1d1596d 100644
--- a/src/SetTime.c
+++ b/src/SetTime.c
@@ -57,6 +57,35 @@ static const CCTK_REAL alpha_array[6] = {
1.0/2.0,
};
+/* RK65 coefficients */
+static const CCTK_REAL alpha_array65[8] = {
+ 0.0,
+ 1.0/10.0,
+ 2.0/9.0,
+ 3.0/7.0,
+ 3.0/5.0,
+ 4.0/5.0,
+ 1.0,
+ 1.0
+};
+
+/* RK87 coefficients */
+ static const CCTK_REAL alpha_array87[13] = {
+ 0.0,
+ 1.0/18.0,
+ 1.0/12.0,
+ 1.0/8.0,
+ 5.0/16.0,
+ 3.0/8.0,
+ 59.0/400.0,
+ 93.0/200.0,
+ 5490023248.0/9719169821.0,
+ 13.0/20.0,
+ 1201146811.0/1299019798.0,
+ 1.0,
+ 1.0
+};
+
/********************************************************************
********************* External Routines **********************
********************************************************************/
@@ -216,6 +245,22 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
+ ((alpha_array[substep] - 1)
* (* Original_Delta_Time) / cctkGH->cctk_timefac));
}
+ else if (CCTK_EQUALS(ODE_Method,"RK65"))
+ {
+ const int substep = MoL_Intermediate_Steps - (* MoL_Intermediate_Step);
+ cctkGH->cctk_time
+ = ((* Original_Time)
+ + ((alpha_array65[substep] - 1)
+ * (* Original_Delta_Time) / cctkGH->cctk_timefac));
+ }
+ else if (CCTK_EQUALS(ODE_Method,"RK87"))
+ {
+ const int substep = MoL_Intermediate_Steps - (* MoL_Intermediate_Step);
+ cctkGH->cctk_time
+ = ((* Original_Time)
+ + ((alpha_array87[substep] - 1)
+ * (* Original_Delta_Time) / cctkGH->cctk_timefac));
+ }
#ifdef MOLDEBUG
printf("MoL has once more reset t (%d): %f.\n",
*MoL_Intermediate_Step, cctkGH->cctk_time);
@@ -301,6 +346,20 @@ int MoL_ResetDeltaTime(CCTK_ARGUMENTS)
= ((alpha_array[substep + 1] - alpha_array[substep])
* (* Original_Delta_Time));
}
+ else if (CCTK_EQUALS(ODE_Method,"RK65"))
+ {
+ const int substep = MoL_Intermediate_Steps - (* MoL_Intermediate_Step);
+ cctkGH->cctk_delta_time
+ = ((alpha_array65[substep + 1] - alpha_array65[substep])
+ * (* Original_Delta_Time));
+ }
+ else if (CCTK_EQUALS(ODE_Method,"RK87"))
+ {
+ const int substep = MoL_Intermediate_Steps - (* MoL_Intermediate_Step);
+ cctkGH->cctk_delta_time
+ = ((alpha_array87[substep + 1] - alpha_array87[substep])
+ * (* Original_Delta_Time));
+ }
#ifdef MOLDEBUG
printf("MoL has once more reset dt (%d): %f.\n",
*MoL_Intermediate_Step,