OrdinaryDiffEqTsit5
The Tsitouras 5/4 method is the recommended default solver for most non-stiff differential equation problems. This method provides an excellent balance of efficiency, reliability, and accuracy.
Key Properties
Tsit5 offers:
- Fifth-order accuracy with embedded fourth-order error estimation
- Excellent efficiency at default tolerances (1e-6 to 1e-3)
- FSAL (First Same As Last) property for computational efficiency
- High-quality interpolation for dense output
- Robust performance across a wide range of problem types
- Optimized coefficients for minimal error in practical applications
When to Use Tsit5
Tsit5 is recommended for:
- Most non-stiff problems as the first choice solver
- Default and higher tolerances (1e-3 to 1e-6)
- General-purpose integration when problem characteristics are unknown
- Educational and research applications as a reliable baseline
- Real-time applications requiring predictable performance
- Problems where simplicity and reliability are preferred over maximum efficiency
Solver Selection Guide
Primary recommendation
Tsit5
: Use as the default choice for non-stiff problems at standard tolerances
Automatic switching
AutoTsit5
: Automatically switches to a stiff solver when stiffness is detected, making it robust for problems of unknown character
When to Consider Alternatives
Consider other solvers when:
- Higher accuracy needed: Use Verner methods (Vern6, Vern7, Vern8, Vern9) for tolerances below 1e-6
- Higher tolerances: Use BS3 or OwrenZen3 for tolerances above 1e-3
- Robust error control needed: Use BS5 when Tsit5 struggles with error estimation
- Equation is stiff: Use implicit methods (SDIRK, BDF) or semi-implicit methods (Rosenbrock) for stiff problems
- Special properties required: Use specialized methods (SSP, symplectic, etc.) for specific problem types
Installation
To be able to access the solvers in OrdinaryDiffEqTsit5
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqTsit5")
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 OrdinaryDiffEqTsit5
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, Tsit5())
Full list of solvers
OrdinaryDiffEqTsit5.Tsit5
— TypeTsit5(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Tsitouras 5/4 Runge-Kutta method. (free 4th order interpolant). Recommended for most non-stiff problems. Good default choice for unknown stiffness. Highly efficient and generic. Very good performance for most non-stiff ODEs. Recommended as default method for unknown stiffness problems.
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{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} }
OrdinaryDiffEqTsit5.AutoTsit5
— FunctionAutomatic switching algorithm that can switch between the (non-stiff) Tsit5()
and stiff_alg
.
AutoTsit5(stiff_alg; kwargs...)
This method is equivalent to AutoAlgSwitch(Tsit5(), stiff_alg; kwargs...)
. To gain access to stiff algorithms you might have to install additional libraries, such as OrdinaryDiffEqRosenbrock
.