ODEInterfaceDiffEq.jl
ODEInterfaceDiffEq.jl provides a wrapper for classic Fortran ODE solvers from the ODEInterface.jl package into the SciML interface. These include the Hairer and Wanner solvers.
Note that this package is not automatically included with DifferentialEquations.jl. To use the following algorithms, you must install and use ODEInterfaceDiffEq.jl:
using Pkg
Pkg.add("ODEInterfaceDiffEq")
import ODEInterfaceDiffEqThese methods can be used independently of the rest of DifferentialEquations.jl.
ODE Solver APIs
ODEInterfaceDiffEq.dopri5 — Type
dopri5()Hairer's classic implementation of the Dormand-Prince 4/5 method.
An explicit Runge-Kutta method of order 5(4) for non-stiff problems. This is one of the most widely used explicit methods for solving ODEs and is generally efficient for problems where evaluations of the right-hand side are relatively cheap.
Solver Properties
- Order: 5
- Adaptive: Yes
- Type: Explicit Runge-Kutta
- Suitable for: Non-stiff problems
ODEInterfaceDiffEq.dop853 — Type
dop853()Explicit Runge-Kutta 8(5,3) method by Dormand-Prince.
A high-order explicit Runge-Kutta method for non-stiff problems requiring high accuracy. Uses adaptive step size control and is particularly effective when tight tolerances are needed.
Solver Properties
- Order: 8
- Adaptive: Yes
- Type: Explicit Runge-Kutta
- Suitable for: Non-stiff problems requiring high accuracy
ODEInterfaceDiffEq.odex — Type
odex()GBS (Gragg-Bulirsch-Stoer) extrapolation algorithm based on the midpoint rule.
An extrapolation method that can achieve very high accuracy for smooth non-stiff problems. Particularly efficient for problems where high accuracy is required and the solution is smooth.
Solver Properties
- Order: up to 12
- Adaptive: Yes
- Type: Extrapolation
- Suitable for: Smooth non-stiff problems requiring high accuracy
ODEInterfaceDiffEq.seulex — Type
seulex(; jac_lower = nothing, jac_upper = nothing)Extrapolation algorithm based on the linearly implicit Euler method.
An implicit extrapolation method suitable for stiff and differential-algebraic problems. Can handle mildly stiff to stiff systems efficiently with automatic stiffness detection.
Arguments
jac_lower: Lower bandwidth of the Jacobian for banded matricesjac_upper: Upper bandwidth of the Jacobian for banded matrices
Solver Properties
- Order: up to 12
- Adaptive: Yes
- Type: Implicit extrapolation
- Suitable for: Stiff problems and DAEs
ODEInterfaceDiffEq.radau — Type
radau(; jac_lower = nothing, jac_upper = nothing)Implicit Runge-Kutta (Radau IIA) method of variable order between 5 and 13.
A high-order implicit Runge-Kutta method particularly well-suited for stiff ODEs and DAEs. The Radau IIA methods have excellent stability properties and can handle very stiff problems efficiently.
Arguments
jac_lower: Lower bandwidth of the Jacobian for banded matricesjac_upper: Upper bandwidth of the Jacobian for banded matrices
Solver Properties
- Order: 5 to 13 (variable)
- Adaptive: Yes
- Type: Implicit Runge-Kutta (Radau IIA)
- Suitable for: Stiff problems and DAEs
ODEInterfaceDiffEq.radau5 — Type
radau5(; jac_lower = nothing, jac_upper = nothing)Implicit Runge-Kutta method (Radau IIA) of order 5.
A robust implicit method for stiff ODEs and DAEs. This is the fixed order 5 version of the Radau method, often preferred when a specific order is desired or for comparison purposes.
Arguments
jac_lower: Lower bandwidth of the Jacobian for banded matricesjac_upper: Upper bandwidth of the Jacobian for banded matrices
Solver Properties
- Order: 5
- Adaptive: Yes
- Type: Implicit Runge-Kutta (Radau IIA)
- Suitable for: Stiff problems and DAEs
ODEInterfaceDiffEq.rodas — Type
rodas(; jac_lower = nothing, jac_upper = nothing)Rosenbrock method of order 4(3).
A Rosenbrock-type method for stiff ODEs and DAEs. Rosenbrock methods are semi-implicit, using only one LU decomposition per step, which can make them more efficient than fully implicit methods for moderately stiff problems.
Arguments
jac_lower: Lower bandwidth of the Jacobian for banded matricesjac_upper: Upper bandwidth of the Jacobian for banded matrices
Solver Properties
- Order: 4
- Adaptive: Yes
- Type: Rosenbrock
- Suitable for: Stiff problems and DAEs
ODEInterfaceDiffEq.ddeabm — Type
ddeabm()Adams-Bashforth-Moulton predictor-corrector method with variable order (1 to 12).
A multistep method based on Adams formulas, using a predictor-corrector approach. Efficient for smooth non-stiff problems, particularly when many function evaluations are expensive, as it reuses previous function evaluations.
Solver Properties
- Order: 1 to 12 (variable)
- Adaptive: Yes
- Type: Adams-Bashforth-Moulton multistep
- Suitable for: Smooth non-stiff problems
ODEInterfaceDiffEq.ddebdf — Type
ddebdf(; jac_lower = nothing, jac_upper = nothing)Backward Differentiation Formula (BDF) with variable order (1 to 5).
A multistep method particularly effective for stiff ODEs and DAEs. BDF methods are the industry standard for stiff problems in many scientific computing applications.
Arguments
jac_lower: Lower bandwidth of the Jacobian for banded matricesjac_upper: Upper bandwidth of the Jacobian for banded matrices
Solver Properties
- Order: 1 to 5 (variable)
- Adaptive: Yes
- Type: BDF multistep
- Suitable for: Stiff problems and DAEs
BVP Solver APIs
Missing docstring for ODEInterfaceDiffEq.BVPM2. Check Documenter's build log for details.