Class UKF¶
Defined in File ukf.h
Inheritance Relationships¶
Base Type¶
public lupnt::Filter(Class Filter)
Class Documentation¶
-
class UKF : public lupnt::Filter¶
Unscented Kalman Filter (UKF): predict/update via the unscented transform (sigma-point propagation) instead of linearization.
An alternative to
EKFfor nonlinear dynamics/measurement models where Jacobians are unavailable or linearization error is too large; constructed and driven the same way as otherFiltersubclasses (SetDynamicsFunction/SetProcessNoiseFunction/SetMeasurementFunction, thenPredict/Update).Public Functions
-
UKF() = default¶
-
inline void Initialize(const Real t0, const State &x0, const MatXd &P0)¶
Initialize the UKF’s epoch, state, and covariance (including prior/posterior copies), set the state dimension
n_x_fromx0, and (re)compute the unscented transform weights viaInitializeUkfParams.Called once by application setup code before the first
Predict/Updatecall.- Parameters:
t0 – Initial epoch [s]
x0 – Initial state estimate, size [n_x]
P0 – Initial state covariance, size [n_x x n_x]
-
MatX ComputeSigmaPoints(const State &state, const MatXd &cov)¶
Generate the
2*n_x_ + 1unscented-transform sigma points for a given meanstateand covariancecov: the central pointstate, plusstate +/- / chol((n_x_ + lambda_) * cov)columns.Called by `Predict` (on the prior `x_`/`P_`, propagated through `f_dyn_`) and by `Update` (on the predicted `x_`/`P_`, passed through `f_meas_`). @param state Mean state about which to spread the sigma points, size [n_x] @param cov Covariance matrix, size [n_x x n_x] @return Sigma points, size [n_x x n_sigma_] (n_sigma_ = 2*n_x_ + 1)
-
State UnscentedTransform(const MatX &sigma_points, MatXd &cov_out, const MatXd *Q = nullptr)¶
Recombine a set of (already-propagated) sigma points into a mean and covariance using the unscented-transform weights
w_m_/w_c_, optionally adding a process noise covarianceQ.Called by
Predictto recombine the propagated dynamics sigma points (withQ= the process noise fromf_proc_) into the predicted state mean/covariance.- Parameters:
sigma_points – Propagated sigma points, size [n_x_ x n_sigma_]
cov_out – Output covariance of the recombined state, size [n_x_ x n_x_]
Q – Optional process noise covariance to add, size [n_x_ x n_x_]
- Returns:
Mean of the sigma points, size [n_x_]
-
inline void SetAlpha(double alpha)¶
Set the unscented-transform spread parameter
alpha(typically small, e.g. 1e-3) and recompute the sigma-point weights.
-
inline void SetBeta(double beta)¶
Set the unscented-transform secondary scaling parameter
beta(2 is optimal for Gaussian distributions) and recompute the sigma-point weights.
-
inline void SetKappa(double kappa)¶
Set the unscented-transform secondary scaling parameter
kappa(typically 0 or3 - n_x) and recompute the sigma-point weights.
-
virtual void Predict(Real t, const State *u = nullptr) override¶
UKF predict step: generate sigma points from
(x_, P_)viaComputeSigmaPoints, propagate each throughf_dyn_(no Jacobian needed), then recombine viaUnscentedTransform(adding process noiseQ_fromf_proc_) to form the predictedx_prior_/P_prior_, which also become the newx_/P_.
-
virtual void Update(const VecX &z_true) override¶
UKF update step: generate sigma points from the predicted
(x_, P_), transform each throughf_meas_to form the measurement sigma points, recombine them (with weightsw_m_/w_c_) into the predicted measurement mean and innovation covarianceS_, compute the state-measurement cross-covariance, form the Kalman gainK_ = cross_cov * S_^-1, and apply the correctiondx_ = K_ * dy_to updatex_post_/P_post_(which also become the newx_/P_).
-
inline VecXd GetMeasurementResidual()¶
Get the measurement residual
dy = z_true - meas_meanfrom the lastUpdatecall, size [n_z].
-
inline MatXd GetKalmanGain()¶
Get the Kalman gain
K_ = cross_cov * S_^-1from the lastUpdatecall, size [n_x x n_z].
-
inline MatXd GetMeasurementNoiseCov()¶
Get the measurement noise covariance
R_(from the first sigma point’sf_meas_evaluation) from the lastUpdatecall, size [n_z x n_z].
-
inline MatXd GetMeasurementJacobian()¶
Get
H_(unused by the UKF’s sigma-point update; retained for interface compatibility withKalmanFilter).
-
inline int GetMeasurementSize()¶
Get the current measurement dimension (number of rows of
H_; seeGetMeasurementJacobian).
-
inline MatXd GetProcessNoise()¶
Get
Q_(unused by the UKF’s sigma-point predict step, which adds process noise directly inUnscentedTransform; retained for interface compatibility withKalmanFilter).
-
inline MatXd GetStateJacobian()¶
Get
F_(unused by the UKF’s sigma-point predict step, which has no explicit state-transition Jacobian; retained for interface compatibility withKalmanFilter).
-
inline MatXd GetInnovationCov()¶
Get the innovation covariance
S_from the lastUpdatecall, size [n_z x n_z].
-
inline MatXd GetMeasurementCov()¶
Same as
GetMeasurementNoiseCov; get the measurement noise covarianceR_.
-
inline VecXd GetStateCorrection()¶
Get the state correction
dx_ = K_ * dy_from the lastUpdatecall, size [n_x].
-
UKF() = default¶