diff options
author | Erik Schnetter <schnetter@cct.lsu.edu> | 2011-01-03 23:04:03 -0500 |
---|---|---|
committer | Barry Wardell <barry.wardell@gmail.com> | 2011-12-14 18:25:53 +0000 |
commit | d707da354918f3003420f26c2c2cc31fe1540913 (patch) | |
tree | fd90463042e5b4f1aff3540259c3519dc784d73f /Carpet/Carpet/src/CallFunction.cc | |
parent | 33bba04684049f70919c71e3b5ca5efcce2d99bd (diff) |
Carpet: Set times of past timelevels correctly
Redesign the way in which the times of the past timelevels are set. If cctk_delta_time changes during initialisation, re-set the times of the past timelevels accordingly.
Diffstat (limited to 'Carpet/Carpet/src/CallFunction.cc')
-rw-r--r-- | Carpet/Carpet/src/CallFunction.cc | 30 |
1 files changed, 30 insertions, 0 deletions
diff --git a/Carpet/Carpet/src/CallFunction.cc b/Carpet/Carpet/src/CallFunction.cc index ffd695c45..73c34d73c 100644 --- a/Carpet/Carpet/src/CallFunction.cc +++ b/Carpet/Carpet/src/CallFunction.cc @@ -1,6 +1,7 @@ #include <algorithm> #include <cassert> #include <cstdlib> +#include <cstring> #include <map> #include <string> @@ -342,6 +343,9 @@ namespace Carpet { } Timer & timer = * ti->second; + // Save the time step size + CCTK_REAL const saved_cctk_delta_time = cctkGH->cctk_delta_time; + user_timer.start(); timer.start(); int const res = CCTK_CallFunction (function, attribute, data); @@ -349,6 +353,32 @@ namespace Carpet { timer.stop(); user_timer.stop(); + // Manage the time step size. If the time step size changes + // during initialisation, assume it is thorn Time, and update + // the time hierarchy. If it changes during evolution, assume it + // is MoL, and do nothing. + if (cctkGH->cctk_iteration == 0 and + cctkGH->cctk_delta_time != saved_cctk_delta_time) + { + // The user changed cctk_delta_time during initialisation -- + // update our internals and the time hierarchy + delta_time = cctkGH->cctk_delta_time * timereflevelfact / mglevelfact; + for (int ml=0; ml<mglevels; ++ml) { + for (int rl=0; rl<reflevels; ++rl) { + // Update the time delta + CCTK_REAL const dt = + delta_time / timereffacts.AT(rl) * ipow(mgfact, ml); + tt->set_delta(ml,rl,dt); + CCTK_REAL const t0 = tt->get_time(ml,rl,0); + // Update the times of the past timelevels + for (int tl=1; tl<timelevels; ++tl) { + CCTK_REAL const t = t0 - tl * dt; + tt->set_time(ml,rl,tl,t); + } + } + } + } + } CallAfterRoutines (cctkGH, function, attribute, data); } |