OrdinaryDiffEqVerner

Verner methods are high-order explicit Runge-Kutta methods designed for high-accuracy integration of non-stiff differential equations. These are the preferred solvers when very low tolerances are required.

Key Properties

Verner methods provide:

  • High-order accuracy (6th through 9th order) for precise integration
  • Excellent efficiency at low tolerances (1e-8 to 1e-15)
  • Robust error estimation with embedded error control
  • Dense output capability with high-quality interpolation

When to Use Verner Methods

Verner methods are recommended for:

  • High-accuracy requirements with tolerances between 1e-8 and 1e-15
  • Smooth non-stiff problems where high precision is critical
  • Long-time integration where error accumulation must be minimized
  • Problems requiring dense output with high interpolation accuracy
  • Orbit computation, molecular dynamics, and other precision-critical applications

Solver Selection Guide

Medium-low tolerance (1e-6 to 1e-8)

  • Vern6: Sixth-order method, good balance of efficiency and accuracy
  • AutoVern6: Automatic switching version for mixed stiffness

Low tolerance (1e-8 to 1e-12) with Float64

  • Vern7: Seventh-order method, excellent for most high-precision needs
  • Vern8: Eighth-order method, best efficiency at very low tolerances
  • AutoVern7, AutoVern8: Automatic switching versions

Very low tolerance (<1e-12)

  • Vern9: Ninth-order method for extreme precision requirements
  • Recommended with BigFloat for tolerances below 1e-15
  • AutoVern9: Automatic switching version for mixed problems

Performance Notes

  • Vern6: Most efficient for tolerances around 1e-6 to 1e-8
  • Vern7: Sweet spot for tolerances around 1e-8 to 1e-10
  • Vern8: Best for tolerances around 1e-10 to 1e-12
  • Vern9: For tolerances below 1e-12, especially with arbitrary precision

The Auto* variants automatically switch to stiff solvers when stiffness is detected, making them robust for problems of unknown character.

Installation

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

using Pkg
Pkg.add("OrdinaryDiffEqVerner")

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 OrdinaryDiffEqVerner

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

Full list of solvers

OrdinaryDiffEqVerner.Vern6Type
Vern6(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False(),
        lazy = true)

Explicit Runge-Kutta Method. Verner's most efficient 6/5 method (lazy 6th 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{verner2010numerically, title={Numerically optimal Runge–Kutta pairs with interpolants}, author={Verner, James H}, journal={Numerical Algorithms}, volume={53}, number={2-3}, pages={383–396}, year={2010}, publisher={Springer} }

source
OrdinaryDiffEqVerner.Vern7Type
Vern7(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False(),
        lazy = true)

Explicit Runge-Kutta Method. Verner's most efficient 7/6 method (lazy 7th order interpolant). Good for problems requiring high accuracy. Slightly more computationally expensive than Tsit5. Performance best when parameter vector remains unchanged. Recommended for high-accuracy non-stiff 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.
  • lazy: determines if the lazy interpolant is used.

References

@article{verner2010numerically, title={Numerically optimal Runge–Kutta pairs with interpolants}, author={Verner, James H}, journal={Numerical Algorithms}, volume={53}, number={2-3}, pages={383–396}, year={2010}, publisher={Springer} }

source
OrdinaryDiffEqVerner.Vern8Type
Vern8(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False(),
        lazy = true)

Explicit Runge-Kutta Method. Verner's most efficient 8/7 method (lazy 8th 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{verner2010numerically, title={Numerically optimal Runge–Kutta pairs with interpolants}, author={Verner, James H}, journal={Numerical Algorithms}, volume={53}, number={2-3}, pages={383–396}, year={2010}, publisher={Springer} }

source
OrdinaryDiffEqVerner.Vern9Type
Vern9(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!,
        thread = OrdinaryDiffEq.False(),
        lazy = true)

Explicit Runge-Kutta Method. Verner's most efficient 9/8 method (lazy 9th 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{verner2010numerically, title={Numerically optimal Runge–Kutta pairs with interpolants}, author={Verner, James H}, journal={Numerical Algorithms}, volume={53}, number={2-3}, pages={383–396}, year={2010}, publisher={Springer} }

source
OrdinaryDiffEqVerner.AutoVern6Function

Automatic switching algorithm that can switch between the (non-stiff) Vern6() and stiff_alg.

AutoVern6(stiff_alg; kwargs...)

This method is equivalent to AutoAlgSwitch(Vern6(), stiff_alg; kwargs...). To gain access to stiff algorithms you might have to install additional libraries, such as OrdinaryDiffEqRosenbrock.

source
OrdinaryDiffEqVerner.AutoVern7Function

Automatic switching algorithm that can switch between the (non-stiff) Vern7() and stiff_alg.

AutoVern7(stiff_alg; kwargs...)

This method is equivalent to AutoAlgSwitch(Vern7(), stiff_alg; kwargs...). To gain access to stiff algorithms you might have to install additional libraries, such as OrdinaryDiffEqRosenbrock.

source
OrdinaryDiffEqVerner.AutoVern8Function

Automatic switching algorithm that can switch between the (non-stiff) Vern8() and stiff_alg.

AutoVern8(stiff_alg; kwargs...)

This method is equivalent to AutoAlgSwitch(Vern8(), stiff_alg; kwargs...). To gain access to stiff algorithms you might have to install additional libraries, such as OrdinaryDiffEqRosenbrock.

source
OrdinaryDiffEqVerner.AutoVern9Function

Automatic switching algorithm that can switch between the (non-stiff) Vern9() and stiff_alg.

AutoVern9(stiff_alg; kwargs...)

This method is equivalent to AutoAlgSwitch(Vern9(), stiff_alg; kwargs...). To gain access to stiff algorithms you might have to install additional libraries, such as OrdinaryDiffEqRosenbrock.

source