aboutsummaryrefslogtreecommitdiff
path: root/src/SetTime.c
diff options
context:
space:
mode:
Diffstat (limited to 'src/SetTime.c')
-rw-r--r--src/SetTime.c44
1 files changed, 40 insertions, 4 deletions
diff --git a/src/SetTime.c b/src/SetTime.c
index 1d1596d..3e9ae3f 100644
--- a/src/SetTime.c
+++ b/src/SetTime.c
@@ -47,8 +47,8 @@ int MoL_ResetDeltaTime(CCTK_ARGUMENTS);
********************* Local Data *****************************
********************************************************************/
-/* RK45 coefficients */
-static const CCTK_REAL alpha_array[6] = {
+/* RK45 Fehlberg coefficients */
+static const CCTK_REAL alpha_array_F[6] = {
0.0,
1.0/4.0,
3.0/8.0,
@@ -57,6 +57,16 @@ static const CCTK_REAL alpha_array[6] = {
1.0/2.0,
};
+/* RK45 Cash-Karp coefficients */
+static const CCTK_REAL alpha_array_CK[6] = {
+ 0.0,
+ 1.0/5.0,
+ 3.0/10.0,
+ 3.0/5.0,
+ 1.0,
+ 7.0/8.0,
+};
+
/* RK65 coefficients */
static const CCTK_REAL alpha_array65[8] = {
0.0,
@@ -237,9 +247,22 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
0.5*(*Original_Delta_Time)/cctkGH->cctk_timefac;
}
}
- else if (CCTK_EQUALS(ODE_Method,"RK45"))
+ else if (CCTK_EQUALS(ODE_Method,"RK45") || CCTK_EQUALS(ODE_Method,"RK45CK"))
{
const int substep = MoL_Intermediate_Steps - (* MoL_Intermediate_Step);
+ const CCTK_REAL * alpha_array;
+ if (CCTK_EQUALS(ODE_Method, "RK45"))
+ {
+ alpha_array = alpha_array_F;
+ }
+ else if (CCTK_EQUALS(ODE_Method, "RK45CK"))
+ {
+ alpha_array = alpha_array_CK;
+ }
+ else
+ {
+ CCTK_WARN (0, "internal error");
+ }
cctkGH->cctk_time
= ((* Original_Time)
+ ((alpha_array[substep] - 1)
@@ -339,9 +362,22 @@ int MoL_ResetDeltaTime(CCTK_ARGUMENTS)
cctkGH->cctk_delta_time = 2.0/3.0*(*Original_Delta_Time);
}
}
- else if (CCTK_EQUALS(ODE_Method,"RK45"))
+ else if (CCTK_EQUALS(ODE_Method,"RK45") || CCTK_EQUALS(ODE_Method,"RK45CK"))
{
const int substep = MoL_Intermediate_Steps - (* MoL_Intermediate_Step);
+ const CCTK_REAL * alpha_array;
+ if (CCTK_EQUALS(ODE_Method, "RK45"))
+ {
+ alpha_array = alpha_array_F;
+ }
+ else if (CCTK_EQUALS(ODE_Method, "RK45CK"))
+ {
+ alpha_array = alpha_array_CK;
+ }
+ else
+ {
+ CCTK_WARN (0, "internal error");
+ }
cctkGH->cctk_delta_time
= ((alpha_array[substep + 1] - alpha_array[substep])
* (* Original_Delta_Time));