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

Manage coefficients of a (extended) Butcher tableau. More...

Public Types

using value_type = real_type
 

Public Member Functions

 ButcherTableau ()=default
 No memory allocation. More...
 
 ButcherTableau (unsigned s, unsigned order, const real_type *a, const real_type *b, const real_type *c)
 Construct a classic non-embedded tableau. More...
 
 ButcherTableau (unsigned s, unsigned embedded_order, unsigned order, const real_type *a, const real_type *b, const real_type *bt, const real_type *c)
 Construct an embedded tableau. More...
 
 ButcherTableau (unsigned s, real_type *data)
 Construct from ARKode standard format. More...
 
real_type a (unsigned i, unsigned j) const
 Read the a_ij coefficients. More...
 
real_type c (unsigned i) const
 Read the c_i coefficients. More...
 
real_type b (unsigned j) const
 Read the b_j coefficients. More...
 
real_type bt (unsigned j) const
 Read the embedded bt_j coefficients. More...
 
real_type d (unsigned j) const
 Return the coefficients for the error estimate Equivalent to b(j)-bt(j) More...
 
unsigned num_stages () const
 The number of stages s. More...
 
unsigned order () const
 global order of accuracy for the method represented by b More...
 
unsigned embedded_order () const
 global order of accuracy for the embedded method represented by bt More...
 
bool isEmbedded () const
 True if the method has an embedding. More...
 
bool isImplicit () const
 True if an element on or above the diagonal in a is non-zero. More...
 
bool isFsal () const
 

Detailed Description

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

Manage coefficients of a (extended) Butcher tableau.

The goal of this class is to represent a Butcher tableau for the use in Runge Kutta type ODE integrators. See https://en.wikipedia.org/wiki/Runge%E2%80%93Kutta_methods for an introduction. The coefficients of the tableau are easily constructible and accessible. Furthermore, we provide utilities like the number of stages, whether the tableau is embedded or not and the order of the method.

Currently available are

Explicit methods

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)

Implicit methods

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
real_typetype of the coefficients
See also
RungeKutta, ERKStep, ARKStep

Member Typedef Documentation

◆ value_type

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

Constructor & Destructor Documentation

◆ ButcherTableau() [1/4]

template<class real_type >
dg::ButcherTableau< real_type >::ButcherTableau ( )
default

No memory allocation.

◆ ButcherTableau() [2/4]

template<class real_type >
dg::ButcherTableau< real_type >::ButcherTableau ( unsigned  s,
unsigned  order,
const real_type *  a,
const real_type *  b,
const real_type *  c 
)
inline

Construct a classic non-embedded tableau.

Parameters
snumber of stages
order(global) order of the resulting method
apointer to s*s real numbers interpreted as a_{ij}=a[i*s+j]
bpointer to s real numbers interpreted as b_j=b[j]
cpointer to s real numbers interpreted as c_i=c[i]
Note
This constructor initializes the embedded coefficients bt=b which makes the embedded method equal to the actual method

◆ ButcherTableau() [3/4]

template<class real_type >
dg::ButcherTableau< real_type >::ButcherTableau ( unsigned  s,
unsigned  embedded_order,
unsigned  order,
const real_type *  a,
const real_type *  b,
const real_type *  bt,
const real_type *  c 
)
inline

Construct an embedded tableau.

Parameters
snumber of stages
embedded_order(global) order of the embedded method (corresponding to bt)
order(global) order of the method (corresponding to b)
apointer to s*s real numbers interpreted as a_{ij}=a[i*s+j]
bpointer to s real numbers interpreted as b_j=b[j]
btpointer to s real numbers interpreted as bt_j=bt[j]
cpointer to s real numbers interpreted as c_i=c[i]

◆ ButcherTableau() [4/4]

template<class real_type >
dg::ButcherTableau< real_type >::ButcherTableau ( unsigned  s,
real_type *  data 
)
inline

Construct from ARKode standard format.

Construct embedded method from single array

Parameters
snumber of stages
dataArray of (s+1)*(s+2) real numbers, interpreted as: c[i]=data[i*(s+1)], a_ij=data[i*(s+1)+j+1], order=data[s*(s+1)], b_j=data[s*(s+1)+j+1], embedded_order = data[(s+1)*(s+1)], bt_j=data[(s+1)*(s+1)+j+1]

Member Function Documentation

◆ a()

template<class real_type >
real_type dg::ButcherTableau< real_type >::a ( unsigned  i,
unsigned  j 
) const
inline

Read the a_ij coefficients.

Parameters
irow number 0<=i<s, i>=s results in undefined behaviour
jcol number 0<=j<s, j>=s results in undefined behaviour
Returns
a_ij

◆ b()

template<class real_type >
real_type dg::ButcherTableau< real_type >::b ( unsigned  j) const
inline

Read the b_j coefficients.

Parameters
jcol number 0<=j<s, j>=s results in undefined behaviour
Returns
b_j

◆ bt()

template<class real_type >
real_type dg::ButcherTableau< real_type >::bt ( unsigned  j) const
inline

Read the embedded bt_j coefficients.

Parameters
jcol number 0<=j<s, j>=s results in undefined behaviour
Returns
bt_j

◆ c()

template<class real_type >
real_type dg::ButcherTableau< real_type >::c ( unsigned  i) const
inline

Read the c_i coefficients.

Parameters
irow number 0<=i<s, i>=s results in undefined behaviour
Returns
c_i

◆ d()

template<class real_type >
real_type dg::ButcherTableau< real_type >::d ( unsigned  j) const
inline

Return the coefficients for the error estimate Equivalent to b(j)-bt(j)

Parameters
jcol number 0<=j<s, j>=s results in undefined behaviour
Returns
b(j)-bt(j)

◆ embedded_order()

template<class real_type >
unsigned dg::ButcherTableau< real_type >::embedded_order ( ) const
inline

global order of accuracy for the embedded method represented by bt

◆ isEmbedded()

template<class real_type >
bool dg::ButcherTableau< real_type >::isEmbedded ( ) const
inline

True if the method has an embedding.

◆ isFsal()

template<class real_type >
bool dg::ButcherTableau< real_type >::isFsal ( ) const
inline

True if the method has the "First Same As Last" property: the last stage is evaluated at the same point as the first stage of the next step.

◆ isImplicit()

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

True if an element on or above the diagonal in a is non-zero.

◆ num_stages()

template<class real_type >
unsigned dg::ButcherTableau< real_type >::num_stages ( ) const
inline

The number of stages s.

◆ order()

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

global order of accuracy for the method represented by b


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