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

Convert identifiers to their corresponding dg::ButcherTableau. More...

Public Types

using value_type = real_type
 

Public Member Functions

 ConvertsToButcherTableau (ButcherTableau< real_type > tableau)
 
 ConvertsToButcherTableau (enum tableau_identifier id)
 Create ButcherTableau from dg::tableau_identifier. More...
 
 ConvertsToButcherTableau (std::string name)
 Create ButcherTableau from its name (very useful) More...
 
 ConvertsToButcherTableau (const char *name)
 Create ButcherTableau from its name (very useful) More...
 
 operator ButcherTableau< real_type > () const
 

Detailed Description

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

Convert identifiers to their corresponding dg::ButcherTableau.

This is a helper class to simplify the interfaces of our timestepper functions and classes. The sole purpose is to implicitly convert either a ButcherTableau or one of the following identifiers to an instance of a ButcherTableau.

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).
Parameters
real_typeThe type of the coefficients in the ButcherTableau

Member Typedef Documentation

◆ value_type

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

Constructor & Destructor Documentation

◆ ConvertsToButcherTableau() [1/4]

template<class real_type >
dg::ConvertsToButcherTableau< real_type >::ConvertsToButcherTableau ( ButcherTableau< real_type >  tableau)
inline

Of course a ButcherTableau converts to a ButcherTableau Useful if you constructed your very own coefficients

◆ ConvertsToButcherTableau() [2/4]

template<class real_type >
dg::ConvertsToButcherTableau< real_type >::ConvertsToButcherTableau ( enum tableau_identifier  id)
inline

Create ButcherTableau from dg::tableau_identifier.

The use of this constructor might be a bit awkward because you'll have to write all caps.

Parameters
idthe identifier, for example dg::DORMAND_PRINCE_7_4_5

◆ ConvertsToButcherTableau() [3/4]

template<class real_type >
dg::ConvertsToButcherTableau< real_type >::ConvertsToButcherTableau ( std::string  name)
inline

Create ButcherTableau from its name (very useful)

Note
In some of the links in the Description below you might want to use the search function of your browser to find the indicated method

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).
Parameters
nameThe name of the tableau as stated in the Name column above, as a string, for example "Dormand-Prince-7-4-5"

◆ ConvertsToButcherTableau() [4/4]

template<class real_type >
dg::ConvertsToButcherTableau< real_type >::ConvertsToButcherTableau ( const char *  name)
inline

Create ButcherTableau from its name (very useful)

Note
In some of the links in the Description below you might want to use the search function of your browser to find the indicated method

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).
Parameters
nameThe name of the tableau as stated in the Name column above, as a string, for example "Dormand-Prince-7-4-5"

Member Function Documentation

◆ operator ButcherTableau< real_type >()

template<class real_type >
dg::ConvertsToButcherTableau< real_type >::operator ButcherTableau< real_type > ( ) const
inline

Convert to ButcherTableau

which means an object can be directly assigned to a ButcherTableau


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