Controller API

Starting with OrdinaryDiffEq v7, the adaptive step-size controller is a first-class object. Each algorithm picks one through default_controller, and users override it by passing controller = ... to solve / init. Tuning knobs like qmin, qmax, gamma, qsteady_min/max, beta1, beta2, and failfactor live on the controller — not on integrator.opts — and the per-step state (q11, qold, dtacc, erracc, …) lives on the controller's cache, which is owned by integrator.controller_cache.

This page documents:

  1. The built-in controllers users can pick from.
  2. The CommonControllerOptions bundle of standard step-size knobs that every controller composes.
  3. The accessor functions the integrator uses to read knobs (get_qmin, get_qmax, …).
  4. The minimal contract you implement when adding a new controller.

Picking a controller

solve(prob, alg; controller = MyController()) overrides the default chosen by default_controller(QT, alg). The built-in choices fall into two families.

Generic controllers

These contain the full step-size logic themselves and work with any adaptive algorithm.

OrdinaryDiffEqCore.IControllerType
IController()

The standard (integral) controller is the most basic step size controller. This controller is usually the first one introduced in numerical analysis classes but should only be used rarely in practice because of efficiency problems for many problems/algorithms.

Construct an integral (I) step size controller adapting the time step based on the formula

Δtₙ₊₁ = εₙ₊₁^(1/k) * Δtₙ

where k = get_current_adaptive_order(alg, integrator.cache) + 1 and εᵢ is the inverse of the error estimate get_EEst(integrator) scaled by the tolerance (Hairer, Nørsett, Wanner, 2008, Section II.4). The step size factor is multiplied by the safety factor gamma and clipped to the interval [qmin, qmax]. A step will be accepted whenever the estimated error get_EEst(integrator) is less than or equal to unity. Otherwise, the step is rejected and re-tried with the predicted step size.

References

source
OrdinaryDiffEqCore.PIControllerType
PIController(beta1, beta2)

The proportional-integral (PI) controller is a widespread step size controller with improved stability properties compared to the IController. This controller is the default for most algorithms in OrdinaryDiffEq.jl.

Construct a PI step size controller adapting the time step based on the formula

Δtₙ₊₁ = εₙ₊₁^β₁ * εₙ^β₂ * Δtₙ

where εᵢ are inverses of the error estimates scaled by the tolerance (Hairer, Nørsett, Wanner, 2010, Section IV.2). The step size factor is multiplied by the safety factor gamma and clipped to the interval [qmin, qmax]. A step will be accepted whenever the estimated error get_EEst(integrator) is less than or equal to unity. Otherwise, the step is rejected and re-tried with the predicted step size.

Note

The coefficients beta1, beta2 are not scaled by the order of the method, in contrast to the PIDController. For the PIController, this scaling by the order must be done when the controller is constructed.

References

source
OrdinaryDiffEqCore.PIDControllerType
PIDController(beta1, beta2, beta3=zero(beta1);
              limiter=default_dt_factor_limiter,
              accept_safety=0.81)

The proportional-integral-derivative (PID) controller is a generalization of the PIController and can have improved stability and efficiency properties.

Construct a PID step size controller adapting the time step based on the formula

Δtₙ₊₁ = εₙ₊₁^(β₁/k) * εₙ^(β₂/k) * εₙ₋₁^(β₃/ k) * Δtₙ

where k = min(alg_order, alg_adaptive_order) + 1 and εᵢ are inverses of the error estimates scaled by the tolerance (Söderlind, 2003). The step size factor is limited by the limiter with default value

limiter(x) = one(x) + atan(x - one(x))

as proposed by Söderlind and Wang (2006). A step will be accepted whenever the predicted step size change is bigger than accept_safety. Otherwise, the step is rejected and re-tried with the predicted step size.

Some standard controller parameters suggested in the literature are

Controllerbeta1beta2beta3
basic1.000.000
PI420.60-0.200
PI332//3-1//30
PI340.70-0.400
H211PI1//61//60
H312PID1//181//91//18
Note

In contrast to the PIController, the coefficients beta1, beta2, beta3 are scaled by the order of the method. Thus, standard controllers such as PI42 can use the same coefficients beta1, beta2, beta3 for different algorithms.

Note

In contrast to other controllers, the PIDController does not use the keyword arguments qmin, qmax to limit the step size change or the safety factor gamma. These common keyword arguments are replaced by the limiter and accept_safety to guarantee a smooth behavior (Söderlind and Wang, 2006). Because of this, a PIDController behaves different from a PIController, even if beta1, beta2 are adapted accordingly and iszero(beta3).

References

  • Söderlind (2003) Digital Filters in Adaptive Time-Stepping DOI: 10.1145/641876.641877
  • Söderlind, Wang (2006) Adaptive time-stepping and computational stability DOI: 10.1016/j.cam.2005.03.008
  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836 # limiter of the dt factor (before clipping)
source
OrdinaryDiffEqCore.PredictiveControllerType
PredictiveController()

The Gustafsson acceleration algorithm accelerates changes so that way algorithms can more swiftly change to handle quick transients. This algorithm is thus well-suited for stiff solvers where this can be expected, and is the default for algorithms like the (E)SDIRK methods.

(; qmin, qmax, gamma) = controller
qmax = get_current_qmax(integrator, qmax)
niters = integrator.cache.nlsolver.iter
fac = min(gamma,
    (1 + 2 * integrator.cache.nlsolver.maxiters) * gamma /
    (niters + 2 * integrator.cache.nlsolver.maxiters))
expo = 1 / (get_current_adaptive_order(alg, integrator.cache) + 1)
qtmp = fastpower(get_EEst(integrator), expo) / fac
@fastmath q = max(inv(qmax), min(inv(qmin), qtmp))
cache.qold = q
q

where controller and cache are the controller and its cache stored in integrator.controller_cache. niters is the number of Newton iterations which was required in the most recent step of the algorithm. Note that these values are used differently depending on acceptance and rejectance. When the step is accepted, the following logic is applied:

(; dtacc, erracc) = cache
(; qmin, qmax, gamma, qsteady_min, qsteady_max) = controller
qmax = get_current_qmax(integrator, qmax)
if integrator.success_iter > 0
    expo = 1 / (get_current_adaptive_order(alg, integrator.cache) + 1)
    qgus = (dtacc / integrator.dt) * fastpower((get_EEst(integrator)^2) / erracc, expo)
    qgus = max(inv(qmax), min(inv(qmin), qgus / gamma))
    qacc = max(q, qgus)
else
    qacc = q
end
if qsteady_min <= qacc <= qsteady_max
    qacc = one(qacc)
end
cache.dtacc = integrator.dt
cache.erracc = max(1e-2, get_EEst(integrator))
integrator.dt / qacc

When it rejects, it's the same as the IController:

if integrator.success_iter == 0
    integrator.dt *= 0.1
else
    integrator.dt = integrator.dt / cache.qold
end
source

Algorithm-integrated controllers

For algorithms whose step-size logic is tied to internal state (the BDF order, the Nordsieck history, extrapolation tableau choice, …), the controller is a thin wrapper that just owns the knobs; the real update lives in the algorithm-specific stepsize_controller! / step_accept_controller! / step_reject_controller! methods.

OrdinaryDiffEqBDF.BDFControllerType
BDFController(; qmin, qmax, qsteady_min, qsteady_max, gamma, qmax_first_step,
              failfactor)

Step-size controller for the variable-order BDF family (QNDF, FBDF, DFBDF). Composes the standard step-size knobs via CommonControllerOptions; the adaptive logic is integrated into the algorithm itself, so the cache falls through to alg-level dispatch (the same way the legacy DummyControllerCache did) but exposes the knobs as a real, settable controller. Pass it explicitly with solve(prob, alg; controller = BDFController(...)), or rely on the default constructed by default_controller(QT, alg).

source
OrdinaryDiffEqNordsieck.JVODEControllerType
JVODEController(; qmin, qmax, qsteady_min, qsteady_max, gamma,
                qmax_first_step, failfactor)

Step-size controller for the variable-order Nordsieck-form JVODE family. Composes the standard step-size knobs via CommonControllerOptions; the adaptive logic is integrated into the algorithm itself, so the controller cache falls through to alg-level dispatch but the knobs are exposed as a real, settable controller (instead of being hard-wired on the algorithm struct as before).

source
OrdinaryDiffEqExtrapolation.ExtrapolationControllerType
ExtrapolationController(; qmin, qmax, gamma, qsteady_min, qsteady_max,
                        qmax_first_step, failfactor)
ExtrapolationController(::Type{QT}, alg; kwargs...)

Step-size controller for the extrapolation family (ExtrapolationMidpointDeuflhard, ExtrapolationMidpointHairerWanner, ImplicitDeuflhardExtrapolation, ImplicitHairerWannerExtrapolation, ImplicitEulerExtrapolation, ImplicitEulerBarycentricExtrapolation). Composes the standard step-size knobs via CommonControllerOptions; the per-step order-and-stepsize selection (Deuflhard / Hairer–Wanner style "work per step" minimization) lives in the algorithm-specific stepsize_controller_internal! / step_accept_controller! / step_reject_controller! methods, so the controller-cache dispatch falls through to alg-level logic.

source
ImplicitDiscreteSolve.KantorovichTypeControllerType
KantorovichTypeController

Default controller for implicit discrete solvers. Assuming a Newton method is used to solve the nonlinear problem, this controller uses convergence rate estimates to adapt the step size based on an posteriori estimate of the change in the Newton convergence radius between steps as described below.

Given the convergence rate estimate Θ₀ for the first iteration, the step size controller adapts the time step as dtₙ₊₁ = γ * (g(Θbar) / (g(Θ₀)))^(1 / p) dtₙ with g(x) = √(1 + 4x) - 1. p denotes the order of the solver – i.e. the order of the extrapolation algorithm to compute the initial guess for the solve at tₙ given a solution at tₙ₋₁ – and Θbar denotes the target convergence rate. γ is a safety factor.

A factor Θreject controls the rejection of a time step if any Θₖ > Θreject. In this case the first Θₖ violating this criterion is taken and the time step is adapted such that dtₙ₊₁ = γ * (g(Θbar) / (g(Θk)))^(1 / p) dtₙ. This behavior can be changed by setting strict = false. In this case the step is accepted whenever the Newton solver converges.

The controller furthermore limits the growth and shrinkage of the time step by a factor between qmin and qmax.

The baseline algorithm has been derived in Peter Deuflhard's book "Newton Methods for Nonlinear Problems" in Section 5.1.3 (Adaptive pathfollowing algorithms). Please note that some implementation details deviate from the original algorithm.

source

Composite controllers

OrdinaryDiffEqCore.CompositeControllerType
CompositeController(controllers::Tuple)

Controller used by CompositeAlgorithm — the per-step dispatch walks controllers[integrator.cache.current] so each branch of the composite algorithm keeps its own controller and per-step state. The default policy is to construct one of these from the per-branch default_controller entries, but a user can supply any tuple of AbstractControllers of matching length.

source

CommonControllerOptions

Every controller embeds a CommonControllerOptions{T} as its basic field. This holds the seven standard step-size scalars and is the target of resolve_basic.

OrdinaryDiffEqCore.CommonControllerOptionsType
CommonControllerOptions{T}

Type-stable, fully-resolved bundle of the standard step-size knobs every adaptive controller uses. All seven fields have concrete element type T (no Union{Nothing, T}), so reads in the hot path (stepsize_controller!, accessors like get_qmin) are inferable.

The fields are:

  • qmin / qmax: lower / upper bounds on the per-step shrink/grow factor.
  • qmax_first_step: looser qmax applied to the very first accepted step (mirrors Sundials CVODE — the initial dt from automatic step-size selection is approximate, so a much larger growth is allowed once).
  • gamma: safety factor applied to the controller's predicted dt change (used by IController, PIController, PredictiveController, and the algorithm-specific BDF/JVODE controllers; PIDController ignores gamma and limits dt via its limiter / accept_safety instead).
  • qsteady_min / qsteady_max: deadband — if the proposed factor lies inside this interval, dt is held constant.
  • failfactor: post-Newton-failure shrink factor used by post_newton_controller!.

User-supplied overrides flow through controllers as a NamedTuple of keyword arguments (whatever subset the user passed). At setup_controller_cache time, resolve_basic constructs a fresh CommonControllerOptions{QT} by filling in unset entries from the algorithm-specific defaults (qmin_default(alg), qmax_default(alg), …). So a user-constructed BDFController() picks up BDF-tuned defaults (qmax = 5, qsteady_min = 9//10, qsteady_max = 12//10, …) while a user-constructed IController() falls back to the generic defaults (qmin = 1//5, qmax = 10, …). This matches the historical behavior of the per-algorithm step-size knobs that used to live on the OrdinaryDiffEq algorithm structs themselves.

source
OrdinaryDiffEqCore.resolve_basicFunction
resolve_basic(overrides::NamedTuple, alg, ::Type{QT})
resolve_basic(opts::CommonControllerOptions, alg, ::Type{QT})

Return a fully-resolved CommonControllerOptions{QT} by combining the user-supplied overrides (a NamedTuple containing whatever subset of qmin, qmax, qmax_first_step, gamma, qsteady_min, qsteady_max, failfactor the user passed) with the algorithm-specific defaults (qmin_default(alg), qmax_default(alg), …). Each entry is converted to QT. Called from each controller's setup_controller_cache method so that user-constructed controllers (e.g. BDFController(qmax = 20)) get sensible per-algorithm defaults for the knobs they didn't set.

The CommonControllerOptions-input form lets already-resolved controllers be re-resolved against a (possibly new) element type QT without going through the override mechanism — useful when a composite algorithm re-keys controller caches to a different scalar type.

source

Reading knobs from the integrator

Hot-path code (the integrator's rejection logic, the implicit Newton loop, the per-step controller logic itself) reads knobs through these accessor functions rather than directly off the controller struct. That layer of indirection lets CompositeControllerCache delegate to the currently-active sub-cache.

OrdinaryDiffEqCore.get_qminFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_qmaxFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_qmax_first_stepFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_gammaFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_qsteady_minFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_qsteady_maxFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_failfactorFunction
get_qmin(integrator)
get_qmax(integrator)
get_qmax_first_step(integrator)
get_gamma(integrator)
get_qsteady_min(integrator)
get_qsteady_max(integrator)
get_failfactor(integrator)

Read a step-size knob from the integrator's controller. Default dispatch reads integrator.controller_cache.controller.basic.X — i.e. it goes through the CommonControllerOptions embedded on every concrete controller (IController/PIController/PIDController/ PredictiveController/BDFController/JVODEController).

CompositeControllerCache overrides each accessor to delegate to the currently active sub-cache (mirroring how stepsize_controller! and friends dispatch). The transitional DummyControllerCache also provides overrides for the BDF/Nordsieck cases that haven't been migrated yet.

These accessors are what the integrator-level paths (e.g. the isoutofdomain rejection path for qmin, post_newton_controller! for failfactor) call instead of reading integrator.opts.X — the v7 controller refactor moved these knobs off DEOptions and onto the controller object.

source
OrdinaryDiffEqCore.get_current_qmaxFunction
get_current_qmax(integrator, qmax)

Return the effective maximum step size growth factor for the current step. On the first step (before any successful steps), returns qmax_first_step from the controller (defaulting to 10000). This allows a much larger step size increase on the first step since the initial dt from the automatic step size selection algorithm is only approximate.

This mirrors the behavior of Sundials CVODE, which limits hnew/hold to 10 on normal steps but 10^4 on the first step.

See also: https://github.com/SciML/DifferentialEquations.jl/issues/299

source
OrdinaryDiffEqCore.get_EEstFunction
get_EEst(controller_cache)
get_EEst(integrator)

Return the scalar error estimate stored by the controller cache. When passed an integrator, forwards to its controller_cache. The fallback reads the EEst field on the cache; controllers that represent the error estimate differently can override this method.

source
OrdinaryDiffEqCore.set_EEst!Function
set_EEst!(controller_cache, val)
set_EEst!(integrator, val)

Store val as the scalar error estimate on the controller cache. When passed an integrator, forwards to its controller_cache. The fallback writes the EEst field on the cache; controllers that represent the error estimate differently can override this method.

source

Per-algorithm defaults

Algorithms specialize these single-argument functions to set their own default values for the standard step-size knobs. The defaults are consulted by resolve_basic at setup_controller_cache time for any knob the user did not override.

OrdinaryDiffEqCore.qmin_defaultFunction
qmin_default(alg)
qmax_default(alg)

Algorithm-specific default for the lower / upper bound on the per-step shrink/grow factor. Used by resolve_basic to fill in CommonControllerOptions when the user didn't override. Override these for a new algorithm type to change its defaults. The generic fallbacks are qmin = 1 // 5, qmax = 10 for adaptive algorithms (qmin = 0 for non-adaptive).

source
OrdinaryDiffEqCore.qmax_defaultFunction
qmin_default(alg)
qmax_default(alg)

Algorithm-specific default for the lower / upper bound on the per-step shrink/grow factor. Used by resolve_basic to fill in CommonControllerOptions when the user didn't override. Override these for a new algorithm type to change its defaults. The generic fallbacks are qmin = 1 // 5, qmax = 10 for adaptive algorithms (qmin = 0 for non-adaptive).

source
OrdinaryDiffEqCore.qmax_first_step_defaultFunction
qmax_first_step_default(alg)

Algorithm-specific default for the looser qmax applied to the very first accepted step (mirrors Sundials CVODE — the initial dt from the automatic step-size selection is approximate, so a much larger growth is allowed once). Generic fallback is 10000. See https://github.com/SciML/DifferentialEquations.jl/issues/299.

source
OrdinaryDiffEqCore.gamma_defaultFunction
gamma_default(alg)

Algorithm-specific default for the safety factor on the controller's predicted dt change. The generic fallback is 9 // 10 for adaptive algorithms and 0 otherwise.

source
OrdinaryDiffEqCore.qsteady_min_defaultFunction
qsteady_min_default(alg)
qsteady_max_default(alg)

Algorithm-specific defaults for the deadband interval — if the proposed step-size factor lies inside [qsteady_min, qsteady_max], the controller holds dt constant. Generic fallback is 1 for both (no deadband). Adaptive implicit algorithms widen qsteady_max to 6 // 5 to reduce Jacobian recomputation; BDF widens both (9 // 10, 12 // 10).

source
OrdinaryDiffEqCore.qsteady_max_defaultFunction
qsteady_min_default(alg)
qsteady_max_default(alg)

Algorithm-specific defaults for the deadband interval — if the proposed step-size factor lies inside [qsteady_min, qsteady_max], the controller holds dt constant. Generic fallback is 1 for both (no deadband). Adaptive implicit algorithms widen qsteady_max to 6 // 5 to reduce Jacobian recomputation; BDF widens both (9 // 10, 12 // 10).

source
OrdinaryDiffEqCore.beta1_defaultFunction
beta1_default(alg, beta2)
beta2_default(alg)

Algorithm-specific defaults for the proportional / integral gains of the PIController. Defaults are scaled by the algorithm's order (beta2 = 2 // (5 alg_order(alg)), beta1 = 7 // (10 alg_order(alg))) for adaptive algorithms and 0 otherwise. beta1 is computed from beta2, so override beta2_default first if you want both to change.

source
OrdinaryDiffEqCore.beta2_defaultFunction
beta1_default(alg, beta2)
beta2_default(alg)

Algorithm-specific defaults for the proportional / integral gains of the PIController. Defaults are scaled by the algorithm's order (beta2 = 2 // (5 alg_order(alg)), beta1 = 7 // (10 alg_order(alg))) for adaptive algorithms and 0 otherwise. beta1 is computed from beta2, so override beta2_default first if you want both to change.

source
OrdinaryDiffEqCore.default_controllerFunction
default_controller(::Type{QT}, alg)

Return the step-size controller used by alg when the user does not pass one explicitly via solve(prob, alg; controller = ...). QT is the scalar type used internally for q, dt, and the step-size factors — usually Float64. The generic fallback is PIController(QT, alg); override for new algorithm types whose default should differ (e.g. BDF uses BDFController, JVODE uses JVODEController).

source

For example, the BDF family overrides qmax_default(::QNDF) = 5 // 1 and qsteady_max_default(::QNDF) = 12 // 10 to match the historical QNDF(qmax=5, qsteady_max=12//10) defaults.

Implementing a custom controller

To plug in a new step-size strategy, you implement two types — the controller (immutable knob holder) and the controller cache (mutable per-solve state) — and a small set of methods.

Required types

struct MyController{B} <: OrdinaryDiffEqCore.AbstractController
    basic::B  # NamedTuple of overrides (unresolved) or CommonControllerOptions{T} (resolved)
end

mutable struct MyControllerCache{T, E} <: OrdinaryDiffEqCore.AbstractControllerCache
    controller::MyController{OrdinaryDiffEqCore.CommonControllerOptions{T}}
    # ... per-step scratch state owned by the controller ...
    EEst::E   # scalar error estimate; required unless you override `get_EEst` / `set_EEst!`
end

The cache must expose the scalar error estimate via an EEst field, or provide its own get_EEst / set_EEst! overrides.

Required methods

OrdinaryDiffEqCore.AbstractControllerType
AbstractController

Supertype of every step-size controller. A concrete subtype is a small holder for the controller's tuning knobs — typically just a basic field holding a CommonControllerOptions (resolved) or a NamedTuple of user-supplied overrides (unresolved, awaiting resolve_basic) plus any controller-specific scalars (like PIController's beta1 / beta2).

Per-solve mutable state (q11, qold, dtacc, erracc, the scalar error estimate, …) lives on the cache subtype, not on the controller itself. See AbstractControllerCache.

Required methods to implement when adding a new controller: setup_controller_cache, stepsize_controller!, step_accept_controller!, step_reject_controller!. Optional methods: accept_step_controller, post_newton_controller!, reinit_controller!, sync_controllers!, reset_alg_dependent_opts!.

source
OrdinaryDiffEqCore.AbstractControllerCacheType
AbstractControllerCache

Each controller cache is expected to own the scalar error estimate used by the step-size logic and exposed to users as get_EEst(integrator). The accessors get_EEst and set_EEst! read/write that scalar; the default implementations dispatch on a EEst field. Controllers that track multiple error estimates (e.g. BDF methods with EEst1, EEst2) can either still hold a canonical scalar EEst or override the accessors directly.

source
OrdinaryDiffEqCore.setup_controller_cacheFunction
setup_controller_cache(alg, algcache, controller::AbstractController, ::Type{EEstT})::AbstractControllerCache

This function takes a controller together with the time stepping algorithm to construct and initialize the respective cache for the controller. The EEstT type parameter is the element type of the error estimate stored on the returned cache (matches what used to live on get_EEst(integrator)).

source
OrdinaryDiffEqCore.stepsize_controller!Function
stepsize_controller!(integrator, alg)
stepsize_controller!(integrator, controller_cache::AbstractControllerCache, alg)

Update the cache to compute the new order of the time marching algorithm and prepare for an update of the time step. The update can be either due to a rejection or acceptance of the current time step.

source
OrdinaryDiffEqCore.step_accept_controller!Function
step_accept_controller!(integrator, alg, q)
step_accept_controller!(integrator, controller_cache::AbstractControllerCache, alg, q)

This function gets called in case of an accepted time step right after stepsize_controller!. It returns the proposed new time step length. Please note that the time step length might not be applied as is and subject to further modification to e.g. match the next time stop.

source

In summary:

MethodWhat it doesWhen called
setup_controller_cache(alg, alg_cache, controller, ::Type{EEstT})Allocate the controller cache, resolve any unresolved CommonControllerOptions.init, once per solve.
stepsize_controller!(integrator, cache, alg)Compute and return the step-size factor q based on the current error estimate.Every step, after perform_step!.
step_accept_controller!(integrator, cache, alg, q)Return the new dt for the next step when the current one is accepted.Once per accepted step.
step_reject_controller!(integrator, cache, alg)Mutate integrator.dt so the rejected step can be retried.Once per rejected step.

Optional methods

These have library-provided defaults that work for most controllers but are commonly overridden for algorithm-integrated controllers (BDF, JVODE, …) where the algorithm itself owns the step-size logic.

OrdinaryDiffEqCore.accept_step_controllerFunction
accept_step_controller(integrator, alg)::Bool
accept_step_controller(integrator, controller_cache::AbstractControllerCache, alg)::Bool

This function decides whether the current time step should be accepted or rejected. A return value of false corresponds to a rejection.

source
OrdinaryDiffEqCore.post_newton_controller!Function
post_newton_controller!(integrator, alg)
post_newton_controller!(integrator, cache::AbstractControllerCache, alg)

Hook called by implicit-solver perform-step routines after a Newton iteration fails to converge. The default behavior is to shrink integrator.dt by get_failfactor(integrator) so the step is retried with a smaller dt. BDF / Nordsieck controllers override this to also reduce the working order on repeated Newton failures.

source
OrdinaryDiffEqCore.reinit_controller!Function
reinit_controller!(integrator, cache::AbstractControllerCache)

Hook called from reinit!(integrator) so the controller can reset its per-solve state (q11, qold, dtacc, erracc, error history, …) to the same values it would have at init. Default is a no-op.

source
OrdinaryDiffEqCore.sync_controllers!Function
sync_controllers!(cache_dst::AbstractControllerCache, cache_src::AbstractControllerCache)

Copy per-solve state from cache_src to cache_dst. Used by CompositeController when the active branch changes so the new branch's controller starts from the previous branch's state (e.g. copying qold and q11 from a PIController to a fresh IController). Default is a no-op; override when your controller carries scratch state worth preserving across switches.

source
OrdinaryDiffEqCore.reset_alg_dependent_opts!Function
reset_alg_dependent_opts!(cache::AbstractControllerCache, alg_old, alg_new)

Hook called when a CompositeAlgorithm switches branches: gives the controller cache a chance to re-derive any algorithm-specific cached values (e.g. order-scaled beta1 / beta2). Default is a no-op.

source

If no user default, then this will change the default to the defaults for the second algorithm. Except if the user default turns out to be the default for the current alg, then it will change anyway and keep changing afterwards (e.g. adaptive).

source
MethodDefaultWhen to override
accept_step_controller(integrator, cache, alg)::Boolget_EEst(integrator) <= 1Controllers that accept on a different criterion (e.g. PIDController uses dt_factor >= accept_safety).
post_newton_controller!(integrator, cache, alg)Shrinks integrator.dt by get_failfactor(integrator)Implicit-solver controllers that also reduce the BDF order on Newton failure (BDFController, JVODEController).
reinit_controller!(integrator, cache)no-opStateful controllers that need to reset q11, qold, dtacc, etc. on reinit!.
sync_controllers!(cache1, cache2)no-opComposite-alg switching, when state must transfer between sub-controllers.
reset_alg_dependent_opts!(cache, alg1, alg2)no-opRe-derive beta1 / beta2 etc. when a composite algorithm switches between branches.
get_qmin(integrator, cache) (and get_qmax / get_gamma / get_qsteady_min / get_qsteady_max / get_qmax_first_step / get_failfactor)Reads cache.controller.basic.XControllers that don't embed a CommonControllerOptions or store the knobs somewhere other than cache.controller.basic.
get_EEst(cache) / set_EEst!(cache, val)Read/write the EEst field on the cacheCaches that store the error differently (vector, per-stage, …).

Worked example

Here is a stripped-down implementation of an I-controller using the new interface. It composes CommonControllerOptions so users can pass qmin / qmax / gamma as kwargs, and it picks up algorithm-specific defaults automatically.

A note on conventions before the code: the controller protocol carries q as a divisor — the new step size is dt / q, so q > 1 shrinks the step, q < 1 grows it, and the controller clamps q to [inv(qmax), inv(qmin)]. Every built-in controller follows this; new controllers should too so they compose with each other and with CompositeController.

using OrdinaryDiffEqCore
using OrdinaryDiffEqCore: AbstractController, AbstractControllerCache,
                         CommonControllerOptions, resolve_basic,
                         _resolved_QT, get_EEst, get_current_adaptive_order,
                         get_current_qmax, fastpower
import OrdinaryDiffEqCore: setup_controller_cache, stepsize_controller!,
                          step_accept_controller!, step_reject_controller!

# Type-stability trick: `basic` is either a `NamedTuple` of user-supplied
# overrides (before setup_controller_cache resolves it) or a fully-resolved
# `CommonControllerOptions{T}` (after). The Union constraint lets the
# compiler still specialize on both branches.
struct MyIController{B <: Union{NamedTuple, CommonControllerOptions}} <: AbstractController
    basic::B
end

# Keyword form: the user's `qmin = …`, `qmax = …`, … ride along as a
# `NamedTuple` until setup_controller_cache merges them with the
# algorithm's defaults.
MyIController(; kwargs...) = MyIController(NamedTuple(kwargs))

# The cache holds the resolved controller (concrete `CommonControllerOptions{T}`)
# and any per-solve scratch state — for an I-controller, just the scalar EEst.
mutable struct MyIControllerCache{T, E} <: AbstractControllerCache
    controller::MyIController{CommonControllerOptions{T}}
    EEst::E
end

# Called once at `init` time. Picks the scalar type `QT` for the step-size
# factors, then fills in every unset knob from `qmin_default(alg)`,
# `qmax_default(alg)`, … via `resolve_basic`. After this, every read of
# `cache.controller.basic.X` is type-stable.
function setup_controller_cache(alg, alg_cache, controller::MyIController, ::Type{E}) where {E}
    QT = _resolved_QT(controller.basic)
    basic = resolve_basic(controller.basic, alg, QT)
    return MyIControllerCache{QT, E}(MyIController(basic), oneunit(E))
end

# Called after each `perform_step!`. Computes the divisor `q` from the
# scaled error estimate using the standard formula
# `q = EEst^(1/(p+1)) / gamma`, clamped to `[inv(qmax), inv(qmin)]`.
# `get_current_qmax` returns a looser `qmax_first_step` on the very first
# step (CVODE-style) and the regular `qmax` afterwards.
@inline function stepsize_controller!(integrator, cache::MyIControllerCache, alg)
    (; qmin, qmax, gamma) = cache.controller.basic
    qmax = get_current_qmax(integrator, qmax)
    EEst = get_EEst(integrator)             # scaled error estimate from perform_step!
    if iszero(EEst)
        return inv(qmax)                    # grow as fast as allowed
    end
    expo = 1 / (get_current_adaptive_order(alg, integrator.cache) + 1)
    qtmp = fastpower(EEst, expo) / gamma    # raw factor predicted by I controller
    return max(inv(qmax), min(inv(qmin), qtmp))
end

# Called once per accepted step. Returns the new `dt` for the next step.
# The deadband (`qsteady_min ≤ q ≤ qsteady_max`) holds `dt` constant when
# the proposed change is small — useful for avoiding jitter, and required
# by Jacobian-reusing implicit solvers.
function step_accept_controller!(integrator, cache::MyIControllerCache, alg, q)
    (; qsteady_min, qsteady_max) = cache.controller.basic
    if qsteady_min <= q <= qsteady_max
        q = one(q)                          # within deadband — keep dt
    end
    return integrator.dt / q                # new dt for the next attempt
end

# Called once per rejected step. Must mutate `integrator.dt` so the
# step gets retried with a smaller value. Here we just shrink by `qmax`
# (the most aggressive contraction the controller is allowed to make).
step_reject_controller!(integrator, cache::MyIControllerCache, alg) =
    (integrator.dt = integrator.dt / get_current_qmax(integrator, cache.controller.basic.qmax))

You can then pass it through solve:

sol = solve(prob, Tsit5(); controller = MyIController(qmax = 5))

Tsit5() is non-stiff so it would normally get PIController; here we swap it for the custom one. The qmax = 5 rides through as part of the unresolved NamedTuple; resolve_basic fills in qmin, gamma, qsteady_min, qsteady_max, qmax_first_step, and failfactor from qmin_default(::Tsit5), gamma_default(::Tsit5), etc., at init time.

Adding a controller-specific knob

If your controller relies on a knob that doesn't exist on CommonControllerOptions, add it as a new field on your controller struct (alongside basic). The built-in controllers follow this pattern: PIController keeps beta1 / beta2 / qoldinit separate from the standard knobs, and PIDController keeps beta / accept_safety / limiter. Separating them keeps CommonControllerOptions small enough to be useful as the "every adaptive controller has these" surface, while still letting specialized controllers carry their own tuning state. The new field moves through setup_controller_cache the same way basic does — unresolved on the user-constructed controller, then converted to the final type QT and stored on the resolved controller when the cache is built.