ryujin 2.1.1 revision 0eab90fbc6e1ac9f2e0a2e6d16f9f023c13a02f7
List of all members
ryujin::TimeIntegrator< Description, dim, Number > Class Template Referencefinal

#include <source/time_integrator.h>

Inheritance diagram for ryujin::TimeIntegrator< Description, dim, Number >:
Inheritance graph
[legend]
Collaboration diagram for ryujin::TimeIntegrator< Description, dim, Number >:
Collaboration graph
[legend]

Public Types

Typedefs and constexpr constants
using HyperbolicSystem = typename Description::HyperbolicSystem
 
using View = typename Description::template HyperbolicSystemView< dim, Number >
 
using ParabolicSystem = typename Description::ParabolicSystem
 
using StateVector = typename View::StateVector
 

Public Member Functions

Constructor and setup
 TimeIntegrator (const MPIEnsemble &mpi_ensemble, const OfflineData< dim, Number > &offline_data, const HyperbolicModule< Description, dim, Number > &hyperbolic_module, const ParabolicModule< Description, dim, Number > &parabolic_module, const std::string &subsection="/TimeIntegrator")
 
void prepare ()
 

Functions for performing explicit time steps

Number step (StateVector &state_vector, Number t, Number t_final=std::numeric_limits< Number >::max())
 
const auto & time_stepping_scheme () const
 
const auto & efficiency () const
 
Number step_ssprk_22 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_ssprk_33 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_erk_11 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_erk_22 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_erk_33 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_erk_43 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_erk_54 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_strang_ssprk_33_cn (StateVector &state_vector, Number t, Number tau_max)
 
Number step_strang_erk_33_cn (StateVector &state_vector, Number t, Number tau_max)
 
Number step_strang_erk_43_cn (StateVector &state_vector, Number t, Number tau_max)
 
Number step_imex_11 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_imex_22 (StateVector &state_vector, Number t, Number tau_max)
 
Number step_imex_33 (StateVector &state_vector, Number t, Number tau_max)
 

Detailed Description

template<typename Description, int dim, typename Number = double>
class ryujin::TimeIntegrator< Description, dim, Number >

The TimeIntegrator class implements IMEX timestepping strategies based on explicit and diagonally-implicit Runge Kutta schemes.

Definition at line 252 of file time_integrator.h.

Member Typedef Documentation

◆ HyperbolicSystem

template<typename Description , int dim, typename Number = double>
using ryujin::TimeIntegrator< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem

Definition at line 260 of file time_integrator.h.

◆ View

template<typename Description , int dim, typename Number = double>
using ryujin::TimeIntegrator< Description, dim, Number >::View = typename Description::template HyperbolicSystemView<dim, Number>

Definition at line 262 of file time_integrator.h.

◆ ParabolicSystem

template<typename Description , int dim, typename Number = double>
using ryujin::TimeIntegrator< Description, dim, Number >::ParabolicSystem = typename Description::ParabolicSystem

Definition at line 265 of file time_integrator.h.

◆ StateVector

template<typename Description , int dim, typename Number = double>
using ryujin::TimeIntegrator< Description, dim, Number >::StateVector = typename View::StateVector

Definition at line 267 of file time_integrator.h.

Constructor & Destructor Documentation

◆ TimeIntegrator()

template<typename Description , int dim, typename Number >
ryujin::TimeIntegrator< Description, dim, Number >::TimeIntegrator ( const MPIEnsemble mpi_ensemble,
const OfflineData< dim, Number > &  offline_data,
const HyperbolicModule< Description, dim, Number > &  hyperbolic_module,
const ParabolicModule< Description, dim, Number > &  parabolic_module,
const std::string &  subsection = "/TimeIntegrator< Description, dim, Number >" 
)

Constructor.

Definition at line 29 of file time_integrator.template.h.

References ryujin::cruise_control, ryujin::erk_33, and ryujin::strang_erk_33_cn.

Member Function Documentation

◆ prepare()

template<typename Description , int dim, typename Number >
void ryujin::TimeIntegrator< Description, dim, Number >::prepare

Prepare time integration. A call to prepare() allocates temporary storage and is necessary before any of the following time-stepping functions can be called.

Definition at line 82 of file time_integrator.template.h.

References ryujin::erk_11, ryujin::erk_22, ryujin::erk_33, ryujin::erk_43, ryujin::erk_54, ryujin::imex_11, ryujin::imex_22, ryujin::imex_33, ryujin::ssprk_22, ryujin::ssprk_33, ryujin::strang_erk_33_cn, ryujin::strang_erk_43_cn, and ryujin::strang_ssprk_33_cn.

◆ step()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step ( StateVector state_vector,
Number  t,
Number  t_final = std::numeric_limits<Number>::max() 
)

Given a reference to a previous state vector U performs an explicit time step (and store the result in U). The function returns the chosen time step size tau. The time step size tau is selected such that $t + tau <= t_final$.

Note
This function switches between different Runge-Kutta methods depending on chosen runtime parameters.
Depending on chosen run time parameters different CFL adaptation and recovery strategies for invariant domain violations are used.

Definition at line 227 of file time_integrator.template.h.

References ryujin::erk_11, ryujin::erk_22, ryujin::erk_33, ryujin::erk_43, ryujin::erk_54, ryujin::imex_11, ryujin::imex_22, ryujin::imex_33, ryujin::none, ryujin::raise_exception, ryujin::ssprk_22, ryujin::ssprk_33, ryujin::strang_erk_33_cn, ryujin::strang_erk_43_cn, and ryujin::strang_ssprk_33_cn.

◆ time_stepping_scheme()

template<typename Description , int dim, typename Number = double>
const auto & ryujin::TimeIntegrator< Description, dim, Number >::time_stepping_scheme ( ) const
inline

The selected time-stepping scheme.

Definition at line 318 of file time_integrator.h.

◆ efficiency()

template<typename Description , int dim, typename Number = double>
const auto & ryujin::TimeIntegrator< Description, dim, Number >::efficiency ( ) const
inline

The eficiency of the selected time-stepping scheme expressed as the ratio of step size of the combined method to step size of an elementary forward Euler step. For example, SSPRK33 has an efficiency ratio of 1 whereas ERK33 has an efficiency ratio of 3.

Definition at line 326 of file time_integrator.h.

◆ step_ssprk_22()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_ssprk_22 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit second-order strong-stability preserving Runge-Kutta SSPRK(2,2;1/2) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 335 of file time_integrator.template.h.

◆ step_ssprk_33()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_ssprk_33 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit third-order strong-stability preserving Runge-Kutta SSPRK(3,3;1/3) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 364 of file time_integrator.template.h.

◆ step_erk_11()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_erk_11 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit first-order Euler step ERK(1,1;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 400 of file time_integrator.template.h.

◆ step_erk_22()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_erk_22 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit second-order Runge-Kutta ERK(2,2;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 420 of file time_integrator.template.h.

◆ step_erk_33()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_erk_33 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit third-order Runge-Kutta ERK(3,3;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 445 of file time_integrator.template.h.

◆ step_erk_43()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_erk_43 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit 4 stage third-order Runge-Kutta ERK(4,3;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 481 of file time_integrator.template.h.

◆ step_erk_54()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_erk_54 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an explicit 4 stage fourth-order Runge-Kutta ERK(5,4;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 522 of file time_integrator.template.h.

◆ step_strang_ssprk_33_cn()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_strang_ssprk_33_cn ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs a combined explicit implicit Strang split using a third-order Runge-Kutta ERK(3,3;1/3) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 594 of file time_integrator.template.h.

References ryujin::sadd(), and ryujin::Restart::suggested_tau_max.

◆ step_strang_erk_33_cn()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_strang_erk_33_cn ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs a combined explicit implicit Strang split using a third-order Runge-Kutta ERK(3,3;1) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 651 of file time_integrator.template.h.

References ryujin::Restart::suggested_tau_max.

◆ step_strang_erk_43_cn()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_strang_erk_43_cn ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs a combined explicit implicit Strang split using a third-order Runge-Kutta ERK(4,3;1) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max.

Definition at line 715 of file time_integrator.template.h.

References ryujin::Restart::suggested_tau_max.

◆ step_imex_11()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_imex_11 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U, performs an implicit-explicit step IMEX(1,1;1) using a forward euler scheme for the hyperbolic subproblem and backward euler scheme for the parabolic subproblem.

Definition at line 787 of file time_integrator.template.h.

◆ step_imex_22()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_imex_22 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an implicit-explicit IMEX(2,2;1) step using a two stage midpoint rule for the hyperbolic subproblem and a two stage midpoint rule for the parabolic subproblem. The function returns the chosen time step size tau.

Definition at line 813 of file time_integrator.template.h.

◆ step_imex_33()

template<typename Description , int dim, typename Number >
Number ryujin::TimeIntegrator< Description, dim, Number >::step_imex_33 ( StateVector state_vector,
Number  t,
Number  tau_max 
)
protected

Given a reference to a previous state vector U performs an implicit-explicit IMEX(3,3;1) step using a three stage ERK tableau for the hyperbolic subproblem and a three stage DIRK tableau for the parabolic subproblem. The function returns the chosen time step size tau.

Definition at line 852 of file time_integrator.template.h.


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