diff options
author | hawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b> | 2003-04-25 15:37:35 +0000 |
---|---|---|
committer | hawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b> | 2003-04-25 15:37:35 +0000 |
commit | 79bd8857ecf28b9a17255af76b5e08785d9f9a94 (patch) | |
tree | c35b9132d7a1b14ae5bf8d91bccc8be3c391c00b /src/SetTime.c | |
parent | afb053800f476ae76546dc5709f145a5793ff336 (diff) |
The resetting of t and dt is corrected for the generic methods.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@7 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/SetTime.c')
-rw-r--r-- | src/SetTime.c | 98 |
1 files changed, 39 insertions, 59 deletions
diff --git a/src/SetTime.c b/src/SetTime.c index 4e99c22..78af11b 100644 --- a/src/SetTime.c +++ b/src/SetTime.c @@ -66,17 +66,23 @@ int MoL_SetTime(CCTK_ARGUMENTS) DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; + CCTK_REAL beta; + *Original_Time = cctkGH->cctk_time; *Original_Delta_Time = cctkGH->CCTK_DELTA_TIME; cctkGH->cctk_time -= cctkGH->CCTK_DELTA_TIME; - if ( (CCTK_EQUALS(MoL_ODE_Method,"ICN")) || - ( (CCTK_EQUALS(MoL_ODE_Method,"Generic")) && - CCTK_EQUALS(Generic_Type,"ICN") ) ) + if (CCTK_EQUALS(MoL_ODE_Method,"ICN")) { cctkGH->CCTK_DELTA_TIME = 0.5*(*Original_Delta_Time); } + else if (CCTK_EQUALS(MoL_ODE_Method,"Generic")) + { + beta = RKBetaCoefficients[0]; + cctkGH->CCTK_DELTA_TIME = beta*(*Original_Delta_Time); + } + return 0; } @@ -102,81 +108,55 @@ int MoL_ResetTime(CCTK_ARGUMENTS) DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; - - + + CCTK_REAL previous_times[MoL_Intermediate_Steps]; + CCTK_INT alphaindex, i, j; + if (*MoL_Intermediate_Step == 0) { cctkGH->cctk_time = (*Original_Time); cctkGH->CCTK_DELTA_TIME = (*Original_Delta_Time); } - else if (*MoL_Intermediate_Step == 1) + else if (CCTK_EQUALS(MoL_ODE_Method,"ICN")) { - if ( (CCTK_EQUALS(MoL_ODE_Method,"Generic")) - &&(CCTK_EQUALS(Generic_Type,"RK")) ) + if (*MoL_Intermediate_Step == 1) { - if (MoL_Intermediate_Steps == 2) - { - cctkGH->cctk_time = (*Original_Time); - } - /* else if (MoL_Intermediate_Steps == 3) - { - cctkGH->CCTK_DELTA_TIME = 2.0 / 3.0 * (*Original_Delta_Time); - cctkGH->cctk_time = (*Original_Time) - 0.5 * (*Original_Delta_Time); - } - else if (MoL_Intermediate_Steps == 4) - { - cctkGH->CCTK_DELTA_TIME = (*Original_Delta_Time) / 6.0; - cctkGH->cctk_time = (*Original_Time); - }*/ + cctkGH->CCTK_DELTA_TIME = (*Original_Delta_Time); } - else if ( ((CCTK_EQUALS(MoL_ODE_Method,"Generic"))&& - (CCTK_EQUALS(Generic_Type,"ICN"))) || - (CCTK_EQUALS(MoL_ODE_Method,"ICN")) ) + else { - cctkGH->CCTK_DELTA_TIME = *Original_Delta_Time; - cctkGH->cctk_time = (*Original_Time)-0.5*(*Original_Delta_Time); - } - else if ( CCTK_EQUALS(MoL_ODE_Method,"RK2") ) - { - cctkGH->CCTK_DELTA_TIME = 0.5 * (*Original_Delta_Time); - cctkGH->cctk_time = (*Original_Time); + cctkGH->CCTK_DELTA_TIME = 0.5*(*Original_Delta_Time); } + cctkGH->cctk_time = (*Original_Time)-0.5*(*Original_Delta_Time); } - else if (*MoL_Intermediate_Step == 2) + else if (CCTK_EQUALS(MoL_ODE_Method,"Generic")) { - if ( (CCTK_EQUALS(MoL_ODE_Method,"Generic")) - &&(CCTK_EQUALS(Generic_Type,"RK")) ) + cctkGH->CCTK_DELTA_TIME = RKBetaCoefficients[MoL_Intermediate_Steps - + (*MoL_Intermediate_Step)] * + (*Original_Delta_Time); + previous_times[0] = (*Original_Time) - (*Original_Delta_Time); + for (i = MoL_Intermediate_Steps - 1; i > *MoL_Intermediate_Step; i--) { - /* if (MoL_Intermediate_Steps == 3) + for (j = MoL_Intermediate_Steps; j > i; j--) { - cctkGH->CCTK_DELTA_TIME = 0.25 * (*Original_Delta_Time); - cctkGH->cctk_time = (*Original_Time); + alphaindex = (MoL_Intermediate_Steps - i) * + MoL_Intermediate_Steps + j; + previous_times[MoL_Intermediate_Steps - i] = + RKBetaCoefficients[MoL_Intermediate_Steps - j] * + (*Original_Delta_Time) + + RKAlphaCoefficients[alphaindex]; } - else if (MoL_Intermediate_Steps == 4) - { - cctkGH->CCTK_DELTA_TIME = (*Original_Delta_Time); - cctkGH->cctk_time = (*Original_Time); - }*/ - } - else if ( ((CCTK_EQUALS(MoL_ODE_Method,"Generic"))&& - (CCTK_EQUALS(Generic_Type,"ICN"))) || - (CCTK_EQUALS(MoL_ODE_Method,"ICN")) ) - { - cctkGH->CCTK_DELTA_TIME = 0.5*(*Original_Delta_Time); - cctkGH->cctk_time = (*Original_Time)-0.5*(*Original_Delta_Time); } + cctkGH->cctk_time = previous_times[MoL_Intermediate_Steps - + *MoL_Intermediate_Step]; } - else if (*MoL_Intermediate_Step == 3) + else if (CCTK_EQUALS(MoL_ODE_Method,"RK2")) { - /* if ( (CCTK_EQUALS(MoL_ODE_Method,"Generic")) - &&(CCTK_EQUALS(Generic_Type,"RK")) ) + if (*MoL_Intermediate_Step == 1) { - if (MoL_Intermediate_Steps == 4) - { - cctkGH->CCTK_DELTA_TIME = 0.5 * (*Original_Delta_Time); - cctkGH->cctk_time = (*Original_Time) - 0.5 * (*Original_Delta_Time); - } - }*/ + 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, |