Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
DiscreteTime Class Reference

#include <deal.II/base/discrete_time.h>

Public Member Functions

 DiscreteTime (const double start_time, const double end_time, const double desired_start_step_size=0.)
 
double get_current_time () const
 
double get_next_time () const
 
double get_previous_time () const
 
double get_start_time () const
 
double get_end_time () const
 
bool is_at_start () const
 
bool is_at_end () const
 
double get_next_step_size () const
 
double get_previous_step_size () const
 
unsigned int get_step_number () const
 
void set_desired_next_step_size (const double time_step_size)
 
void set_next_step_size (const double time_step_size)
 
void advance_time ()
 
void restart ()
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Private Attributes

double start_time
 
double end_time
 
double current_time
 
double next_time
 
double previous_time
 
double start_step_size
 
unsigned int step_number
 

Detailed Description

This class provides a means to keep track of the simulation time in a time-dependent simulation. It manages stepping forward from a start time \(T_{\text{start}}\) to an end time \(T_{\text{end}}\). It also allows adjusting the time step size during the simulation. This class provides the necessary interface to be incorporated in any time-dependent simulation. The usage of this class is demonstrated in step-19 (and step-83) as well as step-21.

This class provides a number of invariants that are guaranteed to be true at all times.

  • The current simulation time is within the closed interval between the start time and the end time ( \(T_{\text{start}} \le t \le T_{\text{end}}\)).
  • Whenever time is incremented, the step size is positive ( \(dt > 0\)). In other words, time advances in strictly ascending order ( \(m < n \Leftrightarrow t_m < t_n\)).

The model this class follows is that one sets a desired time step length either through the constructor or using set_desired_next_step_size() function. This step size will then be used in all following calls to the advance_time() function, but may be adjusted slightly towards the end of the simulation to ensure that the simulation time hits the end time exactly. The adjustment is useful for the following reasons:

Let's say that you loop over all of the time steps by using a for loop

for (DiscreteTime time(0., 1., 0.3);
time.is_at_end() == false;
time.advance_time())
{
// Insert simulation code here
}
bool is_at_end() const

or, if you like this better, the equivalent while loop:

DiscreteTime time(0., 1., 0.3);
while (time.is_at_end() == false)
{
// Insert simulation code here
time.advance_time();
}

In the above example the time starts at \(T_{\text{start}} = 0\) until \(T_{\text{end}}=1\). Assuming the time step \(dt = 0.3\) is not modified inside the loop, the time is advanced from \(t = 0\) to \(t = 0.3\), \(t = 0.6\), \(t = 0.9\) and finally it reaches the end time at \(t = 1.0\). Here, the final step size needs to be reduced from its desired value of 0.3 to \(dt = 0.1\) in order to ensure that we finish the simulation exactly at the specified end time. In fact, you should assume that not only the last time step length may be adjusted, but also previously ones – for example, this class may take the liberty to spread the decrease in time step size out over several time steps and increment time from \(t=0\), to \(0.3\), \(0.6\), \(0.8\), and finally \(t=T_{\text{end}}=1\) to avoid too large a change in time step size from one step to another.

The other situation in which the time step needs to be adjusted (this time to slightly larger values) is if a time increment falls just short of the final time. Imagine, for example, a similar situation as above, but with different end time:

for (DiscreteTime time(0., 1.21, 0.3);
time.is_at_end() == false;
time.advance_time())
{
// Insert simulation code here
}

Here, the time step from \(t=0.9\) to \(t=1.2\) falls just short of the final time \(T_{\text{end}}=1.21\). Instead of following up with a very small step of length \(dt=0.01\), the class stretches the last time step (or last time steps) slightly to reach the desired end time.

The examples above make clear that the time step size given to this class is only a desired step size. You can query the actual time step size using the get_next_step_size() function.

Details of time-stepping

