OrdinaryDiffEqHighOrderRK
High-order explicit Runge-Kutta methods for non-stiff differential equations requiring high accuracy. These methods provide alternatives to the Verner methods, though OrdinaryDiffEqVerner
methods generally perform better at low tolerances and should be preferred in most cases.
Key Properties
High-order RK methods provide:
- High-order accuracy (7th and 8th order) for precise integration
- Specialized coefficients for specific problem types
- Dense output capabilities for some methods
- Alternative approaches to the more commonly used Verner methods
When to Use High-Order RK Methods
These methods are recommended when:
- Verner methods are not suitable for specific problem characteristics
- Specialized properties are needed (e.g., phase-fitted methods for oscillatory problems)
- Research or comparison purposes require different high-order method families
- Specific applications benefit from particular method properties
Solver Selection Guide
General high-order integration
Specialized cases where these methods may be preferred
TanYam7
: Seventh-order Tanaka-Yamashita methodTsitPap8
: Eighth-order Tsitouras-Papakostas methodDP8
: Eighth-order Dormand-Prince method (Hairer's 8/5/3 implementation)PFRK87
: Phase-fitted eighth-order method for oscillatory problems
Performance Notes
- Verner methods are generally more efficient for most high-accuracy applications
- These methods are provided for specialized use cases and research purposes
- Consider problem-specific properties when choosing between different high-order families
Installation
To be able to access the solvers in OrdinaryDiffEqHighOrderRK
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqHighOrderRK")
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 OrdinaryDiffEqHighOrderRK
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, DP8())
Full list of solvers
OrdinaryDiffEqHighOrderRK.TanYam7
— TypeTanYam7(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Tanaka-Yamashita 7 Runge-Kutta method.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(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
Tanaka M., Muramatsu S., Yamashita S., (1992), On the Optimization of Some Nine-Stage Seventh-order Runge-Kutta Method, Information Processing Society of Japan, 33 (12), pp. 1512-1526.
OrdinaryDiffEqHighOrderRK.TsitPap8
— TypeTsitPap8(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Tsitouras-Papakostas 8/7 Runge-Kutta method.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(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{tsitouras1999cheap, title={Cheap error estimation for Runge–Kutta methods}, author={Tsitouras, Ch and Papakostas, SN}, journal={SIAM Journal on Scientific Computing}, volume={20}, number={6}, pages={2067–2088}, year={1999}, publisher={SIAM}}
OrdinaryDiffEqHighOrderRK.DP8
— TypeDP8(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Hairer's 8/5/3 adaption of the Dormand-Prince Runge-Kutta method.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(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.
OrdinaryDiffEqHighOrderRK.PFRK87
— TypePFRK87(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
omega = 0.0)
Explicit Runge-Kutta Method. Phase-fitted Runge-Kutta of 8th order.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(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{tsitouras2017phase, title={Phase-fitted Runge–Kutta pairs of orders 8 (7)}, author={Tsitouras, Ch and Famelis, I Th and Simos, TE}, journal={Journal of Computational and Applied Mathematics}, volume={321}, pages={226–231}, year={2017}, publisher={Elsevier}}