Shu-Osher fixed-step explicit ODE integrator with Slope Limiter / Filter \( \begin{align} u_0 &= u_n \\ u_i &= \Lambda\Pi \left(\sum_{j=0}^{i-1}\left[ \alpha_{ij} u_j + \Delta t \beta_{ij} f( t_j, u_j)\right]\right)\\ u^{n+1} &= u_s \end{align} \).
More...
template<class ContainerType>
struct dg::ShuOsher< ContainerType >
Shu-Osher fixed-step explicit ODE integrator with Slope Limiter / Filter \( \begin{align} u_0 &= u_n \\ u_i &= \Lambda\Pi \left(\sum_{j=0}^{i-1}\left[ \alpha_{ij} u_j + \Delta t \beta_{ij} f( t_j, u_j)\right]\right)\\ u^{n+1} &= u_s \end{align} \).
where \( \Lambda\Pi\) is the limiter, \( i\in [1,s]\) and s
is the number of stages (i.e. the number of times the right hand side is evaluated.
The method is defined by its (explicit) ShuOsherTableau, given by the coefficients alpha
and beta
, and s
is the number of stages.
- Note
- the original reference for the scheme is Chi-Wang Shu, Stanley Osher, Efficient implementation of essentially non-oscillatory shock-capturing schemes, Journal of Computational Physics, Volume 77, Issue 2, 1988, Pages 439-471
-
This struct can be used to implement the RKDG methods with slope-limiter described in Cockburn, B., Shu, CW. Runge–Kutta Discontinuous Galerkin Methods for Convection-Dominated Problems. Journal of Scientific Computing 16, 173–261 (2001)
You can use one of our predefined methods (only the ones that are marked with "Shu-Osher-Form"):
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 | dg::EXPLICIT_EULER_1_1 | Explicit Euler |
Midpoint-2-2 | dg::MIDPOINT_2_2 | Midpoint method |
Kutta-3-3 | dg::KUTTA_3_3 | Kutta-3-3 |
Runge-Kutta-4-4 | dg::CLASSIC_4_4 | "The" Runge-Kutta method |
SSPRK-2-2 | dg::SSPRK_2_2 | SSPRK (2,2) CFL_eff = 0.5 "Shu-Osher-Form" |
SSPRK-3-2 | dg::SSPRK_3_2 | SSPRK (3,2) CFL_eff = 0.66 "Shu-Osher-Form" |
SSPRK-3-3 | dg::SSPRK_3_3 | SSPRK (3,3) CFL_eff = 0.33 "Shu-Osher-Form" |
SSPRK-5-3 | dg::SSPRK_5_3 | SSPRK (5,3) CFL_eff = 0.5 "Shu-Osher-Form" |
SSPRK-5-4 | dg::SSPRK_5_4 | SSPRK (5,4) CFL_eff = 0.37 "Shu-Osher-Form" |
Heun-Euler-2-1-2 | dg::HEUN_EULER_2_1_2 | Heun-Euler-2-1-2 |
Cavaglieri-3-1-2 (explicit) | dg::CAVAGLIERI_3_1_2 | Low-storage implicit/explicit Runge-Kutta schemes for the simulation of stiff high-dimensional ODE systems IMEXRKCB2 scheme |
Fehlberg-3-2-3 | dg::FEHLBERG_3_2_3 | The original uses the embedding as the solution [Hairer, Noersett, Wanner, Solving ordinary differential Equations I, 1987] |
Fehlberg-4-2-3 | dg::FEHLBERG_4_2_3 | The original uses the embedding as the solution [Hairer, Noersett, Wanner, Solving ordinary differential Equations I, 1987] (fsal) |
Bogacki-Shampine-4-2-3 | dg::BOGACKI_SHAMPINE_4_2_3 | Bogacki-Shampine (fsal) |
Cavaglieri-4-2-3 (explicit) | dg::CAVAGLIERI_4_2_3 | Low-storage implicit/explicit Runge-Kutta schemes for the simulation of stiff high-dimensional ODE systems The SSP scheme IMEXRKCB3c (explicit part) |
ARK-4-2-3 (explicit) | dg::ARK324L2SA_ERK_4_2_3 | ARK-4-2-3 (explicit) |
Zonneveld-5-3-4 | dg::ZONNEVELD_5_3_4 | Zonneveld-5-3-4 |
ARK-6-3-4 (explicit) | dg::ARK436L2SA_ERK_6_3_4 | ARK-6-3-4 (explicit) |
Sayfy_Aburub-6-3-4 | dg::SAYFY_ABURUB_6_3_4 | Sayfy_Aburub_6_3_4 |
Cash_Karp-6-4-5 | dg::CASH_KARP_6_4_5 | Cash-Karp |
Fehlberg-6-4-5 | dg::FEHLBERG_6_4_5 | Runge-Kutta-Fehlberg |
Dormand-Prince-7-4-5 | dg::DORMAND_PRINCE_7_4_5 | Dormand-Prince method (fsal) |
Tsitouras09-7-4-5 | dg::TSITOURAS09_7_4_5 | Tsitouras 5(4) method from 2009 (fsal), The default method in Julia |
Tsitouras11-7-4-5 | dg::TSITOURAS11_7_4_5 | Tsitouras 5(4) method from 2011 (fsal) Further improves Tsitouras09 (Note that in the paper b-bt is printed instead of bt) |
ARK-8-4-5 (explicit) | dg::ARK548L2SA_ERK_8_4_5 | Kennedy and Carpenter (2019) Optimum ARK_2 method (explicit part) |
Verner-9-5-6 | dg::VERNER_9_5_6 | Verner-9-5-6 (fsal) |
Verner-10-6-7 | dg::VERNER_10_6_7 | Verner-10-6-7 |
Fehlberg-13-7-8 | dg::FEHLBERG_13_7_8 | Fehlberg-13-7-8 |
Dormand-Prince-13-7-8 | dg::DORMAND_PRINCE_13_7_8 | [Hairer, Noersett, Wanner, Solving ordinary differential Equations I, 1987] |
Feagin-17-8-10 | dg::FEAGIN_17_8_10 | Feagin (The RK10(8) method) |
- Note
- Uses only
dg::blas1
routines to integrate one step.
- Template Parameters
-
ContainerType | Any 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
◆ container_type
template<class ContainerType >
using dg::ShuOsher< ContainerType >::container_type = ContainerType |
the type of the vector class in use
◆ value_type
template<class ContainerType >
the value type of the time variable (float or double)
◆ ShuOsher() [1/2]
template<class ContainerType >
◆ ShuOsher() [2/2]
template<class ContainerType >
Reserve internal workspace for the integration.
- Parameters
-
tableau | Tableau, name or identifier that ConvertsToShuOsherTableau |
copyable | vector 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) |
◆ construct()
template<class ContainerType >
template<class ... Params>
void dg::ShuOsher< ContainerType >::construct |
( |
Params &&... |
ps | ) |
|
|
inline |
Perfect forward parameters to one of the constructors.
- Template Parameters
-
Params | deduced by the compiler |
- Parameters
-
ps | parameters forwarded to constructors |
◆ copyable()
template<class ContainerType >
const ContainerType & dg::ShuOsher< 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
◆ num_stages()
template<class ContainerType >
unsigned dg::ShuOsher< ContainerType >::num_stages |
( |
| ) |
const |
|
inline |
number of stages of the method given by the current Butcher Tableau
◆ order()
template<class ContainerType >
global order of the method given by the current Butcher Tableau
◆ step()
template<class ContainerType >
template<class ExplicitRHS , class Limiter >
Advance one step.
- Template Parameters
-
ExplicitRHS | The 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 . |
Limiter | The filter or limiter class to use in connection with the time-stepper has a member function operator() of signature void operator()( ContainerType&) The argument is the input vector that the function has to override i.e. y = L( y) translates to L( y). |
- Parameters
-
ode | right hand side subroutine and limiter to use. Typically std::tie( rhs,limiter) |
t0 | start time |
u0 | value 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) |
dt | timestep |
- Note
- on return
rhs(t1, u1)
will be the last call to rhs
(this is useful if ExplicitRHS
holds state, which is then updated to the current timestep)
The documentation for this struct was generated from the following file: