Discontinuous Galerkin Library
#include "dg/algorithm.h"
dg::DIRKStep< ContainerType > Struct Template Reference

Embedded diagonally implicit Runge Kutta time-step with error estimate \( \begin{align} k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i} a_{ij} k_j\right) \\ u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j \\ \tilde u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s \tilde b_j k_j \\ \delta^{n+1} = u^{n+1} - \tilde u^{n+1} = \Delta t\sum_{j=1}^s (b_j-\tilde b_j) k_j \end{align} \). More...

Public Types

using value_type = get_value_type< ContainerType >
 the value type of the time variable (float or double) More...
 
using container_type = ContainerType
 

Public Member Functions

 DIRKStep ()
 No memory allocation. More...
 
 DIRKStep (ConvertsToButcherTableau< value_type > im_tableau, const ContainerType &copyable)
 Construct with a diagonally implicit Butcher Tableau. More...
 
template<class ... Params>
void construct (Params &&...ps)
 Perfect forward parameters to one of the constructors. More...
 
const ContainerType & copyable () const
 Return an object of same size as the object used for construction. More...
 
template<class ImplicitRHS , class Solver >
void step (const std::tuple< ImplicitRHS, Solver > &ode, value_type t0, const ContainerType &u0, value_type &t1, ContainerType &u1, value_type dt, ContainerType &delta)
 Advance one step with error estimate. More...
 
template<class ImplicitRHS , class Solver >
void step (const std::tuple< ImplicitRHS, Solver > &ode, value_type t0, const ContainerType &u0, value_type &t1, ContainerType &u1, value_type dt)
 Advance one step ignoring error estimate and embedded method. More...
 
unsigned order () const
 global order of the method given by the current Butcher Tableau More...
 
unsigned embedded_order () const
 global order of the embedding given by the current Butcher Tableau More...
 
unsigned num_stages () const
 number of stages of the method given by the current Butcher Tableau More...
 

Detailed Description

template<class ContainerType>
struct dg::DIRKStep< ContainerType >

Embedded diagonally implicit Runge Kutta time-step with error estimate \( \begin{align} k_i = f\left( t^n + c_i \Delta t, u^n + \Delta t \sum_{j=1}^{i} a_{ij} k_j\right) \\ u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s b_j k_j \\ \tilde u^{n+1} = u^{n} + \Delta t\sum_{j=1}^s \tilde b_j k_j \\ \delta^{n+1} = u^{n+1} - \tilde u^{n+1} = \Delta t\sum_{j=1}^s (b_j-\tilde b_j) k_j \end{align} \).

You can provide your own coefficients or use one of the methods in the following table:

We follow the naming convention of the ARKode library https://sundials.readthedocs.io/en/latest/arkode/Butcher_link.html (They also provide nice stability plots for their methods) as NAME-S-P-Q or NAME-S-Q, where

  • NAME is the author or name of the method
  • S is the number of stages in the method
  • P is the global order of the embedding
  • Q is the global order of the method
