ODE Solvers

This page documents all problem types, ODE function containers, integrator functions, and time-stepping algorithms provided by ClimaTimeSteppers.jl.

For the mathematical formulation of these algorithms, see the Algorithm Formulations section.

Problem Types

These types define what to solve. ODEProblem is the main entry point; SplitODEProblem is a convenience constructor for multirate problems. Types not exported — access them qualified (e.g. CTS.ODEProblem).

ClimaTimeSteppers.ODEProblemType
ODEProblem(f, u0, tspan, p)

An ordinary differential equation problem $du/dt = f(u, p, t)$.

Arguments

  • f: the ODE function (ClimaODEFunction, IncrementingODEFunction, SplitFunction, etc.)
  • u0: initial state
  • tspan: time span (t_start, t_end)
  • p: parameters passed to f
source
ClimaTimeSteppers.SplitFunctionType
SplitFunction(f1, f2)

A split ODE function representing $du/dt = f_1(u,p,t) + f_2(u,p,t)$.

Used with Multirate algorithms where f1 is the fast component and f2 is the slow component.

Arguments

  • f1: fast component (typically an IncrementingODEFunction)
  • f2: slow component
source
ClimaTimeSteppers.IncrementingODEFunctionType
IncrementingODEFunction{true}(f)

An ODE function that supports the incrementing call form f(du, u, p, t, α=true, β=false), which computes du .= α .* f(u,p,t) .+ β .* du.

Required by LowStorageRungeKutta2N methods.

Arguments

  • f: callable with signature f(du, u, p, t, α, β)

The type parameter iip (true for in-place, false for out-of-place) is specified as a curly-brace parameter: IncrementingODEFunction{true}(f).

source
ClimaTimeSteppers.ODEFunctionType
ODEFunction(f; jac_prototype, Wfact, tgrad)

An ODE function wrapper that optionally carries Jacobian information for implicit solvers. Used as the T_imp! field of ClimaODEFunction.

Arguments

  • f: callable with standard ODE signature f(du, u, p, t)

Keyword Arguments

  • jac_prototype: prototype matrix for the Jacobian (determines sparsity/structure)
  • Wfact: function Wfact(W, u, p, dtγ, t) that computes $W = J \Delta t \gamma - I$
  • tgrad: function tgrad(∂f∂t, u, p, t) for the explicit time derivative
source
ClimaTimeSteppers.ODESolutionType
ODESolution{T, U, P, A}

Solution object returned by solve and solve!.

Fields

  • t::Vector{T}: saved time points
  • u::Vector{U}: saved state values (one per time point)
  • prob: the original ODEProblem
  • alg: the algorithm used

Calling sol(t) returns the saved state nearest to time t (not an interpolation — nearest-neighbor lookup only).

source

ODE Function Types

These types define how the right-hand side is structured. Most users will use ClimaODEFunction for IMEX and Rosenbrock methods, or IncrementingODEFunction for low-storage RK methods.

ClimaTimeSteppers.ClimaODEFunctionType
ClimaODEFunction(; T_imp!, [dss!], [cache!], [cache_imp!])
ClimaODEFunction(; T_exp!, [T_lim!], [T_imp!], [lim!], [dss!], [cache!], [cache_imp!])
ClimaODEFunction(; T_exp_T_lim!, [T_imp!], [lim!], [dss!], [cache!], [cache_imp!])

Container for all tendency and auxiliary functions used by IMEX and Rosenbrock time-stepping algorithms. Tendencies set to nothing are skipped, avoiding unnecessary allocations.

Keyword Arguments

Tendency functions (at least one must be provided):

  • T_exp!(du, u, p, t): explicit tendency (not limited)
  • T_lim!(du, u, p, t): explicit tendency passed through the limiter
  • T_exp_T_lim!(du_exp, du_lim, u, p, t): fused alternative to separate T_exp!/T_lim!
  • T_imp!: implicit tendency — typically an ClimaTimeSteppers.ODEFunction carrying Jacobian info
  • T_imp_subproblem!: optional implicit tendency for a preconditioning Newton solve executed before the main T_imp! solve. The subproblem mutates u in place, and its result becomes the initial guess for the main solve — it does not define an additive tendency. Setting T_imp_subproblem! = T_imp! applies the implicit equation twice (compounding), giving different results from the standard single-solve path.

Auxiliary functions (default to no-ops):

  • lim!(u, p, t, u_ref): limiter applied after incrementing u from u_ref by T_lim!
  • dss!(u, p, t): direct stiffness summation (spectral element continuity)
  • cache!(u, p, t): update the parameter cache p to reflect state u
  • cache_imp!(u, p, t): update cache components needed by T_imp! (defaults to cache!)
  • initialize_subproblem!(u, p, γdt): set up the subproblem before the Newton solve

Internally, T_exp! and T_lim! are merged into a single T_exp_T_lim! at construction time.

source
ClimaTimeSteppers.ForwardEulerODEFunctionType
ForwardEulerODEFunction(f; jac_prototype, Wfact, tgrad)

An ODE function whose call signature is f(un, u, p, t, dt), computing a forward-Euler-style update un .= u .+ dt * tendency(u, p, t).

Arguments

  • f: callable with signature f(un, u, p, t, dt)

Keyword Arguments

  • jac_prototype: prototype matrix for the Jacobian
  • Wfact: function Wfact(W, u, p, dtγ, t) computing $W = J \Delta t \gamma - I$
  • tgrad: function tgrad(∂f∂t, u, p, t) for the explicit time derivative
source

Integrator

The integrator is the stateful object that advances the solution forward in time. Create one with init, advance with step!, and run to completion with solve!. The convenience function solve combines init and solve!.

ClimaTimeSteppers.TimeStepperIntegratorType
TimeStepperIntegrator

The main integrator type for ClimaTimeSteppers. Created by init, advanced by step!, and run to completion by solve!.

Key fields (user-facing)

  • u: current state vector (mutated in place)
  • p: parameters
  • t: current time
  • dt: current timestep (may differ from _dt near tstops)
  • sol: the ODESolution being built
source
ClimaTimeSteppers.initFunction
init(prob, alg; dt, kwargs...)

Create a TimeStepperIntegrator for the given problem and algorithm.

Arguments

Keyword Arguments

  • dt: timestep size (required; must be positive)
  • tstops: additional times at which the integrator must stop
  • saveat: times at which to save the solution (default: only endpoints)
  • save_everystep: save after every step (default: false)
  • callback: a DiscreteCallback or CallbackSet
  • advance_to_tstop: if true, step! advances to the next tstop
  • save_func: function (u, t) -> value applied before saving (default: copy)
  • dtchangeable: allow dt reduction near tstops (default: true)
  • stepstop: stop after this many steps (-1 = unlimited)

Returns

source
ClimaTimeSteppers.step!Function
step!(integrator)
step!(integrator, dt, stop_at_tdt=false)

Advance the integrator. The single-argument form takes one internal step (or advances to the next tstop if advance_to_tstop was set). The two-argument form advances by exactly dt time units.

Arguments

  • dt: time interval to advance (must be positive)
  • stop_at_tdt: if true, add t + dt as a tstop to guarantee exact stopping
source
ClimaTimeSteppers.reinit!Function
reinit!(integrator, [u0]; t0, tf, erase_sol, tstops, saveat, reinit_callbacks)

Reset the integrator to a new initial condition without reallocating caches.

Arguments

  • u0: new initial state (default: original prob.u0)

Keyword Arguments

  • t0, tf: new time span endpoints (default: original prob.tspan)
  • erase_sol: clear saved solution data (default true)
  • tstops, saveat: new stop/save schedules (default: reuse originals)
  • reinit_callbacks: re-run callback initialization (default true)
source
ClimaTimeSteppers.add_tstop!Function
add_tstop!(integrator, t)

Schedule a mandatory stop at time t. The integrator will reduce its timestep if necessary to land exactly on t.

source

Abstract Types

ClimaTimeSteppers.AbstractAlgorithmConstraintType
AbstractAlgorithmConstraint

A mechanism for constraining which operations can be performed by an algorithm for solving ODEs.

For example, an unconstrained algorithm might compute a Runge-Kutta stage by taking linear combinations of tendencies; i.e., by adding quantities of the form dt * tendency(state). On the other hand, a "strong stability preserving" algorithm can only take linear combinations of "incremented states"; i.e., it only adds quantities of the form state + dt * coefficient * tendency(state).

source
ClimaTimeSteppers.SSPType
SSP

Indicates that an algorithm must be "strong stability preserving", which makes it easier to guarantee that the algorithm will preserve monotonicity properties satisfied by the initial state. For example, this ensures that the algorithm will be able to use limiters in a mathematically consistent way.

source

IMEX Algorithms

Implicit-explicit additive Runge-Kutta methods. Use with ClimaODEFunction for problems split into stiff (implicit) and non-stiff (explicit) tendencies. See IMEX ARK formulation for the mathematical details.

ClimaTimeSteppers.IMEXTableauType
IMEXTableau(; a_exp, a_imp, [b_exp], [b_imp], [c_exp], [c_imp])

A pair of Butcher tableaux for IMEX (implicit-explicit) Runge-Kutta methods.

The explicit tableau must be strictly lower triangular and the implicit tableau lower triangular (DIRK). Usually constructed indirectly via an algorithm name such as ARS343.

Keyword Arguments

  • a_exp: explicit Butcher matrix (strictly lower triangular, required)
  • a_imp: implicit Butcher matrix (lower triangular / DIRK, required)
  • b_exp: explicit weights (default: last row of a_exp, i.e. FSAL)
  • b_imp: implicit weights (default: last row of a_imp)
  • c_exp: explicit abscissae (default: row sums of a_exp)
  • c_imp: implicit abscissae (default: row sums of a_imp)
source
ClimaTimeSteppers.IMEXAlgorithmType
IMEXAlgorithm(name, newtons_method, [constraint])
IMEXAlgorithm(tableau, newtons_method, [constraint])

An implicit-explicit Runge-Kutta algorithm.

The preferred constructor takes an algorithm name (e.g. ARS343) which determines the tableau and default constraint automatically. A raw IMEXTableau can also be passed directly.

Arguments

  • name: an algorithm name such as ARS343(), SSP333(), etc.
  • newtons_method: a NewtonsMethod for the implicit stages (nothing if the algorithm is purely explicit)
  • constraint: Unconstrained (default) or SSP

Example

alg = IMEXAlgorithm(ARS343(), NewtonsMethod(; max_iters = 2))
source

ARS family ([3])

ClimaTimeSteppers.ARS111Type
ARS111

An IMEX ARK algorithm from [3], section 2, with 1 implicit stage, 1 explicit stage and 1st order accuracy. Also called IMEX Euler or forward-backward Euler.

source
ClimaTimeSteppers.ARS121Type
ARS121

An IMEX ARK algorithm from [3], section 2, with 1 implicit stage, 2 explicit stages, and 1st order accuracy. Also called IMEX Euler or forward-backward Euler.

source
ClimaTimeSteppers.ARS122Type
ARS122

An IMEX ARK algorithm from [3], section 2, with 1 implicit stage, 2 explicit stages, and 2nd order accuracy. Also called IMEX midpoint.

source

IMKG family ([8])

SSP-IMEX family

ClimaTimeSteppers.SSP222Type
SSP222

An IMEX SSPRK algorithm from [9], with 2 implicit stages, 2 explicit stages, and 2nd order accuracy. Also called SSP2(222) in [20].

source
ClimaTimeSteppers.SSP332Type
SSP332

An IMEX SSPRK algorithm from [9], with 3 implicit stages, 3 explicit stages, and 2nd order accuracy. Also called SSP2(332)a in [20].

source
ClimaTimeSteppers.SSP333Type
SSP333(; β = 1/2 + √3/6)

Family of IMEX SSPRK algorithms parametrized by the value β from [10], Section 3.2, with 3 implicit stages, 3 explicit stages, and 3rd order accuracy. The default value of β results in an SDIRK algorithm, which is also called SSP3(333)c in [20].

source
ClimaTimeSteppers.SSP433Type
SSP433

An IMEX SSPRK algorithm from [9], with 4 implicit stages, 3 explicit stages, and 3rd order accuracy. Also called SSP3(433) in [20].

source

Other ARK methods

ClimaTimeSteppers.DBM453Type
DBM453

An IMEX ARK algorithm from [6], Appendix A, with 4 implicit stages, 5 explicit stages, and 3rd order accuracy.

source
ClimaTimeSteppers.ARK2GKCType
ARK2GKC(; paper_version = false)

An IMEX ARK algorithm from [4] with 2 implicit stages, 3 explicit stages, and 2nd order accuracy. If paper_version = true, the algorithm uses coefficients from the paper. Otherwise, it uses coefficients that make it more stable but less accurate.

source
ClimaTimeSteppers.ARK437L2SA1Type
ARK437L2SA1

An IMEX ARK algorithm from [5], Table 8, with 6 implicit stages, 7 explicit stages, and 4th order accuracy. Written as ARK4(3)7L[2]SA₁ in the paper.

source
ClimaTimeSteppers.ARK548L2SA2Type
ARK548L2SA2

An IMEX ARK algorithm from [5], Table 8, with 7 implicit stages, 8 explicit stages, and 5th order accuracy. Written as ARK5(4)8L[2]SA₂ in the paper.

source

Explicit Algorithms

Explicit Runge-Kutta methods for non-stiff problems.

ClimaTimeSteppers.ExplicitTableauType
ExplicitTableau(; a, [b], [c])

An explicit Butcher tableau. Usually constructed indirectly via an algorithm name such as SSP33ShuOsher.

Keyword Arguments

  • a: Butcher matrix (strictly lower triangular, required)
  • b: weights (default: last row of a, i.e., first same as last, FSAL)
  • c: abscissae (default: row sums of a)
source
ClimaTimeSteppers.ExplicitAlgorithmFunction
ExplicitAlgorithm(name, [constraint])
ExplicitAlgorithm(tableau, [constraint])

An explicit Runge-Kutta algorithm. Shorthand for an IMEXAlgorithm with identical explicit/implicit tableaux and no Newton solver.

Arguments

  • name: an algorithm name such as SSP33ShuOsher() or RK4()
  • constraint: Unconstrained (default) or SSP

Example

alg = ExplicitAlgorithm(SSP33ShuOsher())
source

Low-Storage Runge-Kutta

Memory-efficient explicit methods using only two state-sized arrays. Require an IncrementingODEFunction. See the LSRK formulation.

Multirate Methods

For problems with separated fast and slow timescales, using a SplitODEProblem. See the multirate formulation.

ClimaTimeSteppers.MultirateType
Multirate(fast, slow)

A multirate Runge-Kutta scheme that pairs a slow (outer) algorithm with a fast (inner) algorithm. The problem must be a SplitODEProblem where f1 is the fast tendency and f2 is the slow tendency.

Arguments

Pass fast_dt as a keyword argument to init or solve to set the inner timestep.

Example

using ClimaTimeSteppers
import ClimaTimeSteppers as CTS

prob = CTS.SplitODEProblem(f_fast, f_slow, u0, tspan, p)
alg  = Multirate(LSRK54CarpenterKennedy(), MIS3C())
sol  = CTS.solve(prob, alg; dt = 0.1, fast_dt = 0.01)
source

Multirate Infinitesimal Step (MIS)

Wicker-Skamarock

Rosenbrock Methods

Linearly implicit methods that replace Newton iterations with a single linear solve per stage. See the Rosenbrock formulation.

ClimaTimeSteppers.RosenbrockAlgorithmType
RosenbrockAlgorithm(tableau)

A Rosenbrock-type ODE algorithm. The implicit system at each stage is a single linear solve (no Newton iteration), making Rosenbrock methods cheaper per step than fully implicit IMEX methods when updating the Jacobian is inexpensive.

Arguments

Example

import ClimaTimeSteppers as CTS

alg = CTS.RosenbrockAlgorithm(CTS.tableau(SSPKnoth()))
source