Unscented Kalman Filter
Unscented Kalman filter inference for continuous-discrete nonlinear Gaussian SSMs.
inference_ukf
¶
UKFHyperParams
¶
Bases: NamedTuple
Lightweight container for UKF hyperparameters.
Default values taken from https://github.com/sbitzer/UKF-exposed
emissions_unscented_kalman_filter(params: ParamsCDNLGSSM, t_states: Float[Array, 'num_timesteps 1'], state_means: Float[Array, 'num_timesteps state_dim'], state_covs: Optional[Float[Array, 'num_timesteps state_dim state_dim']] = None, inputs: Optional[Float[Array, 'num_timesteps input_dim']] = None, filter_hyperparams: UKFHyperParams = UKFHyperParams(), warn: bool = True) -> Tuple[Float[Array, 'num_timesteps emission_dim'], Optional[Float[Array, 'num_timesteps emission_dim emission_dim']]]
¶
Compute the emissions corresponding to the UKF linearization of the model.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
params
|
ParamsCDNLGSSM
|
CD-NLGSSM parameters, containing - dynamics RHS drift function - diffusion coefficient matrix L - Brownian covariance matrix Q |
required |
t_states
|
Float[Array, 'num_timesteps 1']
|
continuous-time specific time instants of states |
required |
state_means
|
Float[Array, 'num_timesteps state_dim']
|
state means at time instants t_states, always required |
required |
state_covs
|
Optional[Float[Array, 'num_timesteps state_dim state_dim']]
|
state covariances at time instants t_states, optional - if None, then we assume that the states are point estimates, and simply push through emission function |
None
|
inputs
|
Optional[Float[Array, 'num_timesteps input_dim']]
|
optional array of inputs, of shape (1 + num_timesteps) \times input_dim - The extra input is needed for the initial emission, i.e., it should be at time t_init |
None
|
filter_hyperparams
|
UKFHyperParams
|
hyper-parameters of the filter |
UKFHyperParams()
|
warn
|
bool
|
whether to issue warnings (e.g., about PSD issues) |
True
|
Returns:
| Name | Type | Description |
|---|---|---|
emissions_mean |
Float[Array, 'num_timesteps emission_dim']
|
mean of emissions |
emissions_covariance |
Optional[Float[Array, 'num_timesteps emission_dim emission_dim']]
|
covariance of emissions, if available |
forecast_unscented_kalman_filter(params: ParamsCDNLGSSM, init_forecast: tfd.Distribution, t_init: Float[Array, '1 1'], t_forecast: Float[Array, 'num_timesteps 1'], filter_hyperparams: UKFHyperParams = UKFHyperParams(), inputs: Optional[Float[Array, 'ntime input_dim']] = None, output_fields: Optional[List[str]] = ['forecasted_state_means', 'forecasted_state_covariances'], warn: bool = True) -> GSSMForecast
¶
Run an Unscented Kalman filter to forecast states.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
params
|
ParamsCDNLGSSM
|
CD-NLGSSM parameters, containing - dynamics RHS drift function - diffusion coefficient matrix L - Brownian covariance matrix Q |
required |
init_forecast
|
Distribution
|
initial distribution to forecast with. |
required |
t_init
|
Float[Array, '1 1']
|
time-instant of the initial condition of forecast |
required |
t_forecast
|
Float[Array, 'num_timesteps 1']
|
continuous-time specific time instants to forecast |
required |
filter_hyperparams
|
UKFHyperParams
|
hyper-parameters of the UKF, related to the approximation order |
UKFHyperParams()
|
inputs
|
Optional[Float[Array, 'ntime input_dim']]
|
optional array of inputs. |
None
|
output_fields
|
Optional[List[str]]
|
list of fields to return |
['forecasted_state_means', 'forecasted_state_covariances']
|
warn
|
bool
|
whether to issue warnings (e.g., about PSD issues) |
True
|
Returns:
| Name | Type | Description |
|---|---|---|
forecast |
GSSMForecast
|
forecast object. |
unscented_kalman_filter(params: ParamsCDNLGSSM, emissions: Float[Array, 'ntime emission_dim'], t_emissions: Optional[Float[Array, 'num_timesteps 1']] = None, filter_hyperparams: UKFHyperParams = UKFHyperParams(), inputs: Optional[Float[Array, 'ntime input_dim']] = None, output_fields: Optional[List[str]] = ['filtered_means', 'filtered_covariances', 'predicted_means', 'predicted_covariances'], warn: bool = True) -> PosteriorGSSMFiltered
¶
Run a unscented Kalman filter to produce the marginal likelihood and filtered state estimates.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
params
|
ParamsCDNLGSSM
|
CD-NLGSSM parameters, containing - dynamics RHS drift function - diffusion coefficient matrix L - Brownian covariance matrix Q |
required |
emissions
|
Float[Array, 'ntime emission_dim']
|
array of observations. |
required |
t_emissions
|
Optional[Float[Array, 'num_timesteps 1']]
|
continuous-time specific time instants of observations: if not None, it is an array |
None
|
filter_hyperparams
|
UKFHyperParams
|
hyper-parameters. |
UKFHyperParams()
|
inputs
|
Optional[Float[Array, 'ntime input_dim']]
|
optional array of inputs. |
None
|
output_fields
|
Optional[List[str]]
|
list of fields to include in the output. |
['filtered_means', 'filtered_covariances', 'predicted_means', 'predicted_covariances']
|
warn
|
bool
|
whether to issue warnings (e.g., about PSD issues) |
True
|
Returns:
| Name | Type | Description |
|---|---|---|
filtered_posterior |
PosteriorGSSMFiltered
|
posterior object. |