Since time is marched forward in a discrete manner in our simulations, we need to discuss how we increment time. During time stepping we enter two separate alternating regimes in every step.

  • The snapshot stage (the current stage, the consistent stage): In this part of the algorithm, we are at \(t = t_n\) and all quantities of the simulation (displacements, strains, temperatures, etc.) are up-to-date for \(t = t_n\). In this stage, current time refers to \(t_n\), next time refers to \(t_{n+1}\), previous time refers to \(t_{n-1}\). The other useful notation quantities are the next time step size \(t_{n+1} - t_n\) and previous time step size \(t_n - t_{n-1}\). In this stage, it is a perfect occasion to generate text output using print commands within the user's code. Additionally, post-processed outputs can be prepared here, which can then later be viewed by visualization programs such as Tecplot, Paraview, and VisIt. Additionally, during the snapshot stage, the code can assess the quality of the previous step and decide whether it wants to increase or decrease the time step size. The step size for the next time step can be modified here, by calling set_desired_next_step_size().
  • The update stage (the transition stage, the inconsistent stage): In this section of the program, the internal state of the simulation is getting updated from \(t_n\) to \(t_{n+1}\). All of the variables need to be updated one by one, the step number is incremented, the time is incremented by \(dt = t_{n+1} - t_n\), and time-integration algorithms are used to update the other simulation quantities. In the middle of this stage, some variables have been updated to \(t_{n+1}\) but other variables still represent their value at \(t_n\). Thus, we call this the inconsistent stage, requiring that no post-processing output related to the state variables take place within it. The state variables, namely those related to time, the solution field and any internal variables, are not synchronized and then get updated one by one. In general, the order of updating variables is arbitrary, but some care should be taken if there are interdependencies between them. For example, if some variable such as \(x\) depends on the calculation of another variable such as \(y\), then \(y\) must be updated before \(x\) can be updated.

    The question arises whether time should be incremented before updating state quantities. Multiple possibilities exist, depending on program and formulation requirements, and possibly the programmer's preferences:

    • Time is incremented before the rest of the updates. In this case, even though time is incremented to \(t_{n+1}\), not all variables are updated yet. During this update phase, \(dt\) equals the previous time step size. Previous means that it is referring to the \(dt\) of the advance_time() command that was performed previously. In the following example code, we are assuming that a and b are two state variables that need to be updated in this time step.
      time.advance_time();
      new_a = update_a(a, b, time.get_previous_step_size());
      b = update_b(a, b, time.get_previous_step_size());
      a = new_a;
      Here, the code starts in a consistent state, but once advance_time() is called, the time variable, a, and b are no longer consistent with each other until after the last statement. At that point, the variables are all consistent again.
    • Time is incremented from \(t_n\) to \(t_{n+1}\) after all variables have already been updated for \(t_{n+1}\). During the update stage, \(dt\) is denoted as the next time step size. Next means that \(dt\) of the step corresponds to the advance_time() command that will happen subsequently.
      new_a = update_a(a, b, time.get_next_step_size());
      b = update_b(a, b, time.get_next_step_size());
      a = new_a;
      time.advance_time();
    • Time is incremented in the middle of the other updates: In this case \(dt\) would correspond to next or previous depending of whether it is used before or after the call to advance_time().
      new_a = update_a(a, b, time.get_next_step_size());
      time.advance_time();
      b = update_b(a, b, time.get_previous_step_size());
      a = new_a;

One thing to note is that, during the update phase, \(dt\) is referred to either next or previous time step size, depending on whether advance_time() has been called yet. The notion of current time step size is ill-defined. In fact, in the update stage the definition of every variable depends on whether it has been updated yet or not, hence the name the inconsistent stage.

The following code snippet shows the code sections for the snapshot stage and the update stage in the context of a complete time-dependent simulation. This code follows the coding conventions incorporated in the tutorial examples. Note that even though this example is written in the format of a for loop, it can equivalently be written as a while or do while loop (as shown in step-21).

