diff options
author | yye00 <yye00@5633253d-7678-4964-a54d-f87795f8ee59> | 2004-10-05 20:01:05 +0000 |
---|---|---|
committer | yye00 <yye00@5633253d-7678-4964-a54d-f87795f8ee59> | 2004-10-05 20:01:05 +0000 |
commit | fa567ab69ac4eaa70d703b573d5486267e649039 (patch) | |
tree | 2b43f9694582304ddd1faa5b5706506f5bb4c18a | |
parent | f81ea75b2b50d8f644106ef1110101a7b71ae3af (diff) |
Call CFL setting routines only once per refinement level.
Calculate the grid spacing correctly even for dim>3.
Correctly cast the values before passing them to printf.
Print the newly calculated time step size, not the old time step size.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/Time/trunk@62 5633253d-7678-4964-a54d-f87795f8ee59
-rw-r--r-- | schedule.ccl | 2 | ||||
-rw-r--r-- | src/Simple.c | 43 |
2 files changed, 15 insertions, 30 deletions
diff --git a/schedule.ccl b/schedule.ccl index 12133ce..d70c45a 100644 --- a/schedule.ccl +++ b/schedule.ccl @@ -6,6 +6,7 @@ STORAGE: speedvars, couranttemps schedule Time_Initialise at CCTK_BASEGRID before (Time_Simple, Time_Given) { LANG: C + OPTIONS: level } "Initialise Time variables" if (CCTK_Equals (timestep_method, "courant_static")) @@ -13,6 +14,7 @@ if (CCTK_Equals (timestep_method, "courant_static")) schedule Time_Simple at CCTK_BASEGRID after SpatialCoordinates { LANG: C + OPTIONS: level } "Set timestep based on Courant condition (courant_static)" } else if (CCTK_Equals (timestep_method, "courant_speed")) diff --git a/src/Simple.c b/src/Simple.c index e809b7f..15ac2cf 100644 --- a/src/Simple.c +++ b/src/Simple.c @@ -21,32 +21,18 @@ void Time_Simple(CCTK_ARGUMENTS); void Time_Simple(CCTK_ARGUMENTS) { - DECLARE_CCTK_PARAMETERS - DECLARE_CCTK_ARGUMENTS + DECLARE_CCTK_ARGUMENTS; + DECLARE_CCTK_PARAMETERS; - CCTK_REAL min_spacing=0; + CCTK_REAL min_spacing; + int d; /* Calculate the minimum grid spacing */ - if (cctk_dim>=1) + min_spacing = cctk_delta_space[0]; + for (d=1; d<cctk_dim; ++d) { - min_spacing = cctk_delta_space[0]; - } - - if (cctk_dim>=2) - { - min_spacing = (min_spacing<cctk_delta_space[1] ? - min_spacing : cctk_delta_space[1]); - } - - if (cctk_dim>=3) - { - min_spacing = (min_spacing<cctk_delta_space[2] ? - min_spacing : cctk_delta_space[2]); - } - - if (cctk_dim>=4) - { - CCTK_WARN(0,"Time Step not defined for greater than 4 dimensions"); + min_spacing = (min_spacing<cctk_delta_space[d] + ? min_spacing : cctk_delta_space[d]); } /* Calculate the timestep */ @@ -55,17 +41,14 @@ void Time_Simple(CCTK_ARGUMENTS) if (verbose) { CCTK_VInfo(CCTK_THORNSTRING, - "Using simple courant condition to set timestep"); + "Using a simple Courant condition to set then timestep"); CCTK_VInfo(CCTK_THORNSTRING, - " ... using dtfac of %f",dtfac); + " ... using a dtfac of %g", (double)dtfac); CCTK_VInfo(CCTK_THORNSTRING, - " ... using minimum spacing of %f",min_spacing); + " ... using a minimum spacing of %g", (double)min_spacing); } CCTK_VInfo(CCTK_THORNSTRING, - "Timestep set to %g (courant_static)",CCTK_DELTA_TIME); - + "Timestep set to %g (courant_static)", + (double)(cctkGH->cctk_delta_time/cctkGH->cctk_timefac)); } - - - |