OrdinaryDiffEqExtrapolation
Explicit extrapolation methods that achieve high accuracy through Richardson extrapolation of basic integration schemes. These methods provide adaptive order capabilities and natural parallelism, though they are generally outclassed by modern Runge-Kutta methods for most non-stiff problems.
Key Properties
Extrapolation methods provide:
- Adaptive order capability allowing arbitrarily high orders
- Natural parallelism across different substep sequences
- High accuracy potential for very smooth problems
- Richardson extrapolation to eliminate lower-order error terms
- Automatic stepsize and order control
- Theoretical appeal but often practical limitations
When to Use Extrapolation Methods
These methods are recommended for:
- Very smooth problems where high-order accuracy is beneficial
- Extremely low tolerance requirements where adaptive order helps
- Parallel computing environments that can exploit the natural parallelism
- Research applications exploring adaptive order techniques
- Problems where other high-order methods struggle with accuracy
Important Limitations
- Generally outclassed by modern explicit RK methods (Tsit5, Verner methods)
- Higher computational overhead compared to optimized RK methods
- Best suited for very smooth functions - poor performance on non-smooth problems
- Parallel efficiency gains often don't compensate for increased work
Mathematical Background
Extrapolation methods use sequences of basic integrators (like Euler or midpoint) with different stepsizes, then apply Richardson extrapolation to achieve higher-order accuracy. The adaptive order capability comes from using longer extrapolation sequences.
Solver Selection Guide
Explicit extrapolation methods
AitkenNeville
: Euler extrapolation using Aitken-Neville algorithmExtrapolationMidpointDeuflhard
: Midpoint extrapolation with barycentric coordinatesExtrapolationMidpointHairerWanner
: Midpoint extrapolation following ODEX algorithm
When to consider these methods
- Very low tolerances (< 1e-12) where adaptive order might help
- Extremely smooth problems with analytic solutions
- Parallel computing scenarios with many available cores
- Comparison studies with other high-order methods
Better alternatives for most problems
- For high accuracy: Use Verner methods (Vern7, Vern8, Vern9)
- For general problems: Use Tsit5 or appropriate RK method
- For stiff problems: Consider implicit extrapolation methods
Performance Notes
- Consider stiff extrapolation methods which can perform very well for sufficiently stiff problems
- Test against Verner methods before choosing extrapolation for high accuracy
- Parallelism benefits are problem and hardware dependent
- Most effective on very smooth, well-behaved problems
Installation
To be able to access the solvers in OrdinaryDiffEqExtrapolation
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqExtrapolation")
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 OrdinaryDiffEqExtrapolation
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, ExtrapolationMidpointDeuflhard())
Full list of solvers
OrdinaryDiffEqExtrapolation.AitkenNeville
— TypeAitkenNeville(; max_order::Int = 10,
min_order::Int = 1,
init_order = 3,
thread = OrdinaryDiffEq.False())
Parallelized Explicit Extrapolation Method. Euler extrapolation using Aitken-Neville with the Romberg Sequence.
Keyword Arguments
max_order
: maximum order of the adaptive order algorithm.min_order
: minimum order of the adaptive order algorithm.init_order
: initial order of the adaptive order algorithm.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
@inproceedings{elrod2022parallelizing, title={Parallelizing explicit and implicit extrapolation methods for ordinary differential equations}, author={Elrod, Chris and Ma, Yingbo and Althaus, Konstantin and Rackauckas, Christopher and others}, booktitle={2022 IEEE High Performance Extreme Computing Conference (HPEC)}, pages={1–9}, year={2022}, organization={IEEE}}
OrdinaryDiffEqExtrapolation.ExtrapolationMidpointDeuflhard
— TypeExtrapolationMidpointDeuflhard(; max_order = 10,
min_order = 1,
init_order = 5,
thread = OrdinaryDiffEq.True(),
sequence = :harmonic,
sequence_factor = 2)
Parallelized Explicit Extrapolation Method. Midpoint extrapolation using Barycentric coordinates.
Keyword Arguments
max_order
: maximum order of the adaptive order algorithm.min_order
: minimum order of the adaptive order algorithm.init_order
: initial order of the adaptive order algorithm.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.sequence
: the step-number sequences, also called the subdividing sequence. Possible values are:harmonic
,:romberg
or:bulirsch
.sequence_factor
: denotes which even multiple of sequence to take while evaluating internal discretizations.
References
@inproceedings{elrod2022parallelizing, title={Parallelizing explicit and implicit extrapolation methods for ordinary differential equations}, author={Elrod, Chris and Ma, Yingbo and Althaus, Konstantin and Rackauckas, Christopher and others}, booktitle={2022 IEEE High Performance Extreme Computing Conference (HPEC)}, pages={1–9}, year={2022}, organization={IEEE}}
OrdinaryDiffEqExtrapolation.ExtrapolationMidpointHairerWanner
— TypeExtrapolationMidpointHairerWanner(; max_order = 10,
min_order = 2,
init_order = 5,
thread = OrdinaryDiffEq.True(),
sequence = :harmonic,
sequence_factor = 2)
Parallelized Explicit Extrapolation Method. Midpoint extrapolation using Barycentric coordinates, following Hairer's ODEX in the adaptivity behavior.
Keyword Arguments
max_order
: maximum order of the adaptive order algorithm.min_order
: minimum order of the adaptive order algorithm.init_order
: initial order of the adaptive order algorithm.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.sequence
: the step-number sequences, also called the subdividing sequence. Possible values are:harmonic
,:romberg
or:bulirsch
.sequence_factor
: denotes which even multiple of sequence to take while evaluating internal discretizations.
References
@inproceedings{elrod2022parallelizing, title={Parallelizing explicit and implicit extrapolation methods for ordinary differential equations}, author={Elrod, Chris and Ma, Yingbo and Althaus, Konstantin and Rackauckas, Christopher and others}, booktitle={2022 IEEE High Performance Extreme Computing Conference (HPEC)}, pages={1–9}, year={2022}, organization={IEEE}}