// pre-processing/setup stage {
make_grid();
setup_system();
for (DiscreteTime time(0., 1., 0.1); // } end pre-processing/setup stage
time.is_at_end() == false;
time.advance_time()) // part of the update stage, runs at
// the end of the loop body
{
// snapshot stage {
const double time_of_simulation = time.get_next_time();
const double timestep_size = time.get_next_step_size();
std::cout
<< "Timestep: " << time.get_step_number() << " -- "
<< "Solving for the solution at "
<< "t = " << time_of_simulation << " with "
<< "dt = " << timestep_size << '.' << std::endl;
// } end snapshot stage
// update stage {
assemble_system(time_of_simulation, timestep_size);
solve();
update_solutions();
// } end update stage
// snapshot stage {
output_results(time_of_solution);
// propose a new timestep size if need be
// time.set_desired_next_step_size(...);
// } end snapshot stage
}

The run() function in step-19 shows a very similar example where the call to advance_time() ends the update stage and is followed by generating graphical output with the then-current time.

Definition at line 232 of file discrete_time.h.

Constructor & Destructor Documentation

◆ DiscreteTime()

DiscreteTime::DiscreteTime ( const double start_time,
const double end_time,
const double desired_start_step_size = 0. )

Constructor.

Parameters
[in]start_timeThe time at the start of the simulation.
[in]end_timeThe time at the end of the simulation.
[in]desired_start_step_sizeA desired step size for incrementing time for the first step. It is not guaranteed that this value will be actually used as the size of the first step, as discussed in the introduction.
Precondition
desired_start_step_size must be non-negative.
Note
desired_start_step_size is an optional parameter. If it is not provided or it is specified as zero, it indicates that the desired size for the time step will be calculated at a different location in the code. In this case, the created object cannot increment time until the step size is changed by calling set_desired_next_step_size().

Definition at line 45 of file discrete_time.cc.

Member Function Documentation

◆ get_current_time()

double DiscreteTime::get_current_time ( ) const
inline

Return the current time.

Definition at line 515 of file discrete_time.h.

◆ get_next_time()

double DiscreteTime::get_next_time ( ) const
inline

Return the next time that we would reach if we were to advance the time by one step.

Note
If the simulation is at the end time, this method returns the end time.

Definition at line 523 of file discrete_time.h.

◆ get_previous_time()

double DiscreteTime::get_previous_time ( ) const
inline

Return the time we were at before advance_time() was called last time.

Note
If the simulation is at the start time, this method returns the start time.

Definition at line 531 of file discrete_time.h.

◆ get_start_time()

double DiscreteTime::get_start_time ( ) const
inline

Return the start time.

Definition at line 467 of file discrete_time.h.

◆ get_end_time()

double DiscreteTime::get_end_time ( ) const
inline

Return the end of the time interval. The final time step ends exactly at this point. This exact floating-point equality is very important because it allows us to equality-compare current time with end time and decide whether we have reached the end of the simulation.

Definition at line 475 of file discrete_time.h.

◆ is_at_start()

bool DiscreteTime::is_at_start ( ) const
inline

Return whether no step has taken place yet.

Definition at line 483 of file discrete_time.h.

◆ is_at_end()

bool DiscreteTime::is_at_end ( ) const
inline

Return whether time has reached the end time.

Definition at line 491 of file discrete_time.h.

◆ get_next_step_size()

double DiscreteTime::get_next_step_size ( ) const
inline

Return the size of the step from current time step to the next. As discussed in the introduction to the class, this is the actual time step, and may differ from the desired time step set in the constructor or through the set_desired_next_step_size() function.

Note
If the simulation is at the end time, this method returns zero.

Definition at line 499 of file discrete_time.h.

◆ get_previous_step_size()

double DiscreteTime::get_previous_step_size ( ) const
inline

Return the step size of the previous step.

Note
If the simulation is at the start time, this method returns zero.

Definition at line 507 of file discrete_time.h.

◆ get_step_number()

unsigned int DiscreteTime::get_step_number ( ) const
inline

Return the number of times the simulation time has been incremented. Return zero when the simulation is at the start time.

Definition at line 539 of file discrete_time.h.

◆ set_desired_next_step_size()

void DiscreteTime::set_desired_next_step_size ( const double time_step_size)

Set the desired value of the next time step size. By calling this method, we are indicating the next time advance_time() is called, we would like time_step_size to be used to advance the simulation time. However, if the step is too large such that the next simulation time exceeds the end time, the step size is truncated. Additionally, if the step size is such that the next simulation time approximates the end time (but falls just slightly short of it), the step size is adjusted such that the next simulation time exactly matches the end time.

Definition at line 62 of file discrete_time.cc.

◆ set_next_step_size()

void DiscreteTime::set_next_step_size ( const double time_step_size)

Set the actual value of the next time step size. By calling this method, we are indicating the next time advance_time() is called, time_step_size is to be used to advance the simulation time.

Note
The difference between set_next_step_size() and set_desired_next_step_size() is that the former uses the provided \(dt\) exactly without any adjustment, but produces an error (in debug mode) if \(dt\) is not in the acceptable range. Generally, set_desired_next_step_size() is the preferred method because it can adjust the \(dt\) intelligently, based on \(T_{\text{end}}\).
Precondition
\(0 < dt \le T_{\text{end}} - t\).

Definition at line 70 of file discrete_time.cc.

◆ advance_time()

void DiscreteTime::advance_time ( )

Advance the current time based on the value of the current step. If you want to adjust the next time step size, call the method set_desired_next_step_size() before calling this method. If you call this function repeatedly, the time is increased with the same step size until it reaches the end time. See the documentation of set_desired_next_step_size() for explanation of the rules for automatic adjustment of the step size.

Precondition
Current time must be smaller than the end time. The object cannot advance time if it is already at the end time. This rule is created to avoid the creation of an infinite loop when advance_time() is called inside a loop.
The time step size must be nonzero. If the step size is currently zero, change it by calling set_desired_next_step_size() before calling advance_time().

Definition at line 84 of file discrete_time.cc.

◆ restart()

void DiscreteTime::restart ( )

Set the current time equal to start time and set the step size to the initial step size.

Definition at line 100 of file discrete_time.cc.

◆ memory_consumption()

std::size_t DiscreteTime::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 111 of file discrete_time.cc.

◆ serialize()

template<class Archive >
void DiscreteTime::serialize ( Archive & ar,
const unsigned int version )
inline

Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function is used in step-83.

Definition at line 548 of file discrete_time.h.

Member Data Documentation

◆ start_time

double DiscreteTime::start_time
private

The beginning of the time interval.

Definition at line 419 of file discrete_time.h.

◆ end_time

double DiscreteTime::end_time
private

The end of the time interval.

Definition at line 424 of file discrete_time.h.

◆ current_time

double DiscreteTime::current_time
private

The current time.

Definition at line 429 of file discrete_time.h.

◆ next_time

double DiscreteTime::next_time
private

The time at the next step.

Note
Internally, the next simulation time is stored instead of the current step size. For example, when the method set_desired_next_step_size() is called, it computes the appropriate next simulation time and stores it. When advance_time() is called, the current_time is replaced by next_time. This choice for the internal state allows for simpler code and ensures than when we call advance_time() at the last step, the floating-point value of the time exactly matches the end time.

Definition at line 443 of file discrete_time.h.

◆ previous_time

double DiscreteTime::previous_time
private

The previous time.

Definition at line 448 of file discrete_time.h.

◆ start_step_size

double DiscreteTime::start_step_size
private

The size of the first step.

Definition at line 453 of file discrete_time.h.

◆ step_number

unsigned int DiscreteTime::step_number
private

The step number i.e. the number of times the simulation time ha been incremented.

Definition at line 459 of file discrete_time.h.


The documentation for this class was generated from the following files: