OrdinaryDiffEqSSPRK

Strong Stability Preserving Runge-Kutta (SSPRK) methods are specialized explicit Runge-Kutta methods designed to preserve important stability properties of the underlying spatial discretization when applied to hyperbolic partial differential equations and conservation laws.

Key Properties

SSPRK methods provide:

  • Strong stability preservation for convex functionals (total variation, maximum norm, entropy)
  • Optimal SSP coefficients allowing larger stable timesteps
  • Non-oscillatory behavior crucial for hyperbolic PDEs and conservation laws
  • High-order accuracy while maintaining monotonicity properties
  • Specialized variants for different orders and storage requirements

When to Use SSPRK Methods

SSPRK methods are essential for:

  • Hyperbolic partial differential equations (Euler equations, shallow water, etc.)
  • Conservation laws where preserving physical bounds is critical
  • Discontinuous Galerkin methods and other high-order spatial discretizations
  • Problems requiring monotonicity preservation or total variation stability
  • Shock-capturing schemes where spurious oscillations must be avoided
  • Astrophysical simulations and computational fluid dynamics

SSP Coefficient and CFL Conditions

The SSP coefficient determines the maximum allowable timestep for stability preservation. Use OrdinaryDiffEqCore.ssp_coefficient(alg) to query this value for step size calculations. The timestep must satisfy dt ≤ CFL * dx / max_wavespeed where CFL ≤ SSP coefficient.

Solver Selection Guide

Second-order methods

  • SSPRK22: Two-stage, second-order (SSP coefficient = 1)
  • SSPRKMSVS32: Three-step multistep variant

Third-order methods

  • SSPRK33: Three-stage, third-order, optimal (SSP coefficient = 1)
  • SSPRK53: Five-stage, third-order, higher SSP coefficient
  • SSPRK63, SSPRK73, SSPRK83: More stages for larger SSP coefficients
  • SSPRK43: Four-stage with embedded error estimation
  • SSPRK432: Low-storage variant

Fourth-order methods

  • SSPRK54: Five-stage, fourth-order
  • SSPRK104: Ten-stage, fourth-order, large SSP coefficient

Low-storage variants

  • SSPRK53_2N1, SSPRK53_2N2, SSPRK53_H: Two-register storage schemes

Discontinuous Galerkin optimized

  • KYKSSPRK42: Optimized for DG spatial discretizations
  • KYK2014DGSSPRK_3S2: Specialized DG method

Adaptive methods

  • SSPRK432: Third-order with error control
  • SSPRK932: High-stage adaptive method
  • SSPRKMSVS43: Multistep adaptive variant

Installation

To be able to access the solvers in OrdinaryDiffEqSSPRK, you must first install them use the Julia package manager:

using Pkg
Pkg.add("OrdinaryDiffEqSSPRK")

This will only install the solvers listed at the bottom of this page. If you want to explore other solvers for your problem, you will need to install some of the other libraries listed in the navigation bar on the left.

Example usage

using OrdinaryDiffEqSSPRK

function lorenz!(du, u, p, t)
    du[1] = 10.0 * (u[2] - u[1])
    du[2] = u[1] * (28.0 - u[3]) - u[2]
    du[3] = u[1] * u[2] - (8 / 3) * u[3]
end
u0 = [1.0; 0.0; 0.0]
tspan = (0.0, 100.0)
prob = ODEProblem(lorenz!, u0, tspan)
sol = solve(prob, SSPRK22())

Full list of solvers

OrdinaryDiffEqSSPRK.SSPRK22Type
SSPRK22(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A second-order, two-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Shu, Chi-Wang, and Stanley Osher. Efficient implementation of essentially non-oscillatory shock-capturing schemes. Journal of Computational Physics 77.2 (1988): 439-471. https://doi.org/10.1016/0021-9991(88)90177-5

source
OrdinaryDiffEqSSPRK.SSPRK33Type
SSPRK33(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, three-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Shu, Chi-Wang, and Stanley Osher. Efficient implementation of essentially non-oscillatory shock-capturing schemes. Journal of Computational Physics 77.2 (1988): 439-471. https://doi.org/10.1016/0021-9991(88)90177-5

source
OrdinaryDiffEqSSPRK.SSPRK53Type
SSPRK53(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, five-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Ruuth, Steven. Global optimization of explicit strong-stability-preserving Runge-Kutta methods. Mathematics of Computation 75.253 (2006): 183-207

source
OrdinaryDiffEqSSPRK.KYKSSPRK42Type
KYKSSPRK42(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
             step_limiter! = OrdinaryDiffEq.trivial_limiter!,
             thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Optimal strong-stability-preserving Runge-Kutta time discretizations for discontinuous Galerkin methods

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

@article{kubatko2014optimal, title={Optimal strong-stability-preserving Runge–Kutta time discretizations for discontinuous Galerkin methods}, author={Kubatko, Ethan J and Yeager, Benjamin A and Ketcheson, David I}, journal={Journal of Scientific Computing}, volume={60}, pages={313–344}, year={2014}, publisher={Springer}}

source
OrdinaryDiffEqSSPRK.KYK2014DGSSPRK_3S2Type
KYK2014DGSSPRK_3S2(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
                     step_limiter! = OrdinaryDiffEq.trivial_limiter!,
                     thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Optimal strong-stability-preserving Runge-Kutta time discretizations for discontinuous Galerkin methods

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

@article{kubatko2014optimal, title={Optimal strong-stability-preserving Runge–Kutta time discretizations for discontinuous Galerkin methods}, author={Kubatko, Ethan J and Yeager, Benjamin A and Ketcheson, David I}, journal={Journal of Scientific Computing}, volume={60}, pages={313–344}, year={2014}, publisher={Springer}}

source
OrdinaryDiffEqSSPRK.SSPRK53_2N1Type
SSPRK53_2N1(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, five-stage explicit strong stability preserving (SSP) low-storage method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Higueras and T. Roldán. New third order low-storage SSP explicit Runge–Kutta methods arXiv:1809.04807v1.

source
OrdinaryDiffEqSSPRK.SSPRK53_2N2Type
SSPRK53_2N2(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, five-stage explicit strong stability preserving (SSP) low-storage method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Higueras and T. Roldán. New third order low-storage SSP explicit Runge–Kutta methods arXiv:1809.04807v1.

source
OrdinaryDiffEqSSPRK.SSPRK53_HType
SSPRK53_H(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, five-stage explicit strong stability preserving (SSP) low-storage method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Higueras and T. Roldán. New third order low-storage SSP explicit Runge–Kutta methods arXiv:1809.04807v1.

source
OrdinaryDiffEqSSPRK.SSPRK63Type
SSPRK63(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, six-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Ruuth, Steven. Global optimization of explicit strong-stability-preserving Runge-Kutta methods. Mathematics of Computation 75.253 (2006): 183-207

source
OrdinaryDiffEqSSPRK.SSPRK73Type
SSPRK73(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, seven-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Ruuth, Steven. Global optimization of explicit strong-stability-preserving Runge-Kutta methods. Mathematics of Computation 75.253 (2006): 183-207

source
OrdinaryDiffEqSSPRK.SSPRK83Type
SSPRK83(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, eight-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Ruuth, Steven. Global optimization of explicit strong-stability-preserving Runge-Kutta methods. Mathematics of Computation 75.253 (2006): 183-207

source
OrdinaryDiffEqSSPRK.SSPRK43Type
SSPRK43(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, four-stage explicit strong stability preserving (SSP) method.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Optimal third-order explicit SSP method with four stages discovered by

  • J. F. B. M. Kraaijevanger. "Contractivity of Runge-Kutta methods." In: BIT Numerical Mathematics 31.3 (1991), pp. 482–528. DOI: 10.1007/BF01933264.

Embedded method constructed by

  • Sidafa Conde, Imre Fekete, John N. Shadid. "Embedded error estimation and adaptive step-size control for optimal explicit strong stability preserving Runge–Kutta methods." arXiv: 1806.08693

Efficient implementation (and optimized controller) developed by

  • Hendrik Ranocha, Lisandro Dalcin, Matteo Parsani, David I. Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEqSSPRK.SSPRK432Type
SSPRK432(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, four-stage explicit strong stability preserving (SSP) method.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Gottlieb, Sigal, David I. Ketcheson, and Chi-Wang Shu. Strong stability preserving Runge-Kutta and multistep time discretizations. World Scientific, 2011. Example 6.1

source
OrdinaryDiffEqSSPRK.SSPRKMSVS43Type
SSPRKMSVS43(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, four-step explicit strong stability preserving (SSP) linear multistep method. This method does not come with an error estimator and requires a fixed time step size.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Shu, Chi-Wang. Total-variation-diminishing time discretizations. SIAM Journal on Scientific and Statistical Computing 9, no. 6 (1988): 1073-1084. DOI: 10.1137/0909073

source
OrdinaryDiffEqSSPRK.SSPRKMSVS32Type
SSPRKMSVS32(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A second-order, three-step explicit strong stability preserving (SSP) linear multistep method. This method does not come with an error estimator and requires a fixed time step size.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Shu, Chi-Wang. Total-variation-diminishing time discretizations. SIAM Journal on Scientific and Statistical Computing 9, no. 6 (1988): 1073-1084. DOI: 10.1137/0909073

source
OrdinaryDiffEqSSPRK.SSPRK932Type
SSPRK932(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A third-order, nine-stage explicit strong stability preserving (SSP) method.

Consider using SSPRK43 instead, which uses the same main method and an improved embedded method.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Gottlieb, Sigal, David I. Ketcheson, and Chi-Wang Shu. Strong stability preserving Runge-Kutta and multistep time discretizations. World Scientific, 2011.

source
OrdinaryDiffEqSSPRK.SSPRK54Type
SSPRK54(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A fourth-order, five-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Ruuth, Steven. Global optimization of explicit strong-stability-preserving Runge-Kutta methods. Mathematics of Computation 75.253 (2006): 183-207.

source
OrdinaryDiffEqSSPRK.SSPRK104Type
SSPRK104(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. A fourth-order, ten-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Keyword Arguments

  • stage_limiter!: function of the form limiter!(u, integrator, p, t)
  • step_limiter!: function of the form limiter!(u, integrator, p, t)
  • thread: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Ketcheson, David I. Highly efficient strong stability-preserving Runge–Kutta methods with low-storage implementations. SIAM Journal on Scientific Computing 30.4 (2008): 2113-2136.

source