ODEFilter

class probnum.diffeq.odefilter.ODEFilter(*, steprule, prior_process, information_operator=None, approx_strategy=None, with_smoothing=True, init_routine=None, diffusion_model=None, _reference_coordinates=0)

Bases: ODESolver

Probabilistic ODE solver based on Gaussian filtering and smoothing. This is based on continuous-discrete Gaussian filtering.

Note: this is specific for IVPs and does not apply without further considerations to, e.g., BVPs.

Parameters
  • prior_process (MarkovProcess) – Prior Gauss-Markov process.

  • steprule (StepRule) – Step-size-selection rule.

  • information_operator (Optional[ODEInformationOperator]) – Information operator.

  • approx_strategy (Optional[ApproximationStrategy]) – Approximation strategy to turn an intractable information operator into a tractable one.

  • with_smoothing (Optional[bool]) – To smooth after the solve or not to smooth after the solve.

  • init_routine (Optional[InitializationRoutine]) – Initialization algorithm. Either via fitting the prior to a few steps of a Runge-Kutta method (RungeKuttaInitialization) or via Taylor-mode automatic differentiation (:class:TaylorModeInitialization) 1.

  • diffusion_model (Optional[Diffusion]) – Diffusion model. This determines which kind of calibration is used. We refer to Bosch et al. (2020) 2 for a survey.

  • _reference_coordinates (Optional[int]) – Use this state as a reference state to compute the normalized error estimate. Optional. Default is 0 (which amounts to the usual reference state for ODE solvers). Another reasonable choice could be 1, but use this at your own risk!

References

1

Krämer, N., and Hennig, P.. Stable implementation of probabilistic ODE solvers. 2021.

2

Bosch, N., and Hennig, P., and Tronarp, F.. Calibrated adaptive probabilistic ODE solvers. 2021.

Methods Summary

attempt_step(state, dt)

Gaussian IVP filter step as nonlinear Kalman filtering with zero data.

initialize(ivp)

Returns t0 and y0 (for the solver, which might be different to ivp.y0)

method_callback(state)

Optional callback.

perform_full_step(state, initial_dt)

Perform a full ODE solver step.

postprocess(odesol)

If specified (at initialisation), smooth the filter output.

rvlist_to_odesol(times, rvs)

Create an ODESolution object.

solution_generator(ivp[, stop_at, callbacks])

Generate ODE solver steps.

solve(ivp[, stop_at, callbacks])

Solve an IVP.

Methods Documentation

attempt_step(state, dt)[source]

Gaussian IVP filter step as nonlinear Kalman filtering with zero data.

It goes as follows:

1. The current state \(x(t) \sim \mathcal{N}(m(t), P(t))\) is split into a deterministic component and a noisy component,

\[x(t) = m(t) + p(t), \quad p(t) \sim \mathcal{N}(0, P(t))\]

which is required for accurate calibration and error estimation: in order to only work with the local error, ODE solvers often assume that the error at the previous step is zero.

2. The deterministic component is propagated through dynamics model and measurement model

\[\hat{z}(t+\Delta t)\sim\mathcal{N}(H\Phi(\Delta t)m(t),H Q(\Delta t) H^\top)\]

which is a random variable that estimates the expected local defect \(\dot y - f(y)\).

3. The counterpart of \(\hat{z}\) in the (much more likely) interpretation of an erronous previous state (recall that \(\hat z\) was based on the interpretation of an error-free previous state) is computed as,

\[z(t + \Delta t) \sim \mathcal{N}(\mathbb{E}[\hat{z}(t + \Delta t)], \mathbb{C}[\hat{z}(t+\Delta t)]+H\Phi(\Delta t)P(t)\Phi(\Delta t)^\top H^\top),\]

which acknowledges the covariance \(P(t)\) of the previous state. \(\mathbb{E}\) is the mean, and \(\mathbb{C}\) is the covariance. Both \(z(t + \Delta t)\) and \(\hat z(t + \Delta t)\) give rise to a reasonable diffusion_model estimate. Which one to use is handled by the Diffusion attribute of the solver. At this point already we can compute a local error estimate of the current step.

  1. Depending on the diffusion model, there are two options now:

    4.1. For a piecewise constant diffusion, the covariances are calibrated locally – that is, at each step. In this case we update the predicted covariance and measured covariance with the most recent diffusion estimate.

    4.2. For a constant diffusion, the calibration happens post hoc, and the only step that is carried out here is an assembly of the full predicted random variable (up to now, only its parts were available).

5. With the results of either 4.1. or 4.2. (which both return a predicted RV and a measured RV), we finally compute the Kalman update and return the result. Recall that the error estimate has been computed in the third step.

initialize(ivp)[source]

Returns t0 and y0 (for the solver, which might be different to ivp.y0)

method_callback(state)

Optional callback.

Can be overwritten. Do this as soon as it is clear that the current guess is accepted, but before storing it. No return. For example: tune hyperparameters (sigma).

perform_full_step(state, initial_dt)

Perform a full ODE solver step.

This includes the acceptance/rejection decision as governed by error estimation and steprule.

postprocess(odesol)[source]

If specified (at initialisation), smooth the filter output.

rvlist_to_odesol(times, rvs)[source]

Create an ODESolution object.

solution_generator(ivp, stop_at=None, callbacks=None)

Generate ODE solver steps.

Parameters
solve(ivp, stop_at=None, callbacks=None)

Solve an IVP.

Parameters