ODESolution¶
-
class
probnum.diffeq.
ODESolution
(times, rvs, solver)[source]¶ Bases:
probnum.filtsmooth.filtsmoothposterior.FiltSmoothPosterior
Gaussian IVP filtering solution of an ODE problem
Parameters: - times (array_like) – Times of the discrete-time solution.
- rvs (
list
ofRandomVariable
) – Estimated states (in the state-space model view) of the discrete-time solution. - solver (
GaussianIVPFilter
) – Solver used to compute the discrete-time solution.
See also
GaussianIVPFilter
- ODE solver that behaves like a Gaussian filter.
KalmanPosterior
- Posterior over states after Gaussian filtering/smoothing.
Examples
>>> from probnum.diffeq import logistic, probsolve_ivp >>> from probnum import random_variables as rvs >>> initrv = rvs.Dirac(0.15) >>> ivp = logistic(timespan=[0., 1.5], initrv=initrv, params=(4, 1)) >>> solution = probsolve_ivp(ivp, method="ekf0", step=0.1) >>> # Mean of the discrete-time solution >>> print(solution.y.mean) [[0.15 ] [0.2076198 ] [0.27932997] [0.3649165 ] [0.46054129] [0.55945475] [0.65374523] [0.73686744] [0.8053776 ] [0.85895587] [0.89928283] [0.92882899] [0.95007559] [0.96515825] [0.97577054] [0.9831919 ]] >>> # Times of the discrete-time solution >>> print(solution.t) [0. 0.1 0.2 0.3 0.4 0.5 0.6 0.7 0.8 0.9 1. 1.1 1.2 1.3 1.4 1.5] >>> # Individual entries of the discrete-time solution can be accessed with >>> print(solution[5]) <(1,) Normal with dtype=float64> >>> print(solution[5].mean) [0.55945475] >>> # Evaluate the continuous-time solution at a new time point t=0.65 >>> print(solution(0.65).mean) [0.69875089]
Attributes Summary
dy
Derivatives of the discrete-time solution t
Times of the discrete-time solution y
Probabilistic discrete-time solution Methods Summary
__call__
(t)Evaluate the time-continuous solution at time t. sample
([t, size])Draw samples from the filtering/smoothing posterior. Attributes Documentation
-
t
¶ Times of the discrete-time solution
Type: np.ndarray
-
y
¶ Probabilistic discrete-time solution
Probabilistic discrete-time solution at times \(t_1, ..., t_N\), as a list of random variables. To return means and covariances use
y.mean
andy.cov
.Type: list
ofRandomVariable
Methods Documentation
-
__call__
(t)[source]¶ Evaluate the time-continuous solution at time t.
KalmanPosterior.__call__ does the main algorithmic work to return the posterior for a given location. All that is left to do here is to (1) undo the preconditioning, and (2) to slice the state_rv in order to return only the rv for the function value.
Parameters: t (float) – Location / time at which to evaluate the continuous ODE solution. Returns: Probabilistic estimate of the continuous-time solution at time t
.Return type: RandomVariable
-
sample
(t=None, size=())[source]¶ Draw samples from the filtering/smoothing posterior.
If nothing is specified, a single sample is drawn (supported on self.locations). If locations are specified, the samples are drawn on those locations. If size is specified, more than a single sample is drawn.
Parameters: - locations (array_like, optional) – Locations on which the samples are wanted. Default is none, which implies that self.location is used.
- size (int or tuple of ints, optional) – Indicates how many samples are drawn. Default is an empty tuple, in which case a single sample is returned.
Returns: Drawn samples. If size has shape (A1, …, Z1), locations have shape (L,), and the state space model has shape (A2, …, Z2), the output has shape (A1, …, Z1, L, A2, …, Z2). For example: size=4, len(locations)=4, dim=3 gives shape (4, 4, 3).
Return type: numpy.ndarray