ODE solvers¤
See also How to choose a solver.
Term structure
The type of solver chosen determines how the terms
argument of diffeqsolve
should be laid out. Most of them demand that it should be a single AbstractTerm
. But for example diffrax.SemiImplicitEuler
demands that it be a 2-tuple (AbstractTerm, AbstractTerm)
, to represent the two vector fields that solver uses.
If it is different from this default, then you can find the appropriate structure documented below, and available programmatically under <solver>.term_structure
.
Explicit Runge--Kutta (ERK) methods¤
These methods are suitable for most problems.
diffrax.Euler (AbstractItoSolver)
¤
Euler's method.
1st order explicit Runge--Kutta method. Does not support adaptive step sizing. Uses 1 stage. Uses 1st order local linear interpolation for dense/ts output.
When used to solve SDEs, converges to the Itô solution.
diffrax.Heun (AbstractERK, AbstractStratonovichSolver)
¤
Heun's method.
2nd order explicit Runge--Kutta method. Has an embedded Euler method for adaptive step sizing. Uses 2 stages. Uses 2nd-order Hermite interpolation for dense/ts output.
Also sometimes known as either the "improved Euler method", "modified Euler method" or "explicit trapezoidal rule".
Should not be confused with Heun's third order method, which is a different (higher order) method occasionally also just referred to as "Heun's method".
When used to solve SDEs, converges to the Stratonovich solution.
diffrax.Midpoint (AbstractERK, AbstractStratonovichSolver)
¤
Midpoint method.
2nd order explicit Runge--Kutta method. Has an embedded Euler method for adaptive step sizing. Uses 2 stages. Uses 2nd order Hermite interpolation for dense/ts output.
Also sometimes known as the "modified Euler method".
When used to solve SDEs, converges to the Stratonovich solution.
diffrax.Ralston (AbstractERK, AbstractStratonovichSolver)
¤
Ralston's method.
2nd order explicit Runge--Kutta method. Has an embedded Euler method for adaptive step sizing. Uses 2 stages. Uses 2nd order Hermite interpolation for dense output.
When used to solve SDEs, converges to the Stratonovich solution.
diffrax.Bosh3 (AbstractERK)
¤
Bogacki--Shampine's 3/2 method.
3rd order explicit Runge--Kutta method. Has an embedded 2nd order method for adaptive step sizing. Uses 4 stages with FSAL. Uses 3rd order Hermite interpolation for dense/ts output.
Also sometimes known as "Ralston's third order method".
diffrax.Tsit5 (AbstractERK)
¤
Tsitouras' 5/4 method.
5th order explicit Runge--Kutta method. Has an embedded 4th order method for adaptive step sizing. Uses 7 stages with FSAL. Uses 5th order interpolation for dense/ts output.
Reference
@article{tsitouras2011runge,
title={Runge--Kutta pairs of order 5 (4) satisfying only the first column
simplifying assumption},
author={Tsitouras, Ch},
journal={Computers \& Mathematics with Applications},
volume={62},
number={2},
pages={770--775},
year={2011},
publisher={Elsevier}
}
diffrax.Dopri5 (AbstractERK)
¤
Dormand-Prince's 5/4 method.
5th order Runge--Kutta method. Has an embedded 4th order method for adaptive step sizing. Uses 7 stages with FSAL. Uses 5th order interpolation for dense/ts output.
Reference
The original reference for Dormand--Prince's 5(4) method is:
@article{dormand1980family,
author={Dormand, J. R. and Prince, P. J.},
title={A family of embedded {R}unge--{K}utta formulae},
journal={J. Comp. Appl. Math},
year={1980},
volume={6},
pages={19--26}
}
However (despite the name), the Butcher tableau used here is actually due to Shampine:
@article{shampine1986some,
author={Lawrence F. Shampine},
journal={Mathematics of Computation},
number={173},
pages={135--150},
publisher={American Mathematical Society},
title={Some Practical {R}unge-{K}utta Formulas},
volume={46},
year={1986},
doi={https://doi.org/10.2307/2008219}
}
diffrax.Dopri8 (AbstractERK)
¤
Dormand--Prince's 8/7 method.
8th order Runge--Kutta method. Has an embedded 7th order method for adaptive step sizing. Uses 14 stages with FSAL. Uses 8th order interpolation for dense/ts output.
References
Coefficients from:
@article{prince1981high,
author={Prince, P. J and Dormand, J. R.},
title={High order embedded {R}unge--{K}utta formulae},
journal={J. Comp. Appl. Math},
year={1981},
volume={7},
number={1},
pages={67--75}
}
Dense interpolation from:
@article{bogacki1990interpolating,
author={Bogacki, P. and Shampine, L. F.},
title={Interpolating high-order {R}unge--{K}utta formulas},
journal={Comput. Math. with Appl.},
year={1990},
volume={20},
number={3},
pages={15--24},
doi={https://doi.org/10.1016/0898-1221(90)90027-H}
}
Implicit Runge--Kutta (IRK) methods¤
These methods are suitable for stiff problems.
Each of these takes a root_finder
argument at initialisation, defaulting to a Newton solver, which is used to solve the implicit problem at each step. See the page on root finders.
diffrax.ImplicitEuler (AbstractImplicitSolver, AbstractAdaptiveSolver)
¤
Implicit Euler method.
A-B-L stable 1st order SDIRK method. Has an embedded 2nd order Heun method for adaptive step sizing. Uses 1 stage. Uses a 1st order local linear interpolation for dense/ts output.
diffrax.Kvaerno3 (AbstractESDIRK)
¤
Kvaerno's 3/2 method.
A-L stable stiffly accurate 3rd order ESDIRK method. Has an embedded 2nd order method for adaptive step sizing. Uses 4 stages with FSAL. Uses 3rd order Hermite interpolation for dense/ts output.
Reference
@article{kvaerno2004singly,
title={Singly diagonally implicit Runge--Kutta methods with an explicit first
stage},
author={Kv{\ae}rn{\o}, Anne},
journal={BIT Numerical Mathematics},
volume={44},
number={3},
pages={489--502},
year={2004},
publisher={Springer}
}
diffrax.Kvaerno4 (AbstractESDIRK)
¤
Kvaerno's 4/3 method.
A-L stable stiffly accurate 4th order ESDIRK method. Has an embedded 3rd order method for adaptive step sizing. Uses 5 stages with FSAL. Uses 3rd order Hermite interpolation for dense/ts output.
When solving an ODE over the interval \([t_0, t_1]\), note that this method will make some evaluations slightly past \(t_1\).
Reference
@article{kvaerno2004singly,
title={Singly diagonally implicit Runge--Kutta methods with an explicit first
stage},
author={Kv{\ae}rn{\o}, Anne},
journal={BIT Numerical Mathematics},
volume={44},
number={3},
pages={489--502},
year={2004},
publisher={Springer}
}
diffrax.Kvaerno5 (AbstractESDIRK)
¤
Kvaerno's 5/4 method.
A-L stable stiffly accurate 5th order ESDIRK method. Has an embedded 4th order method for adaptive step sizing. Uses 7 stages with FSAL. Uses 3rd order Hermite interpolation for dense/ts output.
When solving an ODE over the interval \([t_0, t_1]\), note that this method will make some evaluations slightly past \(t_1\).
Reference
@article{kvaerno2004singly,
title={Singly diagonally implicit Runge--Kutta methods with an explicit first
stage},
author={Kv{\ae}rn{\o}, Anne},
journal={BIT Numerical Mathematics},
volume={44},
number={3},
pages={489--502},
year={2004},
publisher={Springer}
}
IMEX methods¤
These "implicit-explicit" methods are suitable for problems of the form \(\frac{\mathrm{d}y}{\mathrm{d}t} = f(t, y(t)) + g(t, y(t))\), where \(f\) is the non-stiff part (explicit integration) and \(g\) is the stiff part (implicit integration).
Term structure
These methods should be called with terms=MultiTerm(explicit_term, implicit_term)
.
diffrax.Sil3 (AbstractRungeKutta, AbstractImplicitSolver)
¤
Whitaker--Kar's fast-slow IMEX method.
3rd order in the explicit (ERK) term; 2nd order in the implicit (EDIRK) term. Uses a 2nd-order embedded Heun method for adaptive step sizing. Uses 4 stages with FSAL. Uses 2nd order Hermite interpolation for dense/ts output.
This should be called with terms=MultiTerm(explicit_term, implicit_term)
.
Reference
@article{whitaker2013implicit,
author={Jeffrey S. Whitaker and Sajal K. Kar},
title={Implicit–Explicit Runge–Kutta Methods for Fast–Slow Wave Problems},
journal={Monthly Weather Review},
year={2013},
publisher={American Meteorological Society},
volume={141},
number={10},
doi={https://doi.org/10.1175/MWR-D-13-00132.1},
pages={3426--3434},
}
diffrax.KenCarp3 (AbstractRungeKutta, AbstractImplicitSolver)
¤
Kennedy--Carpenter's 3/2 IMEX method.
3rd order ERK-ESDIRK implicit-explicit (IMEX) method. The implicit part is stiffly accurate and A-L stable. Has an embedded 2nd order method for adaptive step sizing. Uses 4 stages. Uses 2nd order interpolation for dense/ts output.
This should be called with terms=MultiTerm(explicit_term, implicit_term)
.
Reference
@article{kennedy2003additive,
title={Additive Runge--Kutta schemes for convection-diffusion-reaction
equations},
author={Kennedy, Christopher A and Carpenter, Mark H},
journal={Applied numerical mathematics},
volume={44},
number={1-2},
pages={139--181},
year={2003},
publisher={Elsevier}
}
diffrax.KenCarp4 (AbstractRungeKutta, AbstractImplicitSolver)
¤
Kennedy--Carpenter's 4/3 IMEX method.
4th order ERK-ESDIRK implicit-explicit (IMEX) method. The implicit part is stiffly accurate and A-L stable. Has an embedded 3rd order method for adaptive step sizing. Uses 6 stages. Uses 3rd order interpolation for dense/ts output.
This should be called with terms=MultiTerm(explicit_term, implicit_term)
.
Reference
@article{kennedy2003additive,
title={Additive Runge--Kutta schemes for convection-diffusion-reaction
equations},
author={Kennedy, Christopher A and Carpenter, Mark H},
journal={Applied numerical mathematics},
volume={44},
number={1-2},
pages={139--181},
year={2003},
publisher={Elsevier}
}
diffrax.KenCarp5 (AbstractRungeKutta, AbstractImplicitSolver)
¤
Kennedy--Carpenter's 5/4 IMEX method.
5th order ERK-ESDIRK implicit-explicit (IMEX) method. The implicit part is stiffly accurate and A-L stable. Has an embedded 4th order method for adaptive step sizing. Uses 8 stages. Uses 3rd order interpolation for dense/ts output.
This should be called with terms=MultiTerm(explicit_term, implicit_term)
.
Reference
@article{kennedy2003additive,
title={Additive Runge--Kutta schemes for convection-diffusion-reaction
equations},
author={Kennedy, Christopher A and Carpenter, Mark H},
journal={Applied numerical mathematics},
volume={44},
number={1-2},
pages={139--181},
year={2003},
publisher={Elsevier}
}
Symplectic methods¤
These methods are suitable for problems with symplectic structure; that is to say those ODEs of the form
\(\frac{\mathrm{d}v}{\mathrm{d}t}(t) = f(t, w(t))\)
\(\frac{\mathrm{d}w}{\mathrm{d}t}(t) = g(t, v(t))\)
In particular this includes Hamiltonian systems.
Term and state structure
The state of the system (the initial value of which is given by y0
to diffrax.diffeqsolve
) must be a 2-tuple (of PyTrees). The terms (given by the value of terms
to diffrax.diffeqsolve
) must be a 2-tuple of AbstractTerms
.
Letting v, w = y0
and f, g = terms
, then v
is updated according to f(t, w, args)
and w
is updated according to g(t, v, args)
.
See also this Wikipedia page.
diffrax.SemiImplicitEuler (AbstractSolver)
¤
Semi-implicit Euler's method.
Symplectic method. Does not support adaptive step sizing. Uses 1st order local linear interpolation for dense/ts output.
Reversible methods¤
These methods can be run "in reverse": solving from an initial condition y0
to obtain some terminal value y1
, it is possible to reconstruct y0
from y1
with zero truncation error. (There will still be a small amount of floating point error.) This can be done via SaveAt(solver_state=True)
to save the final solver state, and then passing it as diffeqsolve(..., solver_state=solver_state)
on the backwards-in-time pass.
In addition all symplectic methods are reversible, as are some linear multistep methods. (Below are the non-symplectic reversible solvers.)
diffrax.ReversibleHeun (AbstractAdaptiveSolver, AbstractStratonovichSolver)
¤
Reversible Heun method.
Algebraically reversible 2nd order method. Has an embedded 1st order method for adaptive step sizing. Uses 1st order local linear interpolation for dense/ts output.
When used to solve SDEs, converges to the Stratonovich solution.
Reference
@article{kidger2021efficient,
author={Kidger, Patrick and Foster, James and Li, Xuechen and Lyons, Terry},
title={{E}fficient and {A}ccurate {G}radients for {N}eural {SDE}s},
year={2021},
journal={Advances in Neural Information Processing Systems}
}
Linear multistep methods¤
diffrax.LeapfrogMidpoint (AbstractSolver)
¤
Leapfrog/midpoint method.
2nd order linear multistep method. Uses 1st order local linear interpolation for dense/ts output.
Note that this is referred to as the "leapfrog/midpoint method" as this is the name
used by Shampine in the reference below. It should not be confused with any of the
many other "leapfrog methods" (there are several), or with the "midpoint method"
(which is usually taken to refer to the explicit Runge--Kutta method
diffrax.Midpoint
).
Reference
@article{shampine2009stability,
title={Stability of the leapfrog/midpoint method},
author={L. F. Shampine},
journal={Applied Mathematics and Computation},
volume={208},
number={1},
pages={293-298},
year={2009},
}