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

Additive Runge Kutta (semi-implicit) time-step with error estimate following The ARKode library 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

 ARKStep ()
 No memory allocation. More...
 
 ARKStep (std::string name, const ContainerType &copyable)
 Construct with given name. More...
 
 ARKStep (ConvertsToButcherTableau< value_type > ex_tableau, ConvertsToButcherTableau< value_type > im_tableau, const ContainerType &copyable)
 Construct with two Butcher Tableaus. 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 ExplicitRHS , class ImplicitRHS , class Solver >
void step (const std::tuple< ExplicitRHS, ImplicitRHS, Solver > &ode, value_type t0, const ContainerType &u0, value_type &t1, ContainerType &u1, value_type dt, ContainerType &delta)
 Advance one step. 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::ARKStep< ContainerType >

Additive Runge Kutta (semi-implicit) time-step with error estimate following The ARKode library

Currently, the possible Butcher Tableaus for a fully implicit-explicit scheme are the "Cavaglieri-3-1-2", "Cavaglieri-4-2-3", "ARK-4-2-3", "ARK-6-3-4" and "ARK-8-4-5" combinations.

Note
All currently possible schemes enjoy the FSAL qualitiy in the sense that only s-1 implicit solves and 1 evaluation of the implicit part are needed per step; the Cavaglieri methods do not require evaluations of the implicit part at all
Attention
When you use the ARKStep in combination with the Adaptive time step algorithm pay attention to solve the implicit part with sufficient accuracy. Else, the error propagates into the time controller, which will then choose the timestep as if the implicit part was explicit i.e. far too small. This might have to do with stiffness-leakage [Kennedy and Carpenter, Appl. num. Math., (2003)]: "An essential requirement for the viability of stiff/nonstiff IMEX schemes is that the stiffness remains truely separable. If this were not the case then stiffness would leak out of the stiff terms and stiffen the nonstiff terms. It would manifest itself as a loss in stability or a forced reduction in stepsize of the nonstiff terms. A more expensive fully implicit approach might then be required, and hence, methods that leak substantial stiffness might best be avoided".
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::ARKStep< ContainerType >::container_type = ContainerType

the type of the vector class in use

◆ value_type

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

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

Constructor & Destructor Documentation

◆ ARKStep() [1/3]

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

No memory allocation.

◆ ARKStep() [2/3]

template<class ContainerType >
dg::ARKStep< ContainerType >::ARKStep ( std::string  name,
const ContainerType &  copyable 
)
inline

Construct with given name.

Parameters
nameCurrently, one of "Cavaglieri-3-1-2", "Cavaglieri-4-2-3", "ARK-4-2-3", "ARK-6-3-4" or "ARK-8-4-5"
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)

◆ ARKStep() [3/3]

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

Construct with two Butcher Tableaus.

The two Butcher Tableaus represent the parameters for the explicit and implicit parts respectively. If both the explicit and implicit part of your equations are nontrivial, they must be one of the "ARK-X-X-X" methods listed in ConvertsToButcherTableau. Or you have your own tableaus of course but both tableaus must have the same number of steps.

Parameters
ex_tableauTableau for the explicit part
im_tableauTableau for the implicit part (must have the same number of stages as ex_tableau )
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::ARKStep< 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::ARKStep< 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::ARKStep< ContainerType >::embedded_order ( ) const
inline

global order of the embedding given by the current Butcher Tableau

◆ num_stages()

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

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

◆ order()

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

global order of the method given by the current Butcher Tableau

◆ step()

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

Advance one step.

Template Parameters
ExplicitRHSThe explicit (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' = E(t, y) translates to E(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.
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 <explicit rhs, implicit rhs, solver for the rhs> functor. Typically std::tie(explicit_rhs, implicit_rhs, solver)
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)
Note
the implementation is such that on return the last call is the explicit part ex at the new (t1,u1). This is useful if ex holds state, which is then updated to the new timestep and/or if im changes the state of ex
After a solve we immediately call ex on the solution

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