Name Identifier Description
Euler (implicit) dg::IMPLICIT_EULER_1_1 backward Euler \(a_{11} = 1\)
Midpoint (implicit) dg::IMPLICIT_MIDPOINT_1_2 implicit Midpoint \( a_{11} = 0.5 \)
Trapezoidal-2-2 dg::TRAPEZOIDAL_2_2 Crank-Nicolson method \( a_{11} = 0\ a_{22} = 0.5 \)
SDIRK-2-1-2 dg::SDIRK_2_1_2 generic 2nd order A and L-stable \( a_{ii} = 0.29 \)
Cavaglieri-3-1-2 (implicit) dg::CAVAGLIERI_IMPLICIT_3_1_2 Low-storage implicit/explicit Runge-Kutta schemes for the simulation of stiff high-dimensional ODE systems IMEXRKCB2 scheme \( a_{11} = 0 \)
Billington-3-3-2 dg::BILLINGTON_3_3_2 Billington-3-3-2 \( a_{ii} = 0.29\)
TRBDF2-3-3-2 dg::TRBDF2_3_3_2 TRBDF2-3-3-2 \( a_{11} = 0\)
Sanchez-3-3 dg::SANCHEZ_3_3 symplectic DIRK
Kvaerno-4-2-3 dg::KVAERNO_4_2_3 Kvaerno-4-2-3 \( a_{11} = 0\ a_{ii} = 0.44 \)
SDIRK-4-2-3 dg::SDIRK_4_2_3 Cameron2002 \( a_{ii}=0.25\)
Cavaglieri-4-2-3 (implicit) dg::CAVAGLIERI_IMPLICIT_4_2_3 Low-storage implicit/explicit Runge-Kutta schemes for the simulation of stiff high-dimensional ODE systems The SSP scheme IMEXRKCB3c (implicit part) \( a_{11} = 0 \)
ARK-4-2-3 (implicit) dg::ARK324L2SA_DIRK_4_2_3 ARK-4-2-3 (implicit) \( a_{11} = 0\ a_{ii} = 0.44 \)
Sanchez-3-4 dg::SANCHEZ_3_4 symplectic DIRK
Cash-5-2-4 dg::CASH_5_2_4 Cash-5-2-4 \( a_{ii} = 0.44 \)
Cash-5-3-4 dg::CASH_5_3_4 Cash-5-3-4 \( a_{ii} = 0.44 \)
SDIRK-5-3-4 dg::SDIRK_5_3_4 SDIRK-5-3-4 \( a_{ii} = 0.25\)
Kvaerno-5-3-4 dg::KVAERNO_5_3_4 Kvaerno-5-3-4 \( a_{11} = 0\ a_{ii} = 0.44\)
ARK-6-3-4 (implicit) dg::ARK436L2SA_DIRK_6_3_4 ARK-6-3-4 (implicit) \( a_{11} = 0\ a_{ii} = 0.25\)
Sanchez-6-5 dg::SANCHEZ_6_5 symplectic DIRK
Kvaerno-7-4-5 dg::KVAERNO_7_4_5 Kvaerno-7-4-5 \( a_{11} = 0 \ a_{ii} = 0.26\)
ARK-8-4-5 (implicit) dg::ARK548L2SA_DIRK_8_4_5 Kennedy and Carpenter (2019) Optimum ARK_2 method (implicit part) \( a_{11} = 0\ a_{ii} = 0.22\)
Sanchez-7-6 dg::SANCHEZ_7_6 symplectic DIRK
Note
Schemes marked with \( a_{11} = 0\) call the implicit_rhs instead of a solve once per step (the first stage is explicit, typically named EDIRK), all others never call implicit_rhs. Schemes marked with \(a_{ii} = ...\) (lower is better) have all diagonal elements equal (typically named SDIRK for "singly"). Schemes marked with both have equal elements only for i>1 (named ESDIRK).
Template Parameters
ContainerTypeAny class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag. Among others
  • dg::HVec (serial), dg::DVec (cuda / omp), dg::MHVec (mpi + serial) or dg::MDVec (mpi + cuda / omp)
  • std::vector<dg::DVec> (vector of shared device vectors), std::array<double, 4> (array of 4 doubles) or std::map < std::string, dg::DVec> ( a map of named vectors)
  • double (scalar) and other primitive types ...
If there are several ContainerTypes in the argument list, then TensorTraits must exist for all of them
See also
See The dg dispatch system for a detailed explanation of our type dispatch system

Member Typedef Documentation

◆ container_type

template<class ContainerType >
using dg::DIRKStep< ContainerType >::container_type = ContainerType

the type of the vector class in use

◆ value_type

template<class ContainerType >
using dg::DIRKStep< ContainerType >::value_type = get_value_type<ContainerType>

the value type of the time variable (float or double)

Constructor & Destructor Documentation

◆ DIRKStep() [1/2]

template<class ContainerType >
dg::DIRKStep< ContainerType >::DIRKStep ( )
inline

No memory allocation.

◆ DIRKStep() [2/2]

template<class ContainerType >
dg::DIRKStep< ContainerType >::DIRKStep ( ConvertsToButcherTableau< value_type im_tableau,
const ContainerType &  copyable 
)
inline

Construct with a diagonally implicit Butcher Tableau.

The tableau may be one of the implict methods listed in ConvertsToButcherTableau, or you provide your own tableau.

Parameters
im_tableaudiagonally implicit tableau, name or identifier that ConvertsToButcherTableau
copyablevector of the size that is later used in step ( it does not matter what values copyable contains, but its size is important; the step method can only be called with vectors of the same size)

Member Function Documentation

◆ construct()

template<class ContainerType >
template<class ... Params>
void dg::DIRKStep< ContainerType >::construct ( Params &&...  ps)
inline

Perfect forward parameters to one of the constructors.

Template Parameters
Paramsdeduced by the compiler
Parameters
psparameters forwarded to constructors

◆ copyable()

template<class ContainerType >
const ContainerType & dg::DIRKStep< ContainerType >::copyable ( ) const
inline

Return an object of same size as the object used for construction.

Returns
A copyable object; what it contains is undefined, its size is important

◆ embedded_order()

template<class ContainerType >
unsigned dg::DIRKStep< ContainerType >::embedded_order ( ) const
inline

global order of the embedding given by the current Butcher Tableau

◆ num_stages()

template<class ContainerType >
unsigned dg::DIRKStep< ContainerType >::num_stages ( ) const
inline

number of stages of the method given by the current Butcher Tableau

◆ order()

template<class ContainerType >
unsigned dg::DIRKStep< ContainerType >::order ( ) const
inline

global order of the method given by the current Butcher Tableau

◆ step() [1/2]

template<class ContainerType >
template<class ImplicitRHS , class Solver >
void dg::DIRKStep< ContainerType >::step ( const std::tuple< ImplicitRHS, Solver > &  ode,
value_type  t0,
const ContainerType &  u0,
value_type t1,
ContainerType &  u1,
value_type  dt 
)
inline

Advance one step ignoring error estimate and embedded method.

Template Parameters
ImplicitRHSThe implicit (part of the) right hand side is a functor type with no return value (subroutine) of signature void operator()(value_type, const ContainerType&, ContainerType&) The first argument is the time, the second is the input vector, which the functor may not override, and the third is the output, i.e. y' = I(t, y) translates to I(t, y, y'). The two ContainerType arguments never alias each other in calls to the functor. The functor can throw to indicate failure. Exceptions should derive from std::exception.
SolverA functor type for the implicit right hand side. Must solve the equation \( y - \alpha I(y,t) = y^*\) for y for given alpha, t and ys. Alpha is always positive and non-zero. Signature void operator()( value_type alpha, value_type t, ContainerType& y, const ContainerType& ys); The functor can throw. Any Exception should derive from std::exception.
Parameters
odethe <right hand side, solver for the rhs> functors. Typically std::tie(implicit_rhs, solver)
Attention
Contrary to EDIRK methods DIRK and SDIRK methods (all diagonal elements non-zero) never call implicit_rhs
Parameters
t0start time
u0value at t0
t1(write only) end time ( equals t0+dt on return may alias t0)
u1(write only) contains result on return (may alias u0)
dttimestep

◆ step() [2/2]

template<class ContainerType >
template<class ImplicitRHS , class Solver >
void dg::DIRKStep< ContainerType >::step ( const std::tuple< ImplicitRHS, Solver > &  ode,
value_type  t0,
const ContainerType &  u0,
value_type t1,
ContainerType &  u1,
value_type  dt,
ContainerType &  delta 
)
inline

Advance one step with error estimate.

Template Parameters
ImplicitRHSThe implicit (part of the) right hand side is a functor type with no return value (subroutine) of signature void operator()(value_type, const ContainerType&, ContainerType&) The first argument is the time, the second is the input vector, which the functor may not override, and the third is the output, i.e. y' = I(t, y) translates to I(t, y, y'). The two ContainerType arguments never alias each other in calls to the functor. The functor can throw to indicate failure. Exceptions should derive from std::exception.
SolverA functor type for the implicit right hand side. Must solve the equation \( y - \alpha I(y,t) = y^*\) for y for given alpha, t and ys. Alpha is always positive and non-zero. Signature void operator()( value_type alpha, value_type t, ContainerType& y, const ContainerType& ys); The functor can throw. Any Exception should derive from std::exception.
Parameters
odethe <right hand side, solver for the rhs> functors. Typically std::tie(implicit_rhs, solver)
Attention
Contrary to EDIRK methods DIRK and SDIRK methods (all diagonal elements non-zero) never call implicit_rhs
Parameters
t0start time
u0value at t0
t1(write only) end time ( equals t0+dt on return may alias t0)
u1(write only) contains result on return (may alias u0)
dttimestep
deltaContains error estimate (u1 - tilde u1) on return (must have equal size as u0)

The documentation for this struct was generated from the following file: