diff options
author | Erik Schnetter <schnetter@cct.lsu.edu> | 2008-03-01 19:23:10 -0600 |
---|---|---|
committer | Erik Schnetter <schnetter@cct.lsu.edu> | 2008-03-01 19:23:10 -0600 |
commit | a1416c08b9c7dd1ce6ee710d022909caf9dd88e6 (patch) | |
tree | 2c0323f39906997e2bfad4973c3acf942fc720de /Carpet | |
parent | b2f0fee5802fdd1c796bfc3043b0086241cdb121 (diff) |
Carpet: Add new parameter Carpet::refine_timestep
Add a new parameter Carpet::refine_timestep which automatically adjusts
the time step size according to the give spatial and temporal refinement
factors, so that the specified CFL factor Time::dtfac is satisfied on all
refinement levels.
Diffstat (limited to 'Carpet')
-rw-r--r-- | Carpet/Carpet/param.ccl | 4 | ||||
-rw-r--r-- | Carpet/Carpet/schedule.ccl | 14 | ||||
-rw-r--r-- | Carpet/Carpet/src/CarpetBasegrid.cc | 36 | ||||
-rw-r--r-- | Carpet/Carpet/src/carpet.hh | 1 | ||||
-rw-r--r-- | Carpet/Carpet/src/make.code.defn | 1 |
5 files changed, 56 insertions, 0 deletions
diff --git a/Carpet/Carpet/param.ccl b/Carpet/Carpet/param.ccl index ab643e276..b8f2dbe93 100644 --- a/Carpet/Carpet/param.ccl +++ b/Carpet/Carpet/param.ccl @@ -132,6 +132,10 @@ STRING time_refinement_factors "Temporal refinement factors over the coarsest le "^[[:space:]]*\[[[:space:]]*[[:digit:]]+[[:space:]]*(,[[:space:]]*[[:digit:]]+[[:space:]]*)*\][[:space:]]*$" :: "[ <tfact>, ... ]" } "" +BOOLEAN refine_timestep "Correct Time::dtfac for spacings on finer grids" +{ +} "no" + CCTK_INT convergence_level "Convergence level" diff --git a/Carpet/Carpet/schedule.ccl b/Carpet/Carpet/schedule.ccl index b202a4a82..3c9322b5f 100644 --- a/Carpet/Carpet/schedule.ccl +++ b/Carpet/Carpet/schedule.ccl @@ -16,3 +16,17 @@ schedule CarpetParamCheck at PARAMCHECK { LANG: C } "Parameter checking routine" + + + +# Correct time step for finer levels when there are non-trivial +# time refinement factors + +if (refine_timestep) +{ + SCHEDULE CarpetRefineTimeStep AT basegrid AFTER Time_Simple + { + LANG: C + OPTIONS: singlemap + } "Correct time step size for spacing on finer grids" +} diff --git a/Carpet/Carpet/src/CarpetBasegrid.cc b/Carpet/Carpet/src/CarpetBasegrid.cc new file mode 100644 index 000000000..c2dc45f81 --- /dev/null +++ b/Carpet/Carpet/src/CarpetBasegrid.cc @@ -0,0 +1,36 @@ +#include <limits> + +#include "cctk.h" +#include "cctk_Arguments.h" +#include "cctk_Parameters.h" + +#include "carpet.hh" + + + +namespace Carpet { + + using namespace std; + + void CarpetRefineTimeStep (CCTK_ARGUMENTS) + { + DECLARE_CCTK_ARGUMENTS; + DECLARE_CCTK_PARAMETERS; + + // Find the smallest CFL factor for all refinement levels + CCTK_REAL min_cfl = numeric_limits<CCTK_REAL>::max(); + for (int rl = 0; rl < reflevels; ++ rl) { + CCTK_REAL const level_cfl = + timereffacts.AT(rl) / maxval (spacereffacts.AT(rl)); + min_cfl = min (min_cfl, level_cfl); + } + + // Reduce the time step correspondingly + cctkGH->cctk_delta_time *= min_cfl; + CCTK_VInfo (CCTK_THORNSTRING, + "Timestep reduced by a factor of %g to %g", + double (1 / min_cfl), + double (cctkGH->cctk_delta_time/cctkGH->cctk_timefac)); + } + +} // namespace Carpet diff --git a/Carpet/Carpet/src/carpet.hh b/Carpet/Carpet/src/carpet.hh index 271a86086..16dca4228 100644 --- a/Carpet/Carpet/src/carpet.hh +++ b/Carpet/Carpet/src/carpet.hh @@ -22,6 +22,7 @@ namespace Carpet { int CarpetStartup (void); int CarpetMultiModelStartup (void); void CarpetParamCheck (CCTK_ARGUMENTS); + void CarpetRefineTimeStep (CCTK_ARGUMENTS); } // Registered functions diff --git a/Carpet/Carpet/src/make.code.defn b/Carpet/Carpet/src/make.code.defn index 374e3dad4..1f70291cd 100644 --- a/Carpet/Carpet/src/make.code.defn +++ b/Carpet/Carpet/src/make.code.defn @@ -3,6 +3,7 @@ # Source files in this directory SRCS = AllGatherString.cc \ CallFunction.cc \ + CarpetBasegrid.cc \ CarpetParamCheck.cc \ CarpetStartup.cc \ Checksum.cc \ |