aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--doc/documention.tex72
-rw-r--r--interface.ccl2
-rw-r--r--param.ccl15
-rw-r--r--schedule.ccl4
-rw-r--r--src/Courant.c8
5 files changed, 55 insertions, 46 deletions
diff --git a/doc/documention.tex b/doc/documention.tex
index 5629beb..44038b7 100644
--- a/doc/documention.tex
+++ b/doc/documention.tex
@@ -13,51 +13,61 @@
This thorn provides a routines for calculating
the timestep for an evolution based on the spatial Cartesian grid spacing and
a wave speed.
-\section{Comments}
-
-There are currently two methods for calculating the timestep, either the
-simple courant method, or a dynamic courant method. The method is chosen
-using the keyword parameter {\tt time::method}
-The timestep, is passed into the Cactus variable {\tt cctk\_delta\_time}.
-Both timesteps are based on the Courant condition, which states that for
-numerical stability, the chosen timestep should satisfy
-$$
-\Delta t \le \mbox{min}(\Delta x^i)/\mbox{wave speed}/\sqrt(\mbox{dim})
-$$
-The two methods currently implemented are:
+\section{Description}
+
+Thorn {\tt Time} uses one of four methods to decide on the timestep
+to be used for the simulation. The method is chosen using the
+keyword parameter {\tt time::timestep\_method}. (Note: In releases Beta 8 and
+earlier the parameter used was {\tt time::courant\_method}
\begin{itemize}
-\item {\tt time::method = "standard"} (this is the default). The timestep is calculated once
- at the start of the run in the {\tt BASEGRID} timebin, and is then held static
- throughout the run. The algorithm, which uses the parameter {\tt time::dtfac} is
+\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
+ 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}.
$$
\Delta t = \mbox{\tt dtfac} * \mbox{min}(\Delta x^i)
$$
- Note that the parameter {\tt dtfac} should take into account the dimension
- of the space being used, and the wave speed.
-
-\item {\tt time::method = "courant"}. The timestep is calculated dynamically at the
- start of each iteration in the {\tt PRESTEP} timebin. The algorithm is
+ 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
+ dynamic courant type condition, the timestep being set before each
+ timestep 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{maximum wavespeed}/\sqrt(\mbox{dim})
+\Delta t = \mbox{\tt courant\_fac} * \mbox{min}(\Delta x^i)/\mbox{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.
-
-\item {\tt time::method = "courant\_time"}. The timestep is calculated dynamically at the
- start of each iteration in the {\tt PRESTEP} timebin. The algorithm is
+\item{} {\tt time::timestep\_method = ``courant\_time''} 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{minimum time to cross a zone}/\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.
+
+\end{itemize}
-\section{Technical Details}
+In all cases, Thorn {\tt Time} sets the Cactus variable {\tt cctk\_delta\_time}
+which is passed as part of the macro {\tt CCTK\_ARGUMENTS} to thorns called
+by the scheduler.
-If a dynamic {\tt courant} condition is selected, a thorn must set the protected variable
-{\tt courant\_wave\_speed} for the maximum wave speed before Time sets the timestep.
+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})
+$$
-If a dynamic {\tt courant\_time} condition is selected, a thorn must set the protected variable
-{\tt courant\_time} for the minimum time for a wave to cross a zone
- before Time sets the timestep.
\end{itemize}
diff --git a/interface.ccl b/interface.ccl
index db7f042..6ea9cdc 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -8,7 +8,7 @@ protected:
REAL speedvars type=SCALAR
{
courant_wave_speed
- courant_time
+ courant_min_time
} "Speed to use for Courant condition"
private:
diff --git a/param.ccl b/param.ccl
index 79ce823..30fa3ff 100644
--- a/param.ccl
+++ b/param.ccl
@@ -3,19 +3,18 @@
restricted:
-KEYWORD courant_method "Method for calculating timestep"
+KEYWORD timestep_method "Method for calculating timestep"
{
- "none" :: "Use given timestep"
- "standard" :: "Courant condition at BASEGRID"
- "courant" :: "Courant condition at PRESTEP (using wavespeed)"
+ "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)"
-} "standard"
+} "courant_static"
-BOOLEAN courant_outonly "Only output courant timestep?"
+BOOLEAN timestep_outonly "Don't set a dynamic timestep, just output what it would be"
{
} "no"
-
private:
REAL timestep "Absolute value for timestep"
@@ -33,7 +32,7 @@ REAL courant_fac "The courant timestep condition dt = courant_fac*max(delta_spac
0:* :: "Probably only makes sense to be bigger than zero"
} 0.9
-INT outcourant_every "How often to output courant timestep"
+INT outtimestep_every "How often to output courant timestep"
{
0:* :: "Zero means no output"
} 0 \ No newline at end of file
diff --git a/schedule.ccl b/schedule.ccl
index 2399f7f..343bede 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -1,14 +1,14 @@
# Schedule definitions for thorn Time
# $Header$
-if (CCTK_Equals(courant_method,"standard"))
+if (CCTK_Equals(timestep_method,"courant_static"))
{
schedule Time_Simple at CCTK_BASEGRID after SpatialCoordinates
{
LANG: C
} "Set timestep based on Courant condition"
}
-else if (CCTK_Equals(courant_method,"courant") || CCTK_Equals(courant_method,"courant_time"))
+else if (CCTK_Equals(timestep_method,"courant_speed") || CCTK_Equals(timestep_method,"courant_time"))
{
schedule Time_Simple at CCTK_BASEGRID after SpatialCoordinates
{
diff --git a/src/Courant.c b/src/Courant.c
index 81cfa60..eb778d3 100644
--- a/src/Courant.c
+++ b/src/Courant.c
@@ -47,18 +47,18 @@ void Time_Courant(CCTK_ARGUMENTS)
}
/* Calculate the courant timestep */
- if (CCTK_Equals(courant_method,"courant_time"))
+ if (CCTK_Equals(timestep_method,"courant_time"))
{
- *courant_dt = courant_fac*(*courant_time)/sqrt((CCTK_REAL )cctk_dim);
+ *courant_dt = courant_fac*(*courant_min_time)/sqrt((CCTK_REAL )cctk_dim);
}
- else if (CCTK_Equals(courant_method,"courant"))
+ else if (CCTK_Equals(timestep_method,"courant_speed"))
{
*courant_dt = courant_fac/(*courant_wave_speed)/sqrt((CCTK_REAL )cctk_dim);
}
/* Set the Cactus timestep */
- if (! courant_outonly)
+ if (! timestep_outonly)
{
cctkGH->cctk_delta_time = *courant_dt;