aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2003-06-20 15:07:23 +0000
committerhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2003-06-20 15:07:23 +0000
commit8b767a625c0366ee60d55ca483abc793e28a0d68 (patch)
tree4b29433a4a02eba8927d4ce9869bd928a662b2f2 /src
parent00f602b5277ef0e9f7be9743d3484d7b251df739 (diff)
Correct a major screw up with the setting of t/dt that was causing the BSSN_MoL testsuite bh_shift_rad to fail.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@18 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src')
-rw-r--r--src/ICN.c2
-rw-r--r--src/SetTime.c79
2 files changed, 64 insertions, 17 deletions
diff --git a/src/ICN.c b/src/ICN.c
index 3117522..0a1cced 100644
--- a/src/ICN.c
+++ b/src/ICN.c
@@ -89,7 +89,7 @@ void MoL_ICNAdd(CCTK_ARGUMENTS)
CCTK_COMPLEX Complex_Delta_Time = CCTK_Cmplx(CCTK_DELTA_TIME, 0);
#endif
-
+
#ifdef MOLDEBUG
printf("Inside ICN.\nProcessor %d.\nStep %d.\nRefinement %d.\nTimestep %g.\nSpacestep %g.\nTime %g\n",
CCTK_MyProc(cctkGH),
diff --git a/src/SetTime.c b/src/SetTime.c
index 7774102..718ef44 100644
--- a/src/SetTime.c
+++ b/src/SetTime.c
@@ -35,6 +35,8 @@ int MoL_SetTime(CCTK_ARGUMENTS);
int MoL_ResetTime(CCTK_ARGUMENTS);
+int MoL_ResetDeltaTime(CCTK_ARGUMENTS);
+
/********************************************************************
********************* Other Routine Prototypes *********************
********************************************************************/
@@ -93,7 +95,7 @@ int MoL_SetTime(CCTK_ARGUMENTS)
@date Mon May 20 09:49:41 2002
@author Ian Hawke
@desc
- Sets the time and timestep during the MoL evolution loop.
+ Sets the time during the MoL evolution loop.
At the last time all methods should end up with the original
values for time and timestep.
@enddesc
@@ -123,25 +125,13 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
if (*MoL_Intermediate_Step == 0)
{
cctkGH->cctk_time = (*Original_Time);
- cctkGH->cctk_delta_time = (*Original_Delta_Time);
}
else if (CCTK_EQUALS(ODE_Method,"ICN"))
{
- if (*MoL_Intermediate_Step == 1)
- {
- cctkGH->cctk_delta_time = (*Original_Delta_Time);
- }
- else
- {
- cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
- }
cctkGH->cctk_time = (*Original_Time)-0.5*(*Original_Delta_Time)/cctkGH->cctk_timefac;
}
else if (CCTK_EQUALS(ODE_Method,"Generic"))
{
- cctkGH->cctk_delta_time = RKBetaCoefficients[MoL_Intermediate_Steps -
- (*MoL_Intermediate_Step)] *
- (*Original_Delta_Time)/cctkGH->cctk_timefac;
previous_times[0] = (*Original_Time) - (*Original_Delta_Time)/cctkGH->cctk_timefac;
for (i = MoL_Intermediate_Steps - 1; i > *MoL_Intermediate_Step; i--)
{
@@ -162,13 +152,11 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
{
if (*MoL_Intermediate_Step == 1)
{
- cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
cctkGH->cctk_time = (*Original_Time);
}
}
#ifdef MOLDEBUG
- printf("MoL has once more reset t: %f: and dt: %f.\n", cctkGH->cctk_time,
- cctkGH->cctk_delta_time/cctkGH->cctk_timefac);
+ printf("MoL has once more reset t: %f.\n", cctkGH->cctk_time);
#endif
free(previous_times);
@@ -177,6 +165,65 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
return 0;
}
+ /*@@
+ @routine MoL_ResetDeltaTime
+ @date Mon May 20 09:49:41 2002
+ @author Ian Hawke
+ @desc
+ Sets the timestep during the MoL evolution loop.
+ At the last time all methods should end up with the original
+ values for time and timestep.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+int MoL_ResetDeltaTime(CCTK_ARGUMENTS)
+{
+
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+ if (*MoL_Intermediate_Step == 0)
+ {
+ cctkGH->cctk_delta_time = (*Original_Delta_Time);
+ }
+ else if (CCTK_EQUALS(ODE_Method,"ICN"))
+ {
+ if (*MoL_Intermediate_Step == 1)
+ {
+ cctkGH->cctk_delta_time = (*Original_Delta_Time);
+ }
+ else
+ {
+ cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
+ }
+ }
+ else if (CCTK_EQUALS(ODE_Method,"Generic"))
+ {
+ cctkGH->cctk_delta_time = RKBetaCoefficients[MoL_Intermediate_Steps -
+ (*MoL_Intermediate_Step)] *
+ (*Original_Delta_Time)/cctkGH->cctk_timefac;
+ }
+ else if (CCTK_EQUALS(ODE_Method,"RK2"))
+ {
+ if (*MoL_Intermediate_Step == 1)
+ {
+ cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
+ }
+ }
+#ifdef MOLDEBUG
+ printf("MoL has once more reset dt: %f.\n",
+ cctkGH->cctk_delta_time/cctkGH->cctk_timefac);
+#endif
+
+ return 0;
+}
+
/********************************************************************
********************* Local Routines *************************
********************************************************************/