aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorrhaas <rhaas@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2014-07-29 17:06:02 +0000
committerrhaas <rhaas@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2014-07-29 17:06:02 +0000
commit4d32b394d46d3b543f2575390e9d14b3b9025d7b (patch)
tree345be1f3640ce179a93a6e52c33c3ecbf40f268b
parent352240b3264d7655ca3f7d659735fc30b9fffe0b (diff)
check that init_RHS_zero is used for RK?-MR-2-1 methodssvn
these methods assume that the RHS grid function values can be re-used in substeps which is no longer true after MoL revision 225 git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@226 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
-rw-r--r--param.ccl4
-rw-r--r--src/ParamCheck.c30
2 files changed, 24 insertions, 10 deletions
diff --git a/param.ccl b/param.ccl
index c7e134e..4406721 100644
--- a/param.ccl
+++ b/param.ccl
@@ -100,8 +100,8 @@ KEYWORD ODE_Method "The ODE method use by MoL to do time integration"
"RK65" :: "RK65 with error estimation"
"RK87" :: "RK87 with error estimation"
"AB" :: "Adams-Bashforth"
- "RK2-MR-2:1" :: "2nd order 2:1 multirate RK scheme based on RK2 due to Schlegel et al 2009"
- "RK4-MR-2:1" :: "3rd order 2:1 multirate RK scheme based on RK43 due to Schlegel et al 2009"
+ "RK2-MR-2:1" :: "2nd order 2:1 multirate RK scheme based on RK2 due to Schlegel et al 2009. This requires init_RHS_zero='no'."
+ "RK4-MR-2:1" :: "3rd order 2:1 multirate RK scheme based on RK43 due to Schlegel et al 2009. This requires init_RHS_zero='no'."
"RK4-RK2" :: "RK4 as fast method and RK2 as slow method"
} "ICN"
diff --git a/src/ParamCheck.c b/src/ParamCheck.c
index c46343a..a1b313a 100644
--- a/src/ParamCheck.c
+++ b/src/ParamCheck.c
@@ -172,18 +172,32 @@ void MoL_ParamCheck(CCTK_ARGUMENTS)
"number of intermediate steps must be 1");
}
- if ( (CCTK_Equals(ODE_Method, "RK2-MR-2:1"))&&
- ( !((MoL_Intermediate_Steps == 5) && (MoL_Num_Scratch_Levels > 4))) )
+ if ( CCTK_Equals(ODE_Method, "RK2-MR-2:1") )
{
- CCTK_PARAMWARN("When using the multirate 2-1 RK2 evolver the "
- "number of intermediate steps must be 5 and the number of scratch levels at least 5");
+ if ( !((MoL_Intermediate_Steps == 5) && (MoL_Num_Scratch_Levels > 4)) )
+ {
+ CCTK_PARAMWARN("When using the multirate 2-1 RK2 evolver the "
+ "number of intermediate steps must be 5 and the number of scratch levels at least 5");
+ }
+ if (init_RHS_zero)
+ {
+ CCTK_PARAMWARN("When using the multirate 2-1 RK2 evolver the "
+ "parameter MoL::init_RHS_zero must be set to 'no'.");
+ }
}
- if ( (CCTK_Equals(ODE_Method, "RK4-MR-2:1"))&&
- ( !((MoL_Intermediate_Steps == 10) && (MoL_Num_Scratch_Levels > 9))) )
+ if ( CCTK_Equals(ODE_Method, "RK4-MR-2:1") )
{
- CCTK_PARAMWARN("When using the multirate 2-1 RK4 evolver the "
- "number of intermediate steps must be 10 and the number of scratch levels at least 10");
+ if ( !((MoL_Intermediate_Steps == 10) && (MoL_Num_Scratch_Levels > 9)) )
+ {
+ CCTK_PARAMWARN("When using the multirate 2-1 RK4 evolver the "
+ "number of intermediate steps must be 10 and the number of scratch levels at least 10");
+ }
+ if (init_RHS_zero)
+ {
+ CCTK_PARAMWARN("When using the multirate 2-1 RK4 evolver the "
+ "parameter MoL::init_RHS_zero must be set to 'no'.");
+ }
}
if ( (CCTK_Equals(ODE_Method, "RK4-RK2"))&&