aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorallen <allen@5633253d-7678-4964-a54d-f87795f8ee59>2003-06-07 17:21:18 +0000
committerallen <allen@5633253d-7678-4964-a54d-f87795f8ee59>2003-06-07 17:21:18 +0000
commit93c0743acfc174d1455e6325deefd513fb4e5c6d (patch)
tree335cd38fcc3c7681a8b1566b8b2e8312946ac3e7
parent10701d221b2e4ffc413c7ec50bdc35aa5efab4c1 (diff)
Added different bits of documentation, added code for a parameter
which wasn't being used (and changed its name but since it wasn't implemented hopefully noone was using outtimestep_every). Removed two scheduled routines which were not doing anything. Fixes Cactus/1377, Cactus/1378, Cactus/1384, Cactus/1379. Does not fix Cactus/1385 yet. git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/Time/trunk@55 5633253d-7678-4964-a54d-f87795f8ee59
-rw-r--r--doc/documentation.tex43
-rw-r--r--param.ccl14
-rw-r--r--schedule.ccl25
-rw-r--r--src/Courant.c7
-rw-r--r--src/Simple.c2
5 files changed, 60 insertions, 31 deletions
diff --git a/doc/documentation.tex b/doc/documentation.tex
index 69a6593..8e9e4b3 100644
--- a/doc/documentation.tex
+++ b/doc/documentation.tex
@@ -35,10 +35,14 @@ keyword parameter {\tt time::timestep\_method}.
\begin{itemize}
-\item{} {\tt time::timestep\_method = ``given''} The timestep is fixed to the
+\item{} {\tt time::timestep\_method = "given"}
+
+ The timestep is fixed to the
value of the parameter {\tt time::timestep}.
-\item{} {\tt time::timestep\_method = ``courant\_static''} This is the default
+\item{} {\tt time::timestep\_method = "courant\_static"}
+
+ This is the default
method, which calculates the timestep once at the start of the
simulation, based on a simple courant type condition using
the spatial gridsizes and the parameter {\tt time::dtfac}.
@@ -48,25 +52,37 @@ $$
Note that it is up to the user to custom {\tt dtfac} to take
into account the dimension of the space being used, and the wave speed.
-\item{} {\tt time::timestep\_method = ``courant\_speed''} This choice implements a
+\item{} {\tt time::timestep\_method = "courant\_speed"}
+
+ This choice implements a
dynamic courant type condition, the timestep being set before each
iteration using the spatial dimension of the grid, the spatial grid sizes, the
parameter {\tt courant\_fac} and the grid variable {\tt courant\_wave\_speed}.
The algorithm used is
$$
-\Delta t = \mbox{\tt courant\_fac} * \mbox{min}(\Delta x^i)/\mbox{courant\_wave\_speed}/\sqrt(\mbox{dim})
+\Delta t = \mbox{\tt courant\_fac} * \mbox{min}(\Delta x^i)/\mbox{\tt courant\_wave\_speed}/\sqrt{\mbox dim}
$$
For this algorithm to be successful, the variable {\tt courant\_wave\_speed}
- must have been set by a thorn to the maximum wave speed on the grid.
+ must have been set by some thorn to the maximum propagation speed on the grid {\it before}
+ this thorn sets the timestep, that is {\tt AT POSTSTEP BEFORE Time\_Courant} (or earlier
+ in the evolution loop). [Note: The name {\tt courant\_wave\_speed} was poorly
+ chosen here, the required speed is the maximum propagation speed on
+ the grid which may be larger than the maximum wave speed (for example
+ with a shock wave in hydrodynamics, also it is possible to have
+ propagation without waves as with a pure advection equation).
+
+\item{} {\tt time::timestep\_method = "courant\_time"}
-\item{} {\tt time::timestep\_method = ``courant\_time''} This choice is similar to the
+ This choice is similar to the
method {\tt courant\_speed} above, in implementing a dynamic timestep.
However the timestep is chosen using
$$
-\Delta t = \mbox{\tt courant\_fac} * \mbox{\tt courant\_min\_time}/\sqrt(\mbox{dim})
+\Delta t = \mbox{\tt courant\_fac} * \mbox{\tt courant\_min\_time}/\sqrt{\mbox dim}
$$
- where the grid variable {\tt courant\_min\_time} must be set by a thorn to
- the minimum time for a wave to cross a gridzone.
+ where the grid variable {\tt courant\_min\_time} must be set by some thorn to
+ the minimum time for a wave to cross a gridzone {\it before}
+ this thorn sets the timestep, that is {\tt AT POSTSTEP BEFORE Time\_Courant} (or earlier
+ in the evolution loop).
\end{itemize}
@@ -78,23 +94,26 @@ Note that for hyperbolic problems, the Courant condition gives a minimum
requirement for stability, namely that the numerical domain of dependency
must encompass the physical domain of dependency, or
$$
-\Delta t \le \mbox{min}(\Delta x^i)/\mbox{wave speed}/\sqrt(\mbox{dim})
+\Delta t \le \mbox{min}(\Delta x^i)/\mbox{wave speed}/\sqrt{\mbox dim}
$$
\section{Examples}
+\noindent
{\bf Fixed Value Timestep}
{\tt
\begin{verbatim}
-time::timestep_method = ``given''
+time::timestep_method = "given"
time::timestep = 0.1
\end{verbatim}
}
+\noindent
{\bf Calculate Static Timestep Based on Grid Spacings}
+\noindent
The following parameters set the timestep to be 0.25
{\tt
@@ -102,7 +121,7 @@ The following parameters set the timestep to be 0.25
grid::dx = 0.5
grid::dy = 1.0
grid::dz = 1.0
-time::timestep_method = ``courant_static''
+time::timestep_method = "courant_static"
time::dtfac = 0.5
\end{verbatim}
}
diff --git a/param.ccl b/param.ccl
index 2f9e89f..d9c1616 100644
--- a/param.ccl
+++ b/param.ccl
@@ -11,9 +11,9 @@ restricted:
KEYWORD timestep_method "Method for calculating timestep"
{
"given" :: "Use given timestep"
- "courant_static" :: "Courant condition at BASEGRID"
- "courant_speed" :: "Courant condition at PRESTEP (using wavespeed)"
- "courant_time" :: "Courant condition at PRESTEP (using min time)"
+ "courant_static" :: "Courant condition at BASEGRID (using dtfac)"
+ "courant_speed" :: "Courant condition at POSTSTEP (using wavespeed and courant_fac)"
+ "courant_time" :: "Courant condition at POSTSTEP (using min time and courant_fac)"
} "courant_static"
BOOLEAN timestep_outonly "Don't set a dynamic timestep, just output what it would be"
@@ -39,11 +39,11 @@ REAL courant_fac "The courant timestep condition dt = courant_fac*max(delta_spac
*:0 :: "For negative timestep"
} 0.9
-INT outtimestep_every "How often to output courant timestep"
+INT timestep_outevery "How often to output courant timestep"
{
- 0:* :: "Zero means no output"
-} 0
+ 1:* :: "Zero means no output"
+} 1
-BOOLEAN verbose "Give information about timestep setting"
+BOOLEAN verbose "Give selective information about timestep setting"
{
} "no" \ No newline at end of file
diff --git a/schedule.ccl b/schedule.ccl
index 927b5a1..12133ce 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -13,21 +13,28 @@ if (CCTK_Equals (timestep_method, "courant_static"))
schedule Time_Simple at CCTK_BASEGRID after SpatialCoordinates
{
LANG: C
- } "Set timestep based on Courant condition"
+ } "Set timestep based on Courant condition (courant_static)"
}
-else if (CCTK_Equals (timestep_method, "courant_speed") ||
- CCTK_Equals (timestep_method, "courant_time"))
+else if (CCTK_Equals (timestep_method, "courant_speed"))
{
schedule Time_Simple at CCTK_BASEGRID after SpatialCoordinates
{
LANG: C
- } "Set timestep based on Courant condition"
+ } "Set timestep based on Courant condition (courant_speed)"
+
+ schedule Time_Courant at CCTK_POSTSTEP
+ {
+ LANG: C
+ } "Reset timestep each iteration"
+}
+else if (CCTK_Equals (timestep_method, "courant_time"))
+
+{
+ schedule Time_Simple at CCTK_BASEGRID after SpatialCoordinates
+ {
+ LANG: C
+ } "Set timestep based on Courant condition (courant_time)"
- schedule Time_Courant at CCTK_POSTINITIAL
- {
- LANG: C
- } "Reset timestep each iteration"
-
schedule Time_Courant at CCTK_POSTSTEP
{
LANG: C
diff --git a/src/Courant.c b/src/Courant.c
index d40fb84..ba6587e 100644
--- a/src/Courant.c
+++ b/src/Courant.c
@@ -74,12 +74,15 @@ void Time_Courant(CCTK_ARGUMENTS)
cctkGH->cctk_delta_time = *courant_dt;
if (verbose)
{
- CCTK_VInfo(CCTK_THORNSTRING,"Time step set to %f",cctkGH->cctk_delta_time);
+ CCTK_VInfo(CCTK_THORNSTRING,"Time step set to %g",cctkGH->cctk_delta_time);
}
}
else
{
cctkGH->cctk_delta_time = dtfac*min_spacing;
- CCTK_VInfo(CCTK_THORNSTRING,"Courant timestep would be %f",*courant_dt);
+ if (cctkGH->cctk_iteration % timestep_outevery == 0)
+ {
+ CCTK_VInfo(CCTK_THORNSTRING,"Courant timestep would be %g",*courant_dt);
+ }
}
}
diff --git a/src/Simple.c b/src/Simple.c
index 39c049c..af22128 100644
--- a/src/Simple.c
+++ b/src/Simple.c
@@ -63,7 +63,7 @@ void Time_Simple(CCTK_ARGUMENTS)
}
CCTK_VInfo(CCTK_THORNSTRING,
- "Timestep set to %f (Simple Courant)",cctkGH->cctk_delta_time);
+ "Timestep set to %g (courant_static)",cctkGH->cctk_delta_time);
}