State Estimators#
The State Estimator uses sensor information and a Prognostics Model to produce an estimate of system state (which can be used to estimate outputs, event_states, and performance metrics). This state estimate can either be used by itself or as input to a Predictor. A state estimator is typically run each time new information is available.
Here’s an example of its use. In this example we use the unscented kalman filter state estimator and the ThrownObject model.
>>> from progpy.models import ThrownObject
>>> from progpy.state_estimators import UnscentedKalmanFilter
>>>
>>> m = ThrownObject()
>>> initial_state = m.initialize()
>>> filt = UnscentedKalmanFilter(m, initial_state)
>>>
>>> load = {} # No load for ThrownObject
>>> new_data = {'x': 1.8} # Observed state
>>> print('Prior: ', filt.x.mean)
>>> filt.estimate(0.1, load, new_data)
>>> print('Posterior: ', filt.x.mean)
See tutorial and examples for more information and additional features.
Included State Estimators#
The following state estimators are included with this package. A new state estimator can be created by subclassing progpy.state_estimators.StateEstimator. See also: state_estimator_template.py
progpy
- class progpy.state_estimators.ParticleFilter(model, x0, **kwargs)#
Estimates state using a Particle Filter (PF) algorithm.
This class defines logic for a PF using a Prognostics Model (see Prognostics Model Package). This filter uses measurement data with noise to estimate the state of the system using a particles. At each step, particles are predicted forward (with noise). Particles are resampled with replacement from the existing particles according to how well the particles match the observed measurements.
The supported configuration parameters (keyword arguments) for UKF construction are described below:
- Parameters
model (PrognosticsModel) – A prognostics model to be used in state estimation See: Prognostics Model Package
x0 (UncertainData, model.StateContainer, or dict) –
Initial (starting) state, with keys defined by model.states
e.g., x = ScalarData({‘abc’: 332.1, ‘def’: 221.003}) given states = [‘abc’, ‘def’]
- Keyword Arguments
t0 (float, optional) – Starting time (s)
dt (float, optional) – Maximum timestep for prediction in seconds. By default, the timestep dt is the difference between the last and current call of .estimate(). Some models are unstable at larger dt. Setting a smaller dt will force the model to take smaller steps; resulting in multiple prediction steps for each estimate step. Default is the parameters[‘dt’] e.g., dt = 1e-2
num_particles (int, optional) – Number of particles in particle filter
resample_fcn (function, optional) – Resampling function ([weights]) -> [indexes] e.g., filterpy.monte_carlo.residual_resample
- class progpy.state_estimators.UnscentedKalmanFilter(model, x0, **kwargs)#
An Unscented Kalman Filter (UKF) for state estimation
This class defines logic for performing an unscented kalman filter with a Prognostics Model (see Prognostics Model Package). This filter uses measurement data with noise to generate a state estimate and covariance matrix.
The supported configuration parameters (keyword arguments) for UKF construction are described below:
- Parameters
model (PrognosticsModel) – A prognostics model to be used in state estimation See: Prognostics Model Package
x0 (UncertainData, model.StateContainer, or dict) –
Initial (starting) state, with keys defined by model.states
e.g., x = ScalarData({‘abc’: 332.1, ‘def’: 221.003}) given states = [‘abc’, ‘def’]
- Keyword Arguments
alpha (float, optional) – UKF Scaling parameter
beta (float, optional) – UKF Scaling parameter
kappa (float, optional) – UKF Scaling parameter
t0 (float, optional) – Starting time (s)
dt (float, optional) – Maximum timestep for prediction in seconds. By default, the timestep dt is the difference between the last and current call of .estimate(). Some models are unstable at larger dt. Setting a smaller dt will force the model to take smaller steps; resulting in multiple prediction steps for each estimate step. Default is the parameters[‘dt’] e.g., dt = 1e-2
- class progpy.state_estimators.KalmanFilter(model, x0, **kwargs)#
A Kalman Filter (KF) for state estimation
This class defines the logic for performing a kalman filter with a LinearModel (see Prognostics Model Package). This filter uses measurement data with noise to generate a state estimate and covariance matrix.
The supported configuration parameters (keyword arguments) for UKF construction are described below:
- Parameters
model (PrognosticsModel) – A prognostics model to be used in state estimation See: Prognostics Model Package
x0 (UncertainData, model.StateContainer, or dict) –
Initial (starting) state, with keys defined by model.states
e.g., x = ScalarData({‘abc’: 332.1, ‘def’: 221.003}) given states = [‘abc’, ‘def’]
- Keyword Arguments
alpha (float, optional) – KF Scaling parameter. An alpha > 1 turns this into a fading memory filter.
t0 (float, optional) – Starting time (s)
dt (float, optional) – Maximum timestep for prediction in seconds. By default, the timestep dt is the difference between the last and current call of .estimate(). Some models are unstable at larger dt. Setting a smaller dt will force the model to take smaller steps; resulting in multiple prediction steps for each estimate step. Default is the parameters[‘dt’] e.g., dt = 1e-2
Q (list[list[float]], optional) – Kalman Process Noise Matrix
R (list[list[float]], optional) – Kalman Measurement Noise Matrix
State Estimator Interface#
- class progpy.state_estimators.StateEstimator(model, x0, **kwargs)#
Interface class for state estimators
Abstract base class for creating state estimators that perform state estimation. Subclasses must implement this interface. Equivalent to “Observers” in NASA’s Matlab Prognostics Algorithm Library
- Parameters
model (PrognosticsModel) – A prognostics model to be used in state estimation See: Prognostics Model Package
x0 (UncertainData, model.StateContainer, or dict) –
Initial (starting) state, with keys defined by model.states
e.g., x = ScalarData({‘abc’: 332.1, ‘def’: 221.003}) given states = [‘abc’, ‘def’]
- Keyword Arguments
t0 (float) – Initial time at which prediction begins, e.g., 0
dt (float) – Maximum timestep for prediction in seconds. By default, the timestep dt is the difference between the last and current call of .estimate(). Some models are unstable at larger dt. Setting a smaller dt will force the model to take smaller steps; resulting in multiple prediction steps for each estimate step. Default is the parameters[‘dt’] e.g., dt = 1e-2
**kwargs – See state-estimator specific documentation for specific keyword arguments.
- abstract estimate(t: float, u, z, **kwargs) None #
Perform one state estimation step (i.e., update the state estimate, filt.x)
- Parameters
t (float) – Current timestamp in seconds (≥ 0.0) e.g., t = 3.4
u (InputContainer) – Measured inputs, with keys defined by model.inputs. e.g., u = m.InputContainer({‘i’:3.2}) given inputs = [‘i’]
z (OutputContainer) – Measured outputs, with keys defined by model.outputs. e.g., z = m.OutputContainer({‘t’:12.4, ‘v’:3.3}) given outputs = [‘t’, ‘v’]
- Keyword Arguments
dt (float, optional) – Maximum timestep for prediction in seconds. By default, the timestep dt is the difference between the last and current call of .estimate(). Some models are unstable at larger dt. Setting a smaller dt will force the model to take smaller steps; resulting in multiple prediction steps for each estimate step. Default is the parameters[‘dt’] e.g., dt = 1e-2
**kwargs – See state-estimator specific documentation for specific keyword arguments.
Note
This method updates the state estimate stored in filt.x, but doesn’t return the updated estimate. Call filt.x to get the updated estimate.
- abstract property x#
The current estimated state.
Example
state = filt.x