GaussianIVPFilter

class probnum.diffeq.GaussianIVPFilter(ivp, prior, measurement_model, with_smoothing, init_implementation, initrv=None)[source]

Bases: probnum.diffeq.ODESolver

ODE solver that uses a Gaussian filter.

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
  • ivp (IVP) – Initial value problem to be solved.

  • prior (Integrator) – Prior distribution.

  • measurement_model (DiscreteGaussian) – ODE measurement model.

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

  • init_implementation (Callable[[Callable, ndarray, float, Integrator, Normal, Optional[Callable]], Normal]) – Initialization algorithm. Either via Scipy (initialize_odefilter_with_rk) or via Taylor-mode AD (initialize_odefilter_with_taylormode). For more convenient construction, consider GaussianIVPFilter.construct_with_rk_init() and GaussianIVPFilter.construct_with_taylormode_init().

  • initrv (Optional[Normal]) – Initial random variable to be used by the filter. Optional. Default is None, which amounts to choosing a standard-normal initial RV. As part of GaussianIVPFilter.initialise() (called in GaussianIVPFilter.solve()), this variable is improved upon with the help of the initialisation algorithm. The influence of this choice on the posterior may vary depending on the initialization strategy, but is almost always weak.

Attributes Summary

prior

Methods Summary

construct_with_rk_init(ivp, prior, …[, …])

Create a Gaussian IVP filter that is initialised via initialize_odefilter_with_rk().

construct_with_taylormode_init(ivp, prior, …)

Create a Gaussian IVP filter that is initialised via initialize_odefilter_with_taylormode().

initialise()

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

method_callback(time, current_guess, …)

Optional callback.

postprocess(odesol)

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

rvlist_to_odesol(times, rvs)

Create an ODESolution object.

solve(steprule)

Solve an IVP.

step(t, t_new, current_rv)

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

string_to_measurement_model(…[, …])

Construct a measurement model \(\mathcal{N}(g(m), R)\) for an ODE.

Attributes Documentation

prior

Methods Documentation

classmethod construct_with_rk_init(ivp, prior, measurement_model, with_smoothing, initrv=None, init_h0=0.01, init_method='DOP853')[source]

Create a Gaussian IVP filter that is initialised via initialize_odefilter_with_rk().

classmethod construct_with_taylormode_init(ivp, prior, measurement_model, with_smoothing, initrv=None)[source]

Create a Gaussian IVP filter that is initialised via initialize_odefilter_with_taylormode().

initialise()[source]

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

method_callback(time, current_guess, current_error)

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).

postprocess(odesol)[source]

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

rvlist_to_odesol(times, rvs)[source]

Create an ODESolution object.

solve(steprule)

Solve an IVP.

Parameters

steprule (StepRule) – Step-size selection rule, e.g. constant steps or adaptive steps.

step(t, t_new, current_rv)[source]

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

static string_to_measurement_model(measurement_model_string, ivp, prior, measurement_noise_covariance=0.0)[source]

Construct a measurement model \(\mathcal{N}(g(m), R)\) for an ODE.

Return a DiscreteGaussian (either a DiscreteEKFComponent or a DiscreteUKFComponent) that provides a tractable approximation of the transition densities based on the local defect of the ODE

\[g(m) = H_1 m(t) - f(t, H_0 m(t))\]

and user-specified measurement noise covariance \(R\). Almost always, the measurement noise covariance is zero.

Compute either type filter, each with a different interpretation of the Jacobian \(J_g\):

  • EKF0 thinks \(J_g(m) = H_1\)

  • EKF1 thinks \(J_g(m) = H_1 - J_f(t, H_0 m(t)) H_0^\top\)

  • UKF thinks: ‘’What is a Jacobian?’’ and uses the unscented transform to compute a tractable approximation of the transition densities.