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

Manage coefficients of Multistep methods. More...

Public Types

using value_type = real_type
 

Public Member Functions

 MultistepTableau ()
 No memory allocation. More...
 
 MultistepTableau (unsigned steps, unsigned order, const std::vector< real_type > &a_v, const std::vector< real_type > &b_v, const std::vector< real_type > &c_v)
 Construct a tableau. More...
 
real_type a (unsigned i)
 Read the a_i coefficients. More...
 
real_type ex (unsigned i)
 Read the explicit (b_i) coefficients. More...
 
real_type im (unsigned i)
 Read the implicit (c_i) coefficients. More...
 
unsigned steps () const
 The number of stages s. More...
 
unsigned order () const
 global order of accuracy for the method More...
 
bool isExplicit () const
 True if any of the explicit coefficients b_i are non-zero. More...
 
bool isImplicit () const
 True if any of the implicit coefficients c_i are non-zero. More...
 

Detailed Description

template<class real_type>
struct dg::MultistepTableau< real_type >

Manage coefficients of Multistep methods.

The general s-step multistep method has the form

\[ y^{n+1} = \sum_{i=0}^{s-1} a_i y^{n-i} + h \sum_{i=0}^{s-1} b_i E( t_{n-i}, y_{n-i}) + h \sum_{i=0}^s c_i I( t_{n+1-i}, y^{n+1-i}) \]

where E is the explicit and I is implicit part. A purely implicit method is one where all \( b_i\) are zero, while an explicit one is one where all \( c_i\) are zero. A tableau thus consists of the three arrays a, b and c the number of steps and the order of the method. Currently available methods are:

ImEx methods

We follow the naming convention as NAME-S-Q

  • NAME is the author or name of the method
  • S is the number of steps in the method
  • Q is the global order of the method
Name Identifier Description
ImEx-Euler-1-1 dg::IMEX_EULER_1_1 Explicit Euler combined with Implicit Euler
Euler dg::IMEX_EULER_1_1 For convenience
ImEx-Koto-2-2 dg::IMEX_KOTO_2_2 Koto T. Front. Math. China 2009, 4(1): 113-129 A stabilized 2nd order scheme with a large region of stability
ImEx-Adams-X-X dg::IMEX_ADAMS_X_X Hundsdorfer and Ruuth, Journal of Computational Physics 225 (2007)
Note
Possible values for X: 2 (C=0.44), 3 (C=0.16)
ImEx-BDF-X-X dg::IMEX_BDF_X_X The family of schems described in Hundsdorfer and Ruuth, Journal of Computational Physics 225 (2007)
The implicit part is a normal BDF scheme https://en.wikipedia.org/wiki/Backward_differentiation_formula while the explicit part equals the Minimal Projecting method by Alfeld, P., Math. Comput. 33.148 1195-1212 (1979) or extrapolated BDF in Hundsdorfer, W., Ruuth, S. J., & Spiteri, R. J. (2003). Monotonicity-preserving linear multistep methods. SIAM Journal on Numerical Analysis, 41(2), 605-623
Note
Possible values for X: 1 (C=1.00), 2 (C=0.63), 3 (C=0.39), 4 (C=0.22), 5 (C=0.09), 6
Note that X=3 is identical to the "Karniadakis" scheme
Karniadakis dg::IMEX_BDF_3_3 The ImEx-BDF-3-3 scheme is identical to the widely used "Karniadakis" scheme Karniadakis, et al. J. Comput. Phys. 97 (1991)
ImEx-TVB-X-X dg::IMEX_TVB_X_X The family of schems described in < Hundsdorfer and Ruuth, Journal of Computational Physics 225 (2007)
The explicit part is a TVB scheme while the implicit part is optimized to maximize damping of high wavelength
Note
Possible values for X: 3 (C=0.54), 4 (C=0.46), 5 (C=0.38)
Note
the CFL coefficient C is given relative to the forward Euler method: \( \Delta t < C \Delta t_{FE}\).
Attention
The coefficient C is the one that ensures the TVD property of the scheme and is not directly related to the stability region of the scheme
Note
ImEx multistep tableaus can be used in ExplicitMultistep, ImplicitMultistep and ImExMultistep

Explicit methods

We follow the naming convention as NAME-S-Q

  • NAME is the author or name of the method
  • S is the number of steps in the method
  • Q is the global order of the method
Name Identifier Description
AB-X-X dg::AB_X_X The family of schemes described in Linear multistep methods as Adams-Bashforth

\[ u^{n+1} = u^n + \Delta t\sum_{j=0}^{s-1} b_j f\left(t^n - j \Delta t, u^{n-j}\right) \]

Note
Possible stages are X: 1, 2,..., 5, the order of the method is the same as its stages
The Adams-Bashforth schemes implemented here need less storage but may have a smaller region of absolute stability than for example an extrapolated BDF method of the same order.
eBDF-X-X dg::eBDF_X_X The family of schemes described in Hundsdorfer, W., Ruuth, S. J., & Spiteri, R. J. (2003). Monotonicity-preserving linear multistep methods. SIAM Journal on Numerical Analysis, 41(2), 605-623 as extrapolated BDF where it is found to be TVB (total variation bound). The schemes also appear as Minimal Projecting scheme described in Alfeld, P., Math. Comput. 33.148 1195-1212 (1979)
Note
Possible stages are X: 1 (C=1), 2 (C=0.63), 3 (C=0.39), 4 (C=0.22), 5 (C=0.09), 6 with the order the same as the number of stages
TVB-X-X dg::TVB_X_X The family of schemes described in S.J. Ruuth and W. Hundsdorfer, High-order linear multistep methods with general monotonicity and boundedness properties, Journal of Computational Physics, Volume 209, Issue 1, 2005 as Total variation Bound. These schemes have larger allowable step sizes than the eBDF family,
Note
Possible values for X are 1 (C=1), 2 (C=0.5), 3 (C=0.54), 4 (C=0.46), 5 (C=0.38) 6 (C=0.33). We highlight that TVB-3-3 has 38% larger allowable stepsize than eBDF-3-3 and TVB-4-4 has 109% larger stepsize than eBDF-4-4 (to ensure the TVB property, not stability).
SSP-X-Y dg::SSP_X_Y The family of schemes described in Gottlieb, S. On high order strong stability preserving runge-kutta and multi step time discretizations. J Sci Comput 25, 105–128 (2005) as Strong Stability preserving. We implement the lowest order schemes for each stage and disregard the remaining schemes in the paper since their CFL conditions are worse than the TVB scheme of the same order.
Note
Possible values for X-Y : 1-1 (C=1), 2-2 (C=0.5), 3-2 (C=0.5), 4-2 (C=0.66), 5-3 (C=0.5), 6-3 (C=0.567).
These schemes are noteworthy because the coefficients b_i are all positive except for the 2-2 method and the "4-2" and "6-3" methods allow slightly larger allowable stepsize but increased storage requirements than TVB of same order (2 and 3).
Note
Total variation bound (TVB) means \( || v^n|| \leq M ||v^0||\) where the norm signifies the total variation semi-norm. Total variation diminishing (TVD) means M=1, and strong stability preserving (SSP) is the same as TVD, TVB schemes converge to the correct entropy solutions of hyperbolic conservation laws
the CFL coefficient C is given relative to the forward Euler method: \( \Delta t < C \Delta t_{FE}\).
Attention
The coefficient C is the one that ensures the TVD property of the scheme and is not directly related to the stability region of the scheme
Template Parameters
real_typetype of the coefficients
See also
ExplicitMultistep, ImplicitMultistep, ImExMultistep

Member Typedef Documentation

◆ value_type

template<class real_type >
using dg::MultistepTableau< real_type >::value_type = real_type

Constructor & Destructor Documentation

◆ MultistepTableau() [1/2]

template<class real_type >
dg::MultistepTableau< real_type >::MultistepTableau ( )
inline

No memory allocation.

◆ MultistepTableau() [2/2]

template<class real_type >
dg::MultistepTableau< real_type >::MultistepTableau ( unsigned  steps,
unsigned  order,
const std::vector< real_type > &  a_v,
const std::vector< real_type > &  b_v,
const std::vector< real_type > &  c_v 
)
inline

Construct a tableau.

Parameters
stepsnumber of stages
order(global) order of the resulting method
a_vs real numbers
b_vs real numbers (can be empty, which then sets them all to 0 constructing an implicit method)
c_vs+1 real numbers (can be empty, which constructs an explicit method by assigning all c_i 0)

Member Function Documentation

◆ a()

template<class real_type >
real_type dg::MultistepTableau< real_type >::a ( unsigned  i)
inline

Read the a_i coefficients.

Parameters
iidx number 0<=i<s, i>=s results in undefined behaviour
Returns
a_i

◆ ex()

template<class real_type >
real_type dg::MultistepTableau< real_type >::ex ( unsigned  i)
inline

Read the explicit (b_i) coefficients.

Parameters
iidx number 0<=i<s, i>=s results in undefined behaviour
Returns
b_i

◆ im()

template<class real_type >
real_type dg::MultistepTableau< real_type >::im ( unsigned  i)
inline

Read the implicit (c_i) coefficients.

Parameters
iidx number 0<=i<s+1, i>=s+1 results in undefined behaviour
Returns
c_i

◆ isExplicit()

template<class real_type >
bool dg::MultistepTableau< real_type >::isExplicit ( ) const
inline

True if any of the explicit coefficients b_i are non-zero.

◆ isImplicit()

template<class real_type >
bool dg::MultistepTableau< real_type >::isImplicit ( ) const
inline

True if any of the implicit coefficients c_i are non-zero.

◆ order()

template<class real_type >
unsigned dg::MultistepTableau< real_type >::order ( ) const
inline

global order of accuracy for the method

◆ steps()

template<class real_type >
unsigned dg::MultistepTableau< real_type >::steps ( ) const
inline

The number of stages s.


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