|
mrs_lib
Various reusable classes, functions and utilities for use in MRS projects
|
Go to the documentation of this file.
19 template <
int n_states,
int n_inputs,
int n_measurements>
24 template <
int n_states,
int n_inputs,
int n_measurements>
26 : m_alpha(alpha), m_kappa(kappa), m_beta(beta), m_Wm(
W_t::Zero()), m_Wc(
W_t::Zero()), m_transition_model(transition_model), m_observation_model(observation_model)
36 template <
int n_states,
int n_inputs,
int n_measurements>
41 m_lambda = m_alpha*m_alpha*(double(n) + m_kappa) -
double(n);
44 m_Wm(0) = m_lambda / (double(n) + m_lambda);
45 m_Wc(0) = m_Wm(0) + (1.0 - m_alpha*m_alpha + m_beta);
48 for (
int i = 1; i < w; i++)
50 m_Wm(i) = (1.0 - m_Wm(0))/(w - 1.0);
59 template <
int n_states,
int n_inputs,
int n_measurements>
74 template <
int n_states,
int n_inputs,
int n_measurements>
77 m_transition_model = transition_model;
84 template <
int n_states,
int n_inputs,
int n_measurements>
87 m_observation_model = observation_model;
93 template <
int n_states,
int n_inputs,
int n_measurements>
97 const P_t Pa = (double(n) + m_lambda)*P;
100 Eigen::LLT<P_t> llt(Pa);
101 if (llt.info() != Eigen::Success)
103 P_t tmp = Pa + (double(n) + m_lambda)*1e-9*P_t::Identity();
105 if (llt.info() != Eigen::Success)
107 ROS_WARN(
"UKF: taking the square root of covariance during sigma point generation failed.");
108 throw square_root_exception();
112 const P_t Pa_sqrt = llt.matrixL();
118 template <
int n_states,
int n_inputs,
int n_measurements>
121 Eigen::ColPivHouseholderQR<Pzz_t> qr(Pzz);
122 if (!qr.isInvertible())
125 Pzz_t tmp = Pzz + 1e-9*Pzz_t::Identity(Pzz.rows(), Pzz.cols());
127 if (!qr.isInvertible())
130 ROS_ERROR(
"UKF: could not compute matrix inversion!!! Fix your covariances (the measurement's is probably too low...)");
131 throw inverse_exception();
133 ROS_WARN(
"UKF: artificially inflating matrix for inverse computation! Check your covariances (the measurement's might be too low...)");
135 Pzz_t ret = qr.inverse();
141 template <
int n_states,
int n_inputs,
int n_measurements>
142 typename UKF<n_states, n_inputs, n_measurements>::K_t UKF<n_states, n_inputs, n_measurements>::computeKalmanGain([[maybe_unused]]
const x_t& x, [[maybe_unused]]
const z_t& inn,
const K_t& Pxz,
const Pzz_t& Pzz)
const
144 const Pzz_t Pzz_inv = computeInverse(Pzz);
145 const K_t K = Pxz * Pzz_inv;
151 template <
int n_states,
int n_inputs,
int n_measurements>
159 const P_t P_sqrt = computePaSqrt(P);
160 const auto xrep = x.replicate(1, n);
163 S.template block<n, n>(0, 1) = xrep + P_sqrt;
166 S.template block<n, n>(0, n+1) = xrep - P_sqrt;
177 template <
int n_states,
int n_inputs,
int n_measurements>
180 const auto& x = sc.x;
181 const auto& P = sc.P;
184 const X_t S = computeSigmas(x, P);
188 for (
int i = 0; i < w; i++)
190 X.col(i) = m_transition_model(S.col(i), u, dt);
199 for (
int i = 0; i < w; i++)
202 x_t tmp = 1.0/w * X.col(i);
209 for (
int i = 0; i < w; i++)
211 ret.P += m_Wc(i) * (X.col(i) - ret.x) * (X.col(i) - ret.x).transpose();
222 template <
int n_states,
int n_inputs,
int n_measurements>
225 const auto& x = sc.x;
226 const auto& P = sc.P;
227 const X_t S = computeSigmas(x, P);
230 Z_t Z_exp(z.rows(), w);
231 for (
int i = 0; i < w; i++)
233 Z_exp.col(i) = m_observation_model(S.col(i));
237 z_t z_exp = z_t::Zero(z.rows());
238 for (
int i = 0; i < w; i++)
240 z_exp += m_Wm(i) * Z_exp.col(i);
244 Pzz_t Pzz = Pzz_t::Zero(z.rows(), z.rows());
245 for (
int i = 0; i < w; i++)
247 Pzz += m_Wc(i) * (Z_exp.col(i) - z_exp) * (Z_exp.col(i) - z_exp).transpose();
252 K_t Pxz = K_t::Zero(n, z.rows());
253 for (
int i = 0; i < w; i++)
255 Pxz += m_Wc(i) * (S.col(i) - x) * (Z_exp.col(i) - z_exp).transpose();
259 const z_t inn = (z - z_exp);
260 const K_t K = computeKalmanGain(sc.x, inn, Pxz, Pzz);
263 if (!K.array().isFinite().all())
265 ROS_ERROR(
"UKF: inverting of Pzz in correction update produced non-finite numbers!!! Fix your covariances (the measurement's is probably too low...)");
272 ret.P = P - K * Pzz * K.transpose();
typename Eigen::Matrix< double, n, w > X_t
State sigma points matrix.
Definition: ukf.h:50
void setConstants(const double alpha, const double kappa, const double beta)
Changes the Unscented Transform parameters.
Definition: ukf.hpp:61
UKF()
Convenience default constructor.
Definition: ukf.hpp:20
Eigen::Matrix< double, p, p > R_t
Measurement noise covariance matrix type .
Definition: kalman_filter.h:39
typename Eigen::Matrix< double, n, p > K_t
Kalman gain matrix.
Definition: ukf.h:53
typename std::function< x_t(const x_t &, const u_t &, double)> transition_model_t
function of the state transition model typedef
Definition: ukf.h:75
virtual statecov_t correct(const statecov_t &sc, const z_t &z, const R_t &R) const override
Implements the state correction step (measurement update).
Definition: ukf.hpp:223
typename Base_class::statecov_t statecov_t
typedef of a helper struct for state and covariance
Definition: ukf.h:73
typename Eigen::Matrix< double, w, 1 > W_t
weights vector (2n+1)*1 typedef
Definition: ukf.h:71
Eigen::Matrix< double, m, 1 > u_t
Input vector type .
Definition: kalman_filter.h:35
typename std::function< z_t(const x_t &)> observation_model_t
function of the observation model typedef
Definition: ukf.h:77
void setObservationModel(const observation_model_t &observation_model)
Changes the observation model function.
Definition: ukf.hpp:85
Eigen::Matrix< double, n, 1 > x_t
State vector type .
Definition: kalman_filter.h:34
All mrs_lib functions, classes, variables and definitions are contained in this namespace.
Definition: attitude_converter.h:29
void setTransitionModel(const transition_model_t &transition_model)
Changes the transition model function.
Definition: ukf.hpp:75
typename Eigen::Matrix< double, p, p > Pzz_t
Pzz helper matrix.
Definition: ukf.h:52
typename Base_class::P_t P_t
state covariance n*n typedef
Definition: ukf.h:65
virtual statecov_t predict(const statecov_t &sc, const u_t &u, const Q_t &Q, double dt) const override
Implements the state prediction step (time update).
Definition: ukf.hpp:178
Eigen::Matrix< double, n, n > Q_t
Process noise covariance matrix type .
Definition: kalman_filter.h:40
Eigen::Matrix< double, p, 1 > z_t
Measurement vector type .
Definition: kalman_filter.h:36
typename Eigen::Matrix< double, p, w > Z_t
Measurement sigma points matrix.
Definition: ukf.h:51
Defines UKF - a class implementing the Unscented Kalman Filter .
Implementation of the Unscented Kalman filter .
Definition: ukf.h:39
is thrown when taking the inverse of a matrix fails during kalman gain calculation
Definition: ukf.h:89