Reference documentation for deal.II version GIT relicensing-422-gb369f187d8 2024-04-17 18:10:02+00:00
\(\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
Public Member Functions | Private Attributes | List of all members
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 ()
 

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-21.

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

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.

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 231 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 44 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 496 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 504 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 512 of file discrete_time.h.

◆ get_start_time()

double DiscreteTime::get_start_time ( ) const
inline

Return the start time.

Definition at line 448 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 456 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 464 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 472 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 480 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 488 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 520 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 61 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 69 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 83 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 99 of file discrete_time.cc.

Member Data Documentation

◆ start_time

double DiscreteTime::start_time
private

The beginning of the time interval.

Definition at line 400 of file discrete_time.h.

◆ end_time

double DiscreteTime::end_time
private

The end of the time interval.

Definition at line 405 of file discrete_time.h.

◆ current_time

double DiscreteTime::current_time
private

The current time.

Definition at line 410 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 424 of file discrete_time.h.

◆ previous_time

double DiscreteTime::previous_time
private

The previous time.

Definition at line 429 of file discrete_time.h.

◆ start_step_size

double DiscreteTime::start_step_size
private

The size of the first step.

Definition at line 434 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 440 of file discrete_time.h.


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