OrdinaryDiffEqLowOrderRK

Low-order explicit Runge-Kutta methods for non-stiff differential equations. Most of the time, you should use Tsit5, which is the most common and efficient low-order RK method. The alternative methods provided here are for special circumstances where Tsit5 is not suitable.

Key Properties

Low-order explicit RK methods offer:

  • Computational efficiency at higher tolerances (>1e-6)
  • Robust error control for difficult non-stiff problems
  • Specialized interpolation properties for applications requiring dense output
  • Lower-order derivatives requirements for non-smooth functions
  • Good performance for specific problem types

When to Use Alternative Low-Order RK Methods

Choose these methods instead of Tsit5 when:

  • ODE function f is not differentiable to 5th order - use lower-order methods (the more discontinuous, the lower the order needed)
  • Heavy use of interpolations - OwrenZen methods have superior interpolation convergence
  • Delay differential equations - OwrenZen methods are most efficient (see SciMLBenchmarks)
  • Very high tolerances (>1e-3) - BS3 is more efficient than Tsit5
  • Quadratic polynomial ODEs - SIR54 is optimized for these systems
  • Educational purposes - simpler methods for understanding algorithms

Solver Selection Guide

Primary recommendation

For most problems, use Tsit5 instead of these methods.

High tolerances (>1e-3)

  • BS3: Third-order Bogacki-Shampine method, most efficient for very high tolerances

Superior interpolation needs

  • OwrenZen3: Third-order with excellent interpolation convergence
  • OwrenZen5: Fifth-order with excellent interpolation, optimal for DDEs
  • OwrenZen4: Fourth-order interpolation-optimized method

Non-smooth or discontinuous ODEs

  • BS3: Third-order for mildly non-smooth functions
  • Heun: Second-order for more discontinuous functions (not generally recommended)
  • Euler: First-order for highly discontinuous problems

Robust error control alternatives

  • BS5: Fifth-order with very robust error estimation
  • DP5: Fifth-order Dormand-Prince method, classical alternative to Tsit5

Specialized applications

  • RK4: Fourth-order with special residual error control, good for DDEs. Note: Uses adaptive timestepping by default - set adaptive=false in solve() for traditional fixed-step RK4
  • SIR54: Fifth-order optimized for ODEs defined by quadratic polynomials (e.g., SIR-type epidemiological models)
  • Stepanov5: Fifth-order method with enhanced stability properties and optimized error constants
  • Ralston: Second-order with optimized error constants

Periodic and oscillatory problems

  • Anas5: Fifth-order optimized for periodic problems with minimal phase error
  • FRK65: Sixth-order zero dissipation method for oscillatory problems

Advanced specialized methods

  • RKO65: Sixth-order optimized method
  • MSRK5, MSRK6: Multi-stage methods for specific applications
  • PSRK4p7q6, PSRK3p5q4, PSRK3p6q5: Pseudo-symplectic methods
  • Alshina2, Alshina3, Alshina6: Methods with optimized parameters

Installation

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

using Pkg
Pkg.add("OrdinaryDiffEqLowOrderRK")

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 OrdinaryDiffEqLowOrderRK

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, BS3())

Full list of solvers

OrdinaryDiffEqLowOrderRK.EulerType
Euler()

Explicit Runge-Kutta Method. The canonical forward Euler method. Fixed timestep only.

Keyword Arguments

References

E. Hairer, S.P. Norsett, G. Wanner, (1993) Solving Ordinary Differential Equations I. Nonstiff Problems. 2nd Edition. Springer Series in Computational Mathematics, Springer-Verlag.

source
OrdinaryDiffEqLowOrderRK.HeunType
Heun(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
       step_limiter! = OrdinaryDiffEq.trivial_limiter!,
       thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. The second order Heun's method. Uses embedded Euler method for adaptivity.

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

E. Hairer, S.P. Norsett, G. Wanner, (1993) Solving Ordinary Differential Equations I. Nonstiff Problems. 2nd Edition. Springer Series in Computational Mathematics, Springer-Verlag.

source
OrdinaryDiffEqLowOrderRK.RalstonType
Ralston(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. The optimized second order midpoint method. Uses embedded Euler method for adaptivity.

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

E. Hairer, S.P. Norsett, G. Wanner, (1993) Solving Ordinary Differential Equations I. Nonstiff Problems. 2nd Edition. Springer Series in Computational Mathematics, Springer-Verlag.

source
OrdinaryDiffEqLowOrderRK.MidpointType
Midpoint(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. The second order midpoint method. Uses embedded Euler method for adaptivity.

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

E. Hairer, S.P. Norsett, G. Wanner, (1993) Solving Ordinary Differential Equations I. Nonstiff Problems. 2nd Edition. Springer Series in Computational Mathematics, Springer-Verlag.

source
OrdinaryDiffEqLowOrderRK.RK4Type
RK4(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
      step_limiter! = OrdinaryDiffEq.trivial_limiter!,
      thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. The canonical Runge-Kutta Order 4 method. Uses a defect control for adaptive stepping using maximum error over the whole interval. Classic fourth-order method. Good for medium accuracy calculations.

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{shampine2005solving, title={Solving ODEs and DDEs with residual control}, author={Shampine, LF}, journal={Applied Numerical Mathematics}, volume={52}, number={1}, pages={113–127}, year={2005}, publisher={Elsevier} }

source
OrdinaryDiffEqLowOrderRK.BS3Type
BS3(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
      step_limiter! = OrdinaryDiffEq.trivial_limiter!,
      thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Bogacki-Shampine 3/2 method. Third-order adaptive method using embedded Euler method for adaptivity. Recommended for non-stiff problems at moderate tolerances.

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{bogacki19893, title={A 3 (2) pair of Runge-Kutta formulas}, author={Bogacki, Przemyslaw and Shampine, Lawrence F}, journal={Applied Mathematics Letters}, volume={2}, number={4}, pages={321–325}, year={1989}, publisher={Elsevier} }

source
OrdinaryDiffEqLowOrderRK.OwrenZen3Type
OwrenZen3(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Owren-Zennaro optimized interpolation 3/2 method (free 3rd order interpolant).

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{owren1992derivation, title={Derivation of efficient, continuous, explicit Runge–Kutta methods}, author={Owren, Brynjulf and Zennaro, Marino}, journal={SIAM journal on scientific and statistical computing}, volume={13}, number={6}, pages={1488–1501}, year={1992}, publisher={SIAM} }

source
OrdinaryDiffEqLowOrderRK.OwrenZen4Type
OwrenZen4(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Owren-Zennaro optimized interpolation 4/3 method (free 4th order interpolant).

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{owren1992derivation, title={Derivation of efficient, continuous, explicit Runge–Kutta methods}, author={Owren, Brynjulf and Zennaro, Marino}, journal={SIAM journal on scientific and statistical computing}, volume={13}, number={6}, pages={1488–1501}, year={1992}, publisher={SIAM} }

source
OrdinaryDiffEqLowOrderRK.OwrenZen5Type
OwrenZen5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Owren-Zennaro optimized interpolation 5/4 method (free 5th order interpolant).

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{owren1992derivation, title={Derivation of efficient, continuous, explicit Runge–Kutta methods}, author={Owren, Brynjulf and Zennaro, Marino}, journal={SIAM journal on scientific and statistical computing}, volume={13}, number={6}, pages={1488–1501}, year={1992}, publisher={SIAM} }

source
OrdinaryDiffEqLowOrderRK.BS5Type
BS5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
      step_limiter! = OrdinaryDiffEq.trivial_limiter!,
      thread = OrdinaryDiffEq.False(),
      lazy = true)

Explicit Runge-Kutta Method. Bogacki-Shampine 5/4 Runge-Kutta method. (lazy 5th order interpolant).

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.
  • lazy: determines if the lazy interpolant is used.

References

@article{bogacki1996efficient, title={An efficient runge-kutta (4, 5) pair}, author={Bogacki, P and Shampine, Lawrence F}, journal={Computers \& Mathematics with Applications}, volume={32}, number={6}, pages={15–28}, year={1996}, publisher={Elsevier} }

source
OrdinaryDiffEqLowOrderRK.DP5Type
DP5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
      step_limiter! = OrdinaryDiffEq.trivial_limiter!,
      thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Dormand-Prince's 5/4 Runge-Kutta method. (free 4th order interpolant).

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{dormand1980family, title={A family of embedded Runge-Kutta formulae}, author={Dormand, John R and Prince, Peter J}, journal={Journal of computational and applied mathematics}, volume={6}, number={1}, pages={19–26}, year={1980}, publisher={Elsevier} }

source
OrdinaryDiffEqLowOrderRK.Anas5Type
Anas5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False(),
        w = 1)

Explicit Runge-Kutta Method. 4th order Runge-Kutta method designed for periodic problems.

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.
  • w: a periodicity estimate, which when accurate the method becomes 5th order

(and is otherwise 4th order with less error for better estimates).

References

@article{anastassi2005optimized, title={An optimized Runge–Kutta method for the solution of orbital problems}, author={Anastassi, ZA and Simos, TE}, journal={Journal of Computational and Applied Mathematics}, volume={175}, number={1}, pages={1–9}, year={2005}, publisher={Elsevier}}

source
OrdinaryDiffEqLowOrderRK.RKO65Type
RKO65(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Tsitouras' Runge-Kutta-Oliver 6 stage 5th order 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

Tsitouras, Ch. "Explicit Runge–Kutta methods for starting integration of Lane–Emden problem." Applied Mathematics and Computation 354 (2019): 353-364. doi: https://doi.org/10.1016/j.amc.2019.02.047

source
OrdinaryDiffEqLowOrderRK.FRK65Type
FRK65(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False(),
        omega = 0.0)

Explicit Runge-Kutta Method. Zero Dissipation Runge-Kutta of 6th order.

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.
  • omega: a periodicity phase estimate,

when accurate this method results in zero numerical dissipation.

References

@article{medvedev2018fitted, title={Fitted modifications of Runge-Kutta pairs of orders 6 (5)}, author={Medvedev, Maxim A and Simos, TE and Tsitouras, Ch}, journal={Mathematical Methods in the Applied Sciences}, volume={41}, number={16}, pages={6184–6194}, year={2018}, publisher={Wiley Online Library}}

source
OrdinaryDiffEqLowOrderRK.RKMType
RKM(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
      step_limiter! = OrdinaryDiffEq.trivial_limiter!,
      thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. Method designed to have good stability properties when applied to pseudospectral discretizations of hyperbolic partial differential equaitons.

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{mead1999optimal, title={Optimal Runge–Kutta methods for first order pseudospectral operators}, author={Mead, JL and Renaut, RA}, journal={Journal of Computational Physics}, volume={152}, number={1}, pages={404–419}, year={1999}, publisher={Elsevier} }

source
OrdinaryDiffEqLowOrderRK.MSRK5Type
MSRK5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 5th order 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

Misha Stepanov - https://arxiv.org/pdf/2202.08443.pdf : Figure 3.

source
OrdinaryDiffEqLowOrderRK.MSRK6Type
MSRK6(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 6th order 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

Misha Stepanov - https://arxiv.org/pdf/2202.08443.pdf : Table4

source
OrdinaryDiffEqLowOrderRK.PSRK4p7q6Type
PSRK4p7q6(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 6-stage Pseudo-Symplectic 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

@article{Aubry1998, author = {A. Aubry and P. Chartier}, journal = {BIT Numer. Math.}, title = {Pseudo-symplectic {R}unge-{K}utta methods}, volume = {38}, PAGES = {439-461}, year = {1998}, }, @article{Capuano2017, title = {Explicit {R}unge–{K}utta schemes for incompressible flow with improved energy-conservation properties}, journal = {J. Comput. Phys.}, volume = {328}, pages = {86-94}, year = {2017}, issn = {0021-9991}, doi = {https://doi.org/10.1016/j.jcp.2016.10.040}, author = {F. Capuano and G. Coppola and L. Rández and L. {de Luca}},}

source
OrdinaryDiffEqLowOrderRK.PSRK3p5q4Type
PSRK3p5q4(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 4-stage Pseudo-Symplectic 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

@article{Aubry1998, author = {A. Aubry and P. Chartier}, journal = {BIT Numer. Math.}, title = {Pseudo-symplectic {R}unge-{K}utta methods}, year = {1998}, }, @article{Capuano2017, title = {Explicit {R}unge–{K}utta schemes for incompressible flow with improved energy-conservation properties}, journal = {J. Comput. Phys.}, year = {2017}, author = {F. Capuano and G. Coppola and L. Rández and L. {de Luca}},}

source
OrdinaryDiffEqLowOrderRK.PSRK3p6q5Type
PSRK3p6q5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 5-stage Pseudo-Symplectic 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

@article{Aubry1998, author = {A. Aubry and P. Chartier}, journal = {BIT Numer. Math.}, title = {Pseudo-symplectic {R}unge-{K}utta methods}, year = {1998}, }, @article{Capuano2017, title = {Explicit {R}unge–{K}utta schemes for incompressible flow with improved energy-conservation properties}, journal = {J. Comput. Phys.}, year = {2017}, author = {F. Capuano and G. Coppola and L. Rández and L. {de Luca}},}

source
OrdinaryDiffEqLowOrderRK.Stepanov5Type
Stepanov5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 5th order 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

@article{Stepanov2021Embedded5, title={Embedded (4, 5) pairs of explicit 7-stage Runge–Kutta methods with FSAL property}, author={Misha Stepanov}, journal={Calcolo}, year={2021}, volume={59} }

source
OrdinaryDiffEqLowOrderRK.SIR54Type
SIR54(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 5th order method suited for SIR-type epidemic models.

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{Kovalnogov2020RungeKuttaPS, title={Runge–Kutta pairs suited for SIR‐type epidemic models}, author={Vladislav N. Kovalnogov and Theodore E. Simos and Ch. Tsitouras}, journal={Mathematical Methods in the Applied Sciences}, year={2020}, volume={44}, pages={5210 - 5216} }

source
OrdinaryDiffEqLowOrderRK.Alshina2Type
Alshina2(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 2nd order, 2-stage Method with optimal parameters.

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{Alshina2008, doi = {10.1134/s0965542508030068}, url = {https://doi.org/10.1134/s0965542508030068}, year = {2008}, month = mar, publisher = {Pleiades Publishing Ltd}, volume = {48}, number = {3}, pages = {395–405}, author = {E. A. Alshina and E. M. Zaks and N. N. Kalitkin}, title = {Optimal first- to sixth-order accurate Runge-Kutta schemes}, journal = {Computational Mathematics and Mathematical Physics} }

source
OrdinaryDiffEqLowOrderRK.Alshina3Type
Alshina3(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 3rd order, 3-stage Method with optimal parameters.

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{Alshina2008, doi = {10.1134/s0965542508030068}, url = {https://doi.org/10.1134/s0965542508030068}, year = {2008}, month = mar, publisher = {Pleiades Publishing Ltd}, volume = {48}, number = {3}, pages = {395–405}, author = {E. A. Alshina and E. M. Zaks and N. N. Kalitkin}, title = {Optimal first- to sixth-order accurate Runge-Kutta schemes}, journal = {Computational Mathematics and Mathematical Physics} }

source
OrdinaryDiffEqLowOrderRK.Alshina6Type
Alshina6(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

Explicit Runge-Kutta Method. 6th order, 7-stage Method with optimal parameters.

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{Alshina2008, doi = {10.1134/s0965542508030068}, url = {https://doi.org/10.1134/s0965542508030068}, year = {2008}, month = mar, publisher = {Pleiades Publishing Ltd}, volume = {48}, number = {3}, pages = {395–405}, author = {E. A. Alshina and E. M. Zaks and N. N. Kalitkin}, title = {Optimal first- to sixth-order accurate Runge-Kutta schemes}, journal = {Computational Mathematics and Mathematical Physics} }

source