OrdinaryDiffEqExponentialRK
Exponential Runge-Kutta methods are specialized integrators for semi-linear differential equations of the form du/dt = Au + f(u,t)
, where A
is a linear operator (often representing diffusion or dispersion) and f
represents nonlinear terms. These methods are particularly effective for stiff linear parts combined with non-stiff nonlinear terms. Important: The nonlinear term f(u,t)
must be non-stiff for these methods to be effective.
Key Properties
Exponential RK methods provide:
- Exact treatment of linear parts using matrix exponential functions
- High-order accuracy for both linear and nonlinear components
- Excellent stability properties for problems with stiff linear operators
- Efficient handling of semi-linear PDEs after spatial discretization
- Reduced timestep restrictions compared to traditional explicit methods
- Preservation of qualitative behavior for many physical systems
When to Use Exponential RK Methods
These methods are recommended for:
- Semi-linear PDEs with stiff diffusion/dispersion and moderate non-stiff nonlinearity
- Reaction-diffusion systems with fast diffusion and slower non-stiff reactions
- Nonlinear Schrödinger equations and other dispersive wave equations with non-stiff nonlinear terms
- Pattern formation problems (Turing patterns, phase field models) where nonlinearity is non-stiff
- Quantum dynamics with linear Hamiltonian and non-stiff nonlinear interactions
- Problems with strong linear damping or oscillatory linear parts combined with non-stiff nonlinear terms
- Spatially discretized PDEs where the linear part dominates stiffness but the nonlinear part remains non-stiff
Mathematical Background
For problems du/dt = Au + f(u,t)
, exponential methods compute the exact solution of the linear part Au
using exp(A*dt)
and treat the nonlinear part f(u,t)
with Runge-Kutta-like stages. This approach is particularly effective when A
represents well-understood physics (diffusion, dispersion, linear oscillations).
Solver Selection Guide
Basic exponential time differencing (ETD)
LawsonEuler
: First-order exponential Euler methodNorsettEuler
/ETD1
: Alternative first-order schemeETDRK2
: Second-order exponential RKETDRK3
: Third-order exponential RKETDRK4
: Fourth-order exponential RK, popular choiceETD2
: Second-order exponential time differencing (in development)
High-order specialized methods
HochOst4
: Fourth-order exponential RK with enhanced stabilityExp4
: Fourth-order EPIRK scheme
Adaptive exponential Rosenbrock
Exprb32
: Third-order adaptive method with error controlExprb43
: Fourth-order adaptive method
EPIRK (Exponential Propagation Iterative RK) methods
EPIRK4s3A
: Fourth-order with stiff order 4EPIRK4s3B
: Alternative fourth-order variantEPIRK5s3
: Fifth-order method (note: marked as broken)EXPRB53s3
: Fifth-order with stiff order 5EPIRK5P1
,EPIRK5P2
: Fifth-order variants
Performance Recommendations
- For most semi-linear problems:
ETDRK4
- For adaptive stepsize:
Exprb43
- For high stiffness in linear part:
EPIRK4s3A
orEPIRK4s3B
- For maximum accuracy:
EXPRB53s3
Implementation Requirements
These methods require:
- Computation of matrix exponentials
exp(A*dt)
and related functions - Krylov subspace methods for large systems (automatic in most cases)
- Proper problem formulation with identified linear and nonlinear parts
Installation
To be able to access the solvers in OrdinaryDiffEqLinear
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqExponentialRK")
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
Installation
To be able to access the solvers in OrdinaryDiffEqExponentialRK
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqExponentialRK")
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 OrdinaryDiffEqExponentialRK
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, EPIRK5s3())
Full list of solvers
OrdinaryDiffEqExponentialRK.LawsonEuler
— TypeLawsonEuler(; krylov = false,
m = 30,
iop = 0)
Semilinear ODE solver First order exponential Euler scheme (fixed timestepping)
Keyword Arguments
krylov
: Determines whether Krylov approximation or operator caching is used, the latter only available for semilinear problems.krylov=true
is much faster for larger systems and is thus recommended whenever there are >100 ODEs.m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, Marlis, and Alexander Ostermann. “Exponential Integrators.” Acta Numerica 19 (2010): 209–286. doi:10.1017/S0962492910000048.
OrdinaryDiffEqExponentialRK.NorsettEuler
— TypeNorsettEuler(; krylov = false,
m = 30,
iop = 0)
Semilinear ODE solver First order exponential-RK scheme. Alias: ETD1
Keyword Arguments
krylov
: Determines whether Krylov approximation or operator caching is used, the latter only available for semilinear problems.krylov=true
is much faster for larger systems and is thus recommended whenever there are >100 ODEs.m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, Marlis, and Alexander Ostermann. “Exponential Integrators.” Acta Numerica 19 (2010): 209–286. doi:10.1017/S0962492910000048.
OrdinaryDiffEqExponentialRK.ETD2
— TypeETD2: Exponential Runge-Kutta Method Second order Exponential Time Differencing method (in development).
OrdinaryDiffEqExponentialRK.ETDRK2
— TypeETDRK2(; krylov = false,
m = 30,
iop = 0)
Semilinear ODE solver 2nd order exponential-RK scheme.
Keyword Arguments
krylov
: Determines whether Krylov approximation or operator caching is used, the latter only available for semilinear problems.krylov=true
is much faster for larger systems and is thus recommended whenever there are >100 ODEs.m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, Marlis, and Alexander Ostermann. “Exponential Integrators.” Acta Numerica 19 (2010): 209–286. doi:10.1017/S0962492910000048.
OrdinaryDiffEqExponentialRK.ETDRK3
— TypeETDRK3(; krylov = false,
m = 30,
iop = 0)
Semilinear ODE solver 3rd order exponential-RK scheme.
Keyword Arguments
krylov
: Determines whether Krylov approximation or operator caching is used, the latter only available for semilinear problems.krylov=true
is much faster for larger systems and is thus recommended whenever there are >100 ODEs.m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, Marlis, and Alexander Ostermann. “Exponential Integrators.” Acta Numerica 19 (2010): 209–286. doi:10.1017/S0962492910000048.
OrdinaryDiffEqExponentialRK.ETDRK4
— TypeETDRK4(; krylov = false,
m = 30,
iop = 0)
Semilinear ODE solver 4th order exponential-RK scheme (fixed timestepping)
Keyword Arguments
krylov
: Determines whether Krylov approximation or operator caching is used, the latter only available for semilinear problems.krylov=true
is much faster for larger systems and is thus recommended whenever there are >100 ODEs.m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, Marlis, and Alexander Ostermann. “Exponential Integrators.” Acta Numerica 19 (2010): 209–286. doi:10.1017/S0962492910000048.
OrdinaryDiffEqExponentialRK.HochOst4
— TypeHochOst4(; krylov = false,
m = 30,
iop = 0)
Semilinear ODE solver 4th order exponential-RK scheme with stiff order 4.
Keyword Arguments
krylov
: Determines whether Krylov approximation or operator caching is used, the latter only available for semilinear problems.krylov=true
is much faster for larger systems and is thus recommended whenever there are >100 ODEs.m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, Marlis, and Alexander Ostermann. “Exponential Integrators.” Acta Numerica 19 (2010): 209–286. doi:10.1017/S0962492910000048.
Adaptive Exponential Rosenbrock Methods
OrdinaryDiffEqExponentialRK.Exprb32
— TypeExprb32(; m = 30,
iop = 0)
Semilinear ODE solver 3rd order adaptive Exponential-Rosenbrock scheme.
Keyword Arguments
m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, M., & Ostermann, A. (2010). Exponential integrators. Acta Numerica, 19, 209-286. (https://doi.org/10.1017/S0962492910000048)
OrdinaryDiffEqExponentialRK.Exprb43
— TypeExprb43(; m = 30,
iop = 0)
Semilinear ODE solver 4th order adaptive Exponential-Rosenbrock scheme.
Keyword Arguments
m
: Controls the size of Krylov subspace.iop
: If not zero, determines the length of the incomplete orthogonalization procedure (IOP). Note that if the linear operator/Jacobian is hermitian, then the Lanczos algorithm will always be used and the IOP setting is ignored.
References
Hochbruck, M., & Ostermann, A. (2010). Exponential integrators. Acta Numerica, 19, 209-286. (https://doi.org/10.1017/S0962492910000048)