Struct ButcherTableau

Source
pub struct ButcherTableau<T: Real, const S: usize, const I: usize = S> {
    pub c: [T; I],
    pub a: [[T; I]; I],
    pub b: [T; S],
    pub bh: Option<[T; S]>,
    pub bi: Option<[[T; I]; I]>,
    pub er: Option<[T; S]>,
}
Expand description

Butcher Tableau structure for Runge-Kutta methods.

A Butcher tableau encodes the coefficients of a Runge-Kutta method and provides the necessary components for solving ordinary differential equations. This implementation includes support for embedded methods (for error estimation) and dense output through interpolation.

§Generic Parameters

  • T: The type of the coefficients, typically a floating-point type (e.g., f32, f64).
  • S: Number of stages in the method.
  • I: Primary Stages plus extra stages for interpolation (default is equal to S).

§Fields

  • c: Node coefficients (time steps within the interval).

  • a: Runge-Kutta matrix coefficients (coupling between stages).

  • b: Weight coefficients for the primary method’s final stage.

  • bh: Weight coefficients for the embedded method (used for error estimation).

  • bi: Weight coefficients for the interpolation method (used for dense output).

  • er: Error estimation coefficients (optional, not all adaptive methods have these).

    These allow approximation at any point within the integration step.

Fields§

§c: [T; I]§a: [[T; I]; I]§b: [T; S]§bh: Option<[T; S]>§bi: Option<[[T; I]; I]>§er: Option<[T; S]>

Implementations§

Source§

impl<T: Real> ButcherTableau<T, 7>

Source

pub fn dopri5() -> Self

Dormand-Prince 5(4) Tableau with dense output interpolation.

§Overview

This provides a 7-stage Runge-Kutta method with:

  • Primary order: 5
  • Embedded order: 4 (for error estimation)
  • Number of stages: 7 primary + 0 additional stages for interpolation
  • Built-in dense output of order 4
§Efficiency

The DOPRI5 method is popular due to its efficient balance between accuracy and computational cost. It is particularly good for non-stiff problems.

§Interpolation
  • The method provides a 4th-order interpolant using the existing 7 stages
  • The interpolant has continuous first derivatives
  • The interpolant uses the values of the stages to construct a polynomial that allows evaluation at any point within the integration step
§Notes
  • The method was developed by Dormand and Prince in 1980
  • It is one of the most widely used Runge-Kutta methods and is implemented in many software packages for solving ODEs
  • The DOPRI5 method is a member of the Dormand-Prince family of embedded Runge-Kutta methods
§References
  • Dormand, J. R. & Prince, P. J. (1980), “A family of embedded Runge-Kutta formulae”, Journal of Computational and Applied Mathematics, 6(1), pp. 19-26
  • Hairer, E., Nørsett, S. P. & Wanner, G. (1993), “Solving Ordinary Differential Equations I: Nonstiff Problems”, Springer Series in Computational Mathematics, Vol. 8, Springer-Verlag
Source§

impl<T: Real> ButcherTableau<T, 12, 16>

Source

pub fn dop853() -> Self

Dormand-Prince 8(5,3) Tableau with dense output interpolation.

§Overview

This provides a 12-stage Runge-Kutta method with:

  • Primary order: 8
  • Embedded order: 5 (for error estimation)
  • Number of stages: 12 primary + 4 additional stages for interpolation
  • Built-in dense output of order 7
§Efficiency

The DOP853 method is a high-order Runge-Kutta method that provides excellent accuracy and efficiency for non-stiff problems. It’s particularly suitable when high precision is required.

§Interpolation
  • The method provides a 7th-order interpolant using the existing stages plus 3 additional stages for dense output
  • The interpolant has continuous derivatives
  • The interpolation is performed through a sophisticated continuous extension that maintains high accuracy throughout the integration step
§Notes
  • The method was developed by Dormand, Prince and others as an extension of their earlier work
  • It is one of the most accurate explicit Runge-Kutta implementations available for solving ODEs
  • The DOP853 method is widely used in scientific computing where high precision is required
§References
  • Hairer, E., Nørsett, S. P. & Wanner, G. (1993), “Solving Ordinary Differential Equations I: Nonstiff Problems”, Springer Series in Computational Mathematics, Vol. 8, Springer-Verlag
  • Dormand, J. R. & Prince, P. J. (1980), “A family of embedded Runge-Kutta formulae”, Journal of Computational and Applied Mathematics, 6(1), pp. 19-26
Source§

impl<T: Real> ButcherTableau<T, 4>

Source

pub fn rk4() -> Self

Classic Runge-Kutta 4th order method (RK4).

§Overview

This provides a 4-stage, explicit Runge-Kutta method with:

  • Primary order: 4
  • Embedded order: None (this implementation does not include an embedded error estimate)
  • Number of stages: 4
§Interpolation
  • This standard RK4 implementation does not provide coefficients for dense output (interpolation).
§Notes
  • RK4 is a widely used, general-purpose method known for its balance of accuracy and computational efficiency for many problems.
  • It is a fixed-step method as presented here, lacking an embedded formula for adaptive step-size control. For adaptive step sizes, a method with an error estimate (e.g., RKF45) would be required.
§Butcher Tableau
0   |
1/2 | 1/2
1/2 | 0   1/2
1   | 0   0   1
----|--------------------
    | 1/6 1/3 1/3 1/6
§References
  • Kutta, W. (1901). “Beitrag zur näherungsweisen Integration totaler Differentialgleichungen”. Zeitschrift für Mathematik und Physik, 46, 435-453.
  • Runge, C. (1895). “Über die numerische Auflösung von Differentialgleichungen”. Mathematische Annalen, 46, 167-178.
Source

pub fn three_eighths() -> Self

Three-Eighths Rule 4th order method.

§Overview

This provides a 4-stage, explicit Runge-Kutta method with:

  • Primary order: 4
  • Embedded order: None (no embedded error estimate)
  • Number of stages: 4
§Notes
  • The primary advantage of this method is that almost all of the error coefficients are smaller than in the classic RK4 method, but it requires slightly more FLOPs (floating-point operations) per time step.
§Butcher Tableau
0   |
1/3 | 1/3
2/3 | -1/3 1
1   | 1   -1   1
----|--------------------
    | 1/8 3/8 3/8 1/8
§References
  • Butcher, J.C. (2008). “Numerical Methods for Ordinary Differential Equations”.
Source§

impl<T: Real> ButcherTableau<T, 2>

Source

pub fn midpoint() -> Self

Midpoint method (2nd order Runge-Kutta).

§Overview

This provides a 2-stage, explicit Runge-Kutta method with:

  • Primary order: 2
  • Embedded order: None
  • Number of stages: 2
§Butcher Tableau
0   |
1/2 | 1/2
----|--------
    | 0   1
§References
  • Butcher, J.C. (2008). “Numerical Methods for Ordinary Differential Equations”.
Source

pub fn heun() -> Self

Heun’s method (2nd order Runge-Kutta).

§Overview

This provides a 2-stage, explicit Runge-Kutta method with:

  • Primary order: 2
  • Embedded order: None
  • Number of stages: 2
§Butcher Tableau
0   |
1   | 1
----|--------
    | 1/2 1/2
§References
  • Heun, K. (1900). “Neue Methoden zur approximativen Integration der Differentialgleichungen einer unabhängigen Veränderlichen”.
Source

pub fn ralston() -> Self

Ralston’s method (2nd order Runge-Kutta).

§Overview

This provides a 2-stage, explicit Runge-Kutta method with:

  • Primary order: 2
  • Embedded order: None
  • Number of stages: 2
§Butcher Tableau
0   |
2/3 | 2/3
----|--------
    | 1/4 3/4
§References
  • Ralston, A. (1962). “Runge-Kutta Methods with Minimum Error Bounds”.
Source§

impl<T: Real> ButcherTableau<T, 1>

Source

pub fn euler() -> Self

Euler’s method (1st order Runge-Kutta).

§Overview

This provides a 1-stage, explicit Runge-Kutta method with:

  • Primary order: 1
  • Embedded order: None
  • Number of stages: 1
§Butcher Tableau
0 |
--|--
  | 1
§References
  • Euler, L. (1768). “Institutionum calculi integralis”.
Source§

impl<T: Real> ButcherTableau<T, 6>

Source

pub fn rkf45() -> Self

Runge-Kutta-Fehlberg 4(5) method (RKF45).

§Overview

This provides a 6-stage, explicit Runge-Kutta method with:

  • Primary order: 5
  • Embedded order: 4 (for error estimation)
  • Number of stages: 6
§Notes
  • RKF45 is a widely used adaptive step size method that provides a good balance between accuracy and computational efficiency.
  • It uses the difference between 4th and 5th order approximations to estimate error.
§Butcher Tableau
0      |
1/4    | 1/4
3/8    | 3/32         9/32
12/13  | 1932/2197    -7200/2197  7296/2197
1      | 439/216      -8          3680/513    -845/4104
1/2    | -8/27        2           -3544/2565  1859/4104   -11/40
-------|---------------------------------------------------------------
       | 16/135       0           6656/12825  28561/56430 -9/50  2/55  (5th)
       | 25/216       0           1408/2565   2197/4104   -1/5   0     (4th)
§References
  • Fehlberg, E. (1969). “Low-order classical Runge-Kutta formulas with step size control and their application to some heat transfer problems”.
Source

pub fn cash_karp() -> Self

Cash-Karp 4(5) method.

§Overview

This provides a 6-stage, explicit Runge-Kutta method with:

  • Primary order: 5
  • Embedded order: 4 (for error estimation)
  • Number of stages: 6
§Notes
  • The Cash-Karp method is a variant of Runge-Kutta methods with embedded error estimation that often provides better accuracy than RKF45 for some problem types.
§Butcher Tableau
0      |
1/5    | 1/5
3/10   | 3/40         9/40
3/5    | 3/10         -9/10       6/5
1      | -11/54       5/2         -70/27      35/27
7/8    | 1631/55296   175/512     575/13824   44275/110592 253/4096
-------|---------------------------------------------------------------
       | 37/378       0           250/621     125/594     0      512/1771  (5th)
       | 2825/27648   0           18575/48384 13525/55296 277/14336 1/4    (4th)
§References
  • Cash, J.R., Karp, A.H. (1990). “A Variable Order Runge-Kutta Method for Initial Value Problems with Rapidly Varying Right-Hand Sides”.
Source§

impl<T: Real> ButcherTableau<T, 1>

Source

pub fn backward_euler() -> Self

Source§

impl<T: Real> ButcherTableau<T, 2>

Source

pub fn trapezoidal() -> Self

Source

pub fn crank_nicolson() -> Self

Source§

impl<T: Real> ButcherTableau<T, 9, 10>

Source

pub fn rkv655e() -> Self

Verner’s RKV(6,5,5)e method with 5th-order interpolation.

A ‘most efficient’ Runge-Kutta (9,6(5)) FSAL pair with rational coefficients.

§Overview

This is a nine-stage FSAL pair of methods of orders p=6 and p=5, with dominant stage order 3. The coefficients are exact rationals. The method is “most efficient” in the sense that the maximum coefficient in b and A is not large, and the propagating formula nearly minimizes the 2-norm of the local truncation error (T_72 ≈ 1.44e-6).

§Interpolation

Additional stages and interpolating weights allow computation of an approximation at any point in the domain of solution of order up to p. These interpolants have continuous derivatives.

The additional node c[9]=1/2 is chosen to achieve order 5 interpolation in ten stages, with maximum 2-norm of the local truncation error Ti_62 ≈ 3.18e-4 (three local maxima on [0,1]).

The remaining two nodes are chosen to minimize the maximum 2-norm of the local truncation error for the order 6 interpolant, Ti_72 ≈ 4.44e-5 (two local maxima on [0,1]).

§Source

Verner’s RKV655e

Source§

impl<T: Real> ButcherTableau<T, 9, 12>

Source

pub fn rkv656e() -> Self

Verner’s RKV(6,5,5)e method with 5th-order interpolation.

A ‘most efficient’ Runge-Kutta (9,6(5)) FSAL pair with rational coefficients.

§Overview

This is a nine-stage FSAL pair of methods of orders p=6 and p=5, with dominant stage order 3. The coefficients are exact rationals. The method is “most efficient” in the sense that the maximum coefficient in b and A is not large, and the propagating formula nearly minimizes the 2-norm of the local truncation error (T_72 ≈ 1.44e-6).

§Interpolation

Additional stages and interpolating weights allow computation of an approximation at any point in the domain of solution of order up to p. These interpolants have continuous derivatives.

The additional node c[9]=1/2 is chosen to achieve order 5 interpolation in ten stages, with maximum 2-norm of the local truncation error Ti_62 ≈ 3.18e-4 (three local maxima on [0,1]).

The remaining two nodes are chosen to minimize the maximum 2-norm of the local truncation error for the order 6 interpolant, Ti_72 ≈ 4.44e-5 (two local maxima on [0,1]).

§Source

Verner’s RKV655e

Source§

impl<T: Real> ButcherTableau<T, 10, 13>

Source

pub fn rkv766e() -> Self

A ‘most efficient’ Runge-Kutta (10:6(7)) pair with 6th-order interpolation.

§Overview

This provides a 10-stage Runge-Kutta method with:

  • Primary order: 6
  • Embedded order: 7 (for error estimation)
  • Number of stages: 10 primary + 3 additional stages for interpolation
  • Dominant stage order: 3
§Efficiency

This method is considered “most efficient” in the sense that for a specified maximum coefficient from b and A, it almost minimizes the local truncation error 2-norm:

T_82 ~ 0.000003389

(Formulas with slightly different nodes may have a slightly smaller error norm, perhaps achieved by having a larger maximum coefficient.)

§Interpolation
  • Node c[11]=1 was chosen to make the interpolants differentiable

  • Node c[12] was chosen as a zero of a quadratic polynomial to enable efficient interpolation

  • The interpolant has continuous derivatives

  • Maximum 2-norm of the local truncation error over [0,1] for the interpolant:

    Ti_72 ~ 0.0000165

  • The interpolation error has four local maximum values on [0,1]

§Source

Verner’s RKV766e

Source§

impl<T: Real> ButcherTableau<T, 10, 16>

Source

pub fn rkv767e() -> Self

A ‘most efficient’ Runge-Kutta (10:7(6)) pair with 7th-order interpolation.

§Overview

This provides a 10-stage Runge-Kutta method with:

  • Primary order: 7
  • Embedded order: 6 (for error estimation)
  • Number of stages: 10 primary + 6 additional stages for interpolation
  • Dominant stage order: 3
§Efficiency

This method is considered “most efficient” in the sense that for a specified maximum coefficient from b and A, it almost minimizes the local truncation error 2-norm:

T_82 ~ 0.000003389

§Interpolation
  • Additional stages allow computation of approximations at any point with order up to 7

  • The interpolant has continuous derivatives

  • Maximum 2-norm of the local truncation error for the 7th-order interpolant:

    Ti_82 ~ 0.000003389

  • The 2-norm is bounded by this value on [0,1] and is nearly monotone increasing

  • The 7th-order interpolant requires 6 additional stages (16 total)

§Notes
  • This method uses the same primary stages as the 6th-order method (rkv766e)
  • The higher-order interpolant requires more additional stages
§Source

Verner’s RKV766e

Source§

impl<T: Real> ButcherTableau<T, 13, 17>

Source

pub fn rkv877e() -> Self

An efficient Runge-Kutta (13:7) pair method with 7th-order interpolation.

§Overview

This provides a 13-stage Runge-Kutta method with:

  • Primary order: 7
  • Embedded order: 7 (used for stepsize control)
  • Number of stages: 13 primary + 4 additional stages for interpolation
  • Dominant stage order: 4
§Interpolation
  • The method provides a 7th-order interpolant using the nodes c[14]=1, c[15]-c[17].

  • The interpolant has continuous derivatives.

  • Maximum 2-norm of the local truncation error over [0,1] for the interpolant:

    Ti_82 ~ 0.0000063833

  • The interpolation error has three local maximum values on [0,1].

§Notes
  • This method uses the same initial 13 stages as the 8th-order method (rkv878e).
  • The interpolant requires fewer additional stages than the 8th-order version.
§Source

Verner’s RKV878e

Source§

impl<T: Real> ButcherTableau<T, 13, 21>

Source

pub fn rkv878e() -> Self

A highly efficient Runge-Kutta (13:8(7)) pair method.

§Overview

This is an efficient 13-stage Runge-Kutta method that provides:

  • Primary order: 8
  • Embedded order: 7 (for error estimation)
  • Number of stages: 13 primary + 8 additional stages for interpolation
  • Dominant stage order: 4
§Efficiency

The method is particularly “efficient” with a local truncation error 2-norm:

T_92 ~ 0.000000011182

which appears to be the minimum possible for this type of method.

§Interpolation
  • The additional stages (c[14] through c[21]) allow computation of approximations at any point in the solution domain with order 8.

  • The interpolant has continuous derivatives.

  • Maximum 2-norm of the local truncation error over [0,1] for the interpolant:

    Ti_92 ~ 0.0000011515

  • The interpolation error has two local maximum values on [0,1].

§Source

Verner’s RKV878e

Source§

impl<T: Real> ButcherTableau<T, 16, 21>

Source

pub fn rkv988e() -> Self

A better efficient Runge-Kutta (16:8(9)) pair with 8th-order interpolation.

§Overview

This provides a 16-stage Runge-Kutta method with:

  • Primary order: 8
  • Embedded order: 9 (for error estimation)
  • Number of stages: 16 primary + 5 additional stages for interpolation
  • Dominant stage order: 5
§Efficiency

This method is particularly efficient with a local truncation error 2-norm:

T_10,2 ~ 0.00000034399

This 2-norm is smaller than previous (16:8,9) pairs and has two local maximum values on [0,1].

§Interpolation
  • Nodes c[17]=1 and c[18]-c[21] were selected to minimize interpolation error.

  • The interpolant has continuous derivatives.

  • Maximum 2-norm of the local truncation error over [0,1] for the interpolant:

    Ti_9,2 ~ 0.000003721

  • The interpolation error has three local maximum values on [0,1].

§Notes
  • Exact coefficients of the method are possible using surds in terms of sqrt(6).
  • This method shares the same primary stages as the 9th-order method (rkv989e).
  • The interpolant requires fewer additional stages than the 9th-order version.
§Source

Verner’s RKV878e

Source§

impl<T: Real> ButcherTableau<T, 16, 26>

Source

pub fn rkv989e() -> Self

A better efficient Runge-Kutta (16:9(8)) pair with 9th-order interpolation.

§Overview

This provides a 16-stage Runge-Kutta method with:

  • Primary order: 9
  • Embedded order: 8 (for error estimation)
  • Number of stages: 16 primary + 10 additional stages for interpolation
  • Dominant stage order: 5
§Efficiency

This method is particularly efficient with a local truncation error 2-norm:

T_10,2 ~ 0.00000034399

This 2-norm is smaller than previous (16:8,9) pairs and has two local maximum values on [0,1].

§Interpolation
  • The method provides a 9th-order interpolant using additional nodes c[17] through c[25].

  • The interpolant has continuous derivatives.

  • Maximum 2-norm of the local truncation error over [0,1] for the interpolant:

    Ti_10,2 ~ 0.0000007830

  • The interpolation error has two local maximum values on [0,1].

§Notes
  • Exact coefficients of the method are possible using surds in terms of sqrt(6).
  • Exact coefficients for the interpolants cannot be conveniently represented as they require c[18] to be a zero of a cubic polynomial with coefficients that are surds in terms of sqrt(6).
§Source

Verner’s RKV878e

Source§

impl<T: Real> ButcherTableau<T, 2>

Source

pub fn gauss_legendre_4() -> Self

Butcher Tableau for the Gauss-Legendre method of order 4.

§Overview

This provides a 2-stage, implicit Runge-Kutta method (Gauss-Legendre) with:

  • Primary order: 4
  • Number of stages: 2
§Notes
  • Gauss-Legendre methods are A-stable and symmetric.
  • They are highly accurate for their number of stages.
  • The c values are the roots of the Legendre polynomial P_s(2x-1) = 0.
  • This implementation includes two sets of b coefficients. The primary b coefficients are used for the solution, and bh can represent alternative coefficients (often related to error estimation or specific properties).
§Butcher Tableau
(1/2 - sqrt(3)/6) |  1/4                  1/4 - sqrt(3)/6
(1/2 + sqrt(3)/6) |  1/4 + sqrt(3)/6      1/4
-------------------|---------------------------------------
                   |  1/2                  1/2
                   |  1/2 + sqrt(3)/2      1/2 - sqrt(3)/2  (bh coefficients)
§References
  • Hairer, E., Nørsett, S. P., & Wanner, G. (1993). Solving Ordinary Differential Equations I: Nonstiff Problems. Springer. (Page 200, Table 4.5)
Source§

impl<T: Real> ButcherTableau<T, 3>

Source

pub fn gauss_legendre_6() -> Self

Butcher Tableau for the Gauss-Legendre method of order 6.

§Overview

This provides a 3-stage, implicit Runge-Kutta method (Gauss-Legendre) with:

  • Primary order: 6
  • Number of stages: 3
§Notes
  • Gauss-Legendre methods are A-stable and symmetric.
  • They are highly accurate for their number of stages.
  • The c values are the roots of the Legendre polynomial P_s(2x-1) = 0.
  • This implementation includes two sets of b coefficients. The primary b coefficients are used for the solution, and bh can represent alternative coefficients.
§Butcher Tableau
(1/2 - sqrt(15)/10) |  5/36                2/9 - sqrt(15)/15    5/36 - sqrt(15)/30
 1/2                |  5/36 + sqrt(15)/24  2/9                  5/36 - sqrt(15)/24
(1/2 + sqrt(15)/10) |  5/36 + sqrt(15)/30  2/9 + sqrt(15)/15    5/36
--------------------|-------------------------------------------------------------
                    |  5/18              4/9                  5/18
                    | -5/6               8/3                 -5/6                (bh coefficients)
§References
  • Hairer, E., Nørsett, S. P., & Wanner, G. (1993). Solving Ordinary Differential Equations I: Nonstiff Problems. Springer. (Page 200, Table 4.5)
Source§

impl<T: Real> ButcherTableau<T, 2>

Source

pub fn lobatto_iiic_2() -> Self

Butcher Tableau for the Lobatto IIIC method of order 2.

§Overview

This provides a 2-stage, implicit Runge-Kutta method (Lobatto IIIC) with:

  • Primary order: 2
  • Number of stages: 2
§Notes
  • Lobatto IIIC methods are L-stable.
  • They are algebraically stable and thus B-stable, making them suitable for stiff problems.
  • This implementation includes two sets of b coefficients. The primary b coefficients are used for the solution, and bh can represent alternative coefficients if needed (though their specific use here as a second row is characteristic of Lobatto IIIC).
§Butcher Tableau
0   |  1/2  -1/2
1   |  1/2   1/2
----|------------
    |  1/2   1/2   (b coefficients)
    |  1     0     (bh coefficients)
§References
  • Hairer, E., & Wanner, G. (1996). Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems. Springer. (Page 80, Table 5.7)
Source§

impl<T: Real> ButcherTableau<T, 3>

Source

pub fn lobatto_iiic_4() -> Self

Butcher Tableau for the Lobatto IIIC method of order 4.

§Overview

This provides a 3-stage, implicit Runge-Kutta method (Lobatto IIIC) with:

  • Primary order: 4
  • Number of stages: 3
§Notes
  • Lobatto IIIC methods are L-stable.
  • They are algebraically stable and thus B-stable, making them suitable for stiff problems.
  • This implementation includes two sets of b coefficients. The primary b coefficients are used for the solution, and bh represents the second row of coefficients from the tableau.
§Butcher Tableau
0   |  1/6  -1/3   1/6
1/2 |  1/6   5/12 -1/12
1   |  1/6   2/3   1/6
----|------------------
    |  1/6   2/3   1/6   (b coefficients)
    | -1/2   2    -1/2   (bh coefficients)
§References
  • Hairer, E., & Wanner, G. (1996). Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems. Springer. (Page 80, Table 5.7)
Source§

impl<T: Real> ButcherTableau<T, 2>

Source

pub fn radau_iia_3() -> Self

Butcher Tableau for the Radau IIA method of order 3.

§Overview

This provides a 2-stage, implicit Runge-Kutta method (Radau IIA) with:

  • Primary order: 3
  • Embedded order: None (this implementation does not include an embedded error estimate)
  • Number of stages: 2
§Interpolation
  • This standard Radau IIA order 3 implementation does not provide coefficients for dense output (interpolation).
§Notes
  • Radau IIA methods are A-stable and L-stable, making them suitable for stiff differential equations.
  • Being implicit, they generally require solving a system of algebraic equations at each step.
§Butcher Tableau
1/3 |  5/12  -1/12
1   |  3/4    1/4
----|---------------
    |  3/4    1/4
§References
  • Hairer, E., Nørsett, S. P., & Wanner, G. (1993). Solving Ordinary Differential Equations I: Nonstiff Problems. Springer.
  • Hairer, E., & Wanner, G. (1996). Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems. Springer.
Source§

impl<T: Real> ButcherTableau<T, 3>

Source

pub fn radau_iia_5() -> Self

Butcher Tableau for the Radau IIA method of order 5.

§Overview

This provides a 3-stage, implicit Runge-Kutta method (Radau IIA) with:

  • Primary order: 5
  • Embedded order: None (this implementation does not include an embedded error estimate)
  • Number of stages: 3
§Interpolation
  • This standard Radau IIA order 5 implementation does not provide coefficients for dense output (interpolation).
§Notes
  • Radau IIA methods are A-stable and L-stable, making them highly suitable for stiff differential equations.
  • Being implicit, they generally require solving a system of algebraic equations at each step.
  • The c values are the roots of P_s(2x-1) - P_{s-1}(2x-1) = 0 where P_s is the s-th Legendre polynomial. For s=3, the c values are (2/5 - sqrt(6)/10), (2/5 + sqrt(6)/10), and 1.
§Butcher Tableau
c1 | a11 a12 a13
c2 | a21 a22 a23
c3 | a31 a32 a33
---|------------
   | b1  b2  b3

Where: c1 = 2/5 - sqrt(6)/10 c2 = 2/5 + sqrt(6)/10 c3 = 1

a11 = 11/45 - 7sqrt(6)/360 a12 = 37/225 - 169sqrt(6)/1800 a13 = -2/225 + sqrt(6)/75

a21 = 37/225 + 169sqrt(6)/1800 a22 = 11/45 + 7sqrt(6)/360 a23 = -2/225 - sqrt(6)/75

a31 = 4/9 - sqrt(6)/36 a32 = 4/9 + sqrt(6)/36 a33 = 1/9

b1 = 4/9 - sqrt(6)/36 b2 = 4/9 + sqrt(6)/36 b3 = 1/9

§References
  • Hairer, E., Nørsett, S. P., & Wanner, G. (1993). Solving Ordinary Differential Equations I: Nonstiff Problems. Springer.
  • Hairer, E., & Wanner, G. (1996). Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems. Springer.
Source§

impl<T: Real> ButcherTableau<T, 2>

Source

pub fn sdirk21() -> Self

SDIRK-2-1: 2-stage, 2nd order SDIRK method with 1st order embedding

§Overview

This provides a 2-stage, singly diagonally implicit Runge-Kutta method with:

  • Primary order: 2
  • Embedded order: 1 (for error estimation)
  • Number of stages: 2
  • A-stable and B-stable
§Notes
  • This is a simple SDIRK method where all diagonal entries are equal (γ = 1)
  • Good for basic adaptive stepping with stiff problems
  • The embedded method provides basic error estimation for step size control
  • Particularly useful as a starting method for more complex stiff systems
§Butcher Tableau
1    | 1    0
0    | -1   1
-----|--------
     | 1/2  1/2  (2nd order)
     | 1    0    (1st order embedding)
§References
  • Hairer, E., Wanner, G. (1996). “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”
Source§

impl<T: Real> ButcherTableau<T, 3>

Source

pub fn esdirk33() -> Self

ESDIRK-3-3: 3-stage, 3rd order ESDIRK method

§Overview

This provides a 3-stage, explicit singly diagonally implicit Runge-Kutta method with:

  • Primary order: 3
  • Number of stages: 3
  • A-stable
§Notes
  • This method pairs with SSPRK(3,3)-Shu-Osher-ERK to make a 3rd order IMEX method
  • Has an explicit first stage (ESDIRK property) making it computationally efficient
  • The first stage being explicit reduces the computational cost per step
  • Suitable for problems with both stiff and non-stiff components
§Butcher Tableau
0      | 0      0      0
1      | 4γ+2β  1-4γ-2β 0
1/2    | α₃₁    γ      β
-------|------------------
       | 1/6    1/6    2/3

where:

  • β = √3/6 + 1/2 ≈ 0.7886751346
  • γ = (-1/8)(√3 + 1) ≈ -0.3416407865
  • α₃₁ = 1/2 - β - γ ≈ 0.0529656519
§References
  • Conde, S., et al. (2017). “Implicit and implicit-explicit strong stability preserving Runge-Kutta methods”
Source§

impl<T: Real> ButcherTableau<T, 4>

Source

pub fn esdirk324l2sa() -> Self

ESDIRK3(2)4L[2]SA: 4-stage, 3rd order ESDIRK method with embedded 2nd order

§Overview

This provides a 4-stage, explicit singly diagonally implicit Runge-Kutta method with:

  • Primary order: 3
  • Embedded order: 2 (for error estimation)
  • Number of stages: 4
  • A-stable and B-stable
§References
  • Kennedy, C.A. and Carpenter, M.H. (2003). “Additive Runge-Kutta schemes for convection-diffusion-reaction equations”
Source§

impl<T: Real> ButcherTableau<T, 4>

Source

pub fn kvaerno423() -> Self

Kvaerno(4,2,3): 4-stage, 3rd order DIRK method with embedded 2nd order

§Overview

This provides a 4-stage, diagonally implicit Runge-Kutta method with:

  • Primary order: 3
  • Embedded order: 2 (for error estimation)
  • Number of stages: 4
  • A-stable
§Notes
  • Developed specifically for stiff differential equations
  • A-stable method suitable for moderately stiff problems
  • The embedded method provides efficient error estimation for adaptive stepping
  • All diagonal elements are identical (γ ≈ 0.4358665215), simplifying LU factorization
  • Good balance between stability and computational efficiency
§Butcher Tableau
0      | 0      0      0      0
.8717  | .4359  .4359  0      0
1      | .4906  .0736  .4359  0  
1      | .3088  1.4906 -1.2352 .4359
-------|-------------------------
b³     | .3088  1.4906 -1.2352 .4359
b²     | .4906  .0736  .4359  0

where γ ≈ 0.4358665215

§References
  • Kvaerno, A. (2004). “Singly diagonally implicit Runge-Kutta methods with an explicit first stage”
Source§

impl<T: Real> ButcherTableau<T, 7>

Source

pub fn kvaerno745() -> Self

Kvaerno(7,4,5): 7-stage, 5th order DIRK method with embedded 4th order

§Overview

This provides a 7-stage, diagonally implicit Runge-Kutta method with:

  • Primary order: 5
  • Embedded order: 4 (for error estimation)
  • Number of stages: 7
  • A-stable and B-stable
§Notes
  • High-order method designed for stiff differential equations requiring high accuracy
  • A-stable and B-stable for excellent stability properties
  • The embedded 4th order method provides accurate error estimation for adaptive control
  • All diagonal elements are identical (γ = 0.26), enabling efficient LU factorization reuse
  • Particularly effective for smooth solutions where high accuracy is required
  • Higher computational cost per step but fewer steps needed due to high order
§Butcher Tableau
0      | 0      0      0      0      0      0      0
.52    | .26    .26    0      0      0      0      0
1.2303 | .13    .8403  .26    0      0      0      0
.8958  | .2237  .4768  -.0647 .26    0      0      0
.4364  | .1665  .1045  .0363  -.1309 .26    0      0
1      | .1386  .0000  -.0425 .0245  .6194  .26    0
1      | .1366  .0000  -.0550 -.0412 .6299  .0696  .26
-------|--------------------------------------------
b⁵     | .1366  .0000  -.0550 -.0412 .6299  .0696  .26
b⁴     | .1386  .0000  -.0425 .0245  .6194  .26    0

where γ = 0.26 exactly

§References
  • Kvaerno, A. (2004). “Singly diagonally implicit Runge-Kutta methods with an explicit first stage”
Source§

impl<T: Real, const S: usize, const I: usize> ButcherTableau<T, S, I>

Source

pub const STAGES: usize = S

Number of stages in the method

Source

pub const EXTRA_STAGES: usize

Number of extra stages for interpolation

Auto Trait Implementations§

§

impl<T, const S: usize, const I: usize> Freeze for ButcherTableau<T, S, I>
where T: Freeze,

§

impl<T, const S: usize, const I: usize> RefUnwindSafe for ButcherTableau<T, S, I>
where T: RefUnwindSafe,

§

impl<T, const S: usize, const I: usize> Send for ButcherTableau<T, S, I>

§

impl<T, const S: usize, const I: usize> Sync for ButcherTableau<T, S, I>

§

impl<T, const S: usize, const I: usize> Unpin for ButcherTableau<T, S, I>
where T: Unpin,

§

impl<T, const S: usize, const I: usize> UnwindSafe for ButcherTableau<T, S, I>
where T: UnwindSafe,

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> Same for T

Source§

type Output = T

Should always be Self
Source§

impl<SS, SP> SupersetOf<SS> for SP
where SS: SubsetOf<SP>,

Source§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
Source§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
Source§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
Source§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.