mrs_lib
Various reusable classes, functions and utilities for use in MRS projects
ukf.hpp
Go to the documentation of this file.
1 // clang: MatousFormat
2 
3 #ifndef UKF_HPP
4 #define UKF_HPP
5 
12 #include <ros/ros.h>
13 #include <mrs_lib/ukf.h>
14 
15 namespace mrs_lib
16 {
17  /* constructor //{ */
18 
19  template <int n_states, int n_inputs, int n_measurements>
21  {
22  }
23 
24  template <int n_states, int n_inputs, int n_measurements>
25  UKF<n_states, n_inputs, n_measurements>::UKF(const transition_model_t& transition_model, const observation_model_t& observation_model, const double alpha, const double kappa, const double beta)
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)
27  {
28  assert(alpha > 0.0);
29  computeWeights();
30  }
31 
32  //}
33 
34  /* computeWeights() //{ */
35 
36  template <int n_states, int n_inputs, int n_measurements>
38  {
39  // initialize lambda
40  /* m_lambda = double(n) * (m_alpha * m_alpha - 1.0); */
41  m_lambda = m_alpha*m_alpha*(double(n) + m_kappa) - double(n);
42 
43  // initialize first terms of the weights
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);
46 
47  // initialize the rest of the weights
48  for (int i = 1; i < w; i++)
49  {
50  m_Wm(i) = (1.0 - m_Wm(0))/(w - 1.0);
51  m_Wc(i) = m_Wm(i);
52  }
53  }
54 
55  //}
56 
57  /* setConstants() //{ */
58 
59  template <int n_states, int n_inputs, int n_measurements>
60  // update the UKF constants
61  void UKF<n_states, n_inputs, n_measurements>::setConstants(const double alpha, const double kappa, const double beta)
62  {
63  m_alpha = alpha;
64  m_kappa = kappa;
65  m_beta = beta;
66 
67  computeWeights();
68  }
69 
70  //}
71 
72  /* setTransitionModel() method //{ */
73 
74  template <int n_states, int n_inputs, int n_measurements>
76  {
77  m_transition_model = transition_model;
78  }
79 
80  //}
81 
82  /* setObservationModel() method //{ */
83 
84  template <int n_states, int n_inputs, int n_measurements>
86  {
87  m_observation_model = observation_model;
88  }
89 
90  //}
91 
92  /* computePaSqrt() method //{ */
93  template <int n_states, int n_inputs, int n_measurements>
95  {
96  // calculate the square root of the covariance matrix
97  const P_t Pa = (double(n) + m_lambda)*P;
98 
99  /* Eigen::SelfAdjointEigenSolver<P_t> es(Pa); */
100  Eigen::LLT<P_t> llt(Pa);
101  if (llt.info() != Eigen::Success)
102  {
103  P_t tmp = Pa + (double(n) + m_lambda)*1e-9*P_t::Identity();
104  llt.compute(tmp);
105  if (llt.info() != Eigen::Success)
106  {
107  ROS_WARN("UKF: taking the square root of covariance during sigma point generation failed.");
108  throw square_root_exception();
109  }
110  }
111 
112  const P_t Pa_sqrt = llt.matrixL();
113  return Pa_sqrt;
114  }
115  //}
116 
117  /* computeInverse() method //{ */
118  template <int n_states, int n_inputs, int n_measurements>
119  typename UKF<n_states, n_inputs, n_measurements>::Pzz_t UKF<n_states, n_inputs, n_measurements>::computeInverse(const Pzz_t& Pzz) const
120  {
121  Eigen::ColPivHouseholderQR<Pzz_t> qr(Pzz);
122  if (!qr.isInvertible())
123  {
124  // add some stuff to the tmp matrix diagonal to make it invertible
125  Pzz_t tmp = Pzz + 1e-9*Pzz_t::Identity(Pzz.rows(), Pzz.cols());
126  qr.compute(tmp);
127  if (!qr.isInvertible())
128  {
129  // never managed to make this happen except for explicitly putting NaNs in the input
130  ROS_ERROR("UKF: could not compute matrix inversion!!! Fix your covariances (the measurement's is probably too low...)");
131  throw inverse_exception();
132  }
133  ROS_WARN("UKF: artificially inflating matrix for inverse computation! Check your covariances (the measurement's might be too low...)");
134  }
135  Pzz_t ret = qr.inverse();
136  return ret;
137  }
138  //}
139 
140  /* computeKalmanGain() method //{ */
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
143  {
144  const Pzz_t Pzz_inv = computeInverse(Pzz);
145  const K_t K = Pxz * Pzz_inv;
146  return K;
147  }
148  //}
149 
150  /* computeSigmas() method //{ */
151  template <int n_states, int n_inputs, int n_measurements>
152  typename UKF<n_states, n_inputs, n_measurements>::X_t UKF<n_states, n_inputs, n_measurements>::computeSigmas(const x_t& x, const P_t& P) const
153  {
154  // calculate sigma points
155  // fill in the middle of the elipsoid
156  X_t S;
157  S.col(0) = x;
158 
159  const P_t P_sqrt = computePaSqrt(P);
160  const auto xrep = x.replicate(1, n);
161 
162  // positive sigma points
163  S.template block<n, n>(0, 1) = xrep + P_sqrt;
164 
165  // negative sigma points
166  S.template block<n, n>(0, n+1) = xrep - P_sqrt;
167 
168  /* std::cout << "x: " << std::endl << x << std::endl; */
169  /* std::cout << "S rowmean: " << std::endl << S.rowwise().mean() << std::endl; */
170 
171  return S;
172  }
173  //}
174 
175  /* predict() method //{ */
176 
177  template <int n_states, int n_inputs, int n_measurements>
179  {
180  const auto& x = sc.x;
181  const auto& P = sc.P;
182  statecov_t ret;
183 
184  const X_t S = computeSigmas(x, P);
185 
186  // propagate sigmas through the transition model
187  X_t X;
188  for (int i = 0; i < w; i++)
189  {
190  X.col(i) = m_transition_model(S.col(i), u, dt);
191  }
192 
193  /* std::cout << "x: " << std::endl << x << std::endl; */
194  /* std::cout << "X rowmean: " << std::endl << X.rowwise().mean() << std::endl; */
195  /* std::cout << "m_Wm sum: " << m_Wm.sum() << std::endl; */
196 
197  // recompute the state vector
198  ret.x = x_t::Zero();
199  for (int i = 0; i < w; i++)
200  {
201  //TODO: WHY DOES THIS SHIT WORK IF I SUBSTITUTE m_Wm(i) FOR 1.0/w ??
202  x_t tmp = 1.0/w * X.col(i);
203  /* x_t tmp = m_Wm(i) * X.col(i); */
204  ret.x += tmp;
205  }
206 
207  // recompute the covariance
208  ret.P = P_t::Zero();
209  for (int i = 0; i < w; i++)
210  {
211  ret.P += m_Wc(i) * (X.col(i) - ret.x) * (X.col(i) - ret.x).transpose();
212  }
213  ret.P += Q;
214 
215  return ret;
216  }
217 
218  //}
219 
220  /* correct() method //{ */
221 
222  template <int n_states, int n_inputs, int n_measurements>
224  {
225  const auto& x = sc.x;
226  const auto& P = sc.P;
227  const X_t S = computeSigmas(x, P);
228 
229  // propagate sigmas through the observation model
230  Z_t Z_exp(z.rows(), w);
231  for (int i = 0; i < w; i++)
232  {
233  Z_exp.col(i) = m_observation_model(S.col(i));
234  }
235 
236  // compute expected measurement
237  z_t z_exp = z_t::Zero(z.rows());
238  for (int i = 0; i < w; i++)
239  {
240  z_exp += m_Wm(i) * Z_exp.col(i);
241  }
242 
243  // compute the covariance of measurement
244  Pzz_t Pzz = Pzz_t::Zero(z.rows(), z.rows());
245  for (int i = 0; i < w; i++)
246  {
247  Pzz += m_Wc(i) * (Z_exp.col(i) - z_exp) * (Z_exp.col(i) - z_exp).transpose();
248  }
249  Pzz += R;
250 
251  // compute cross covariance
252  K_t Pxz = K_t::Zero(n, z.rows());
253  for (int i = 0; i < w; i++)
254  {
255  Pxz += m_Wc(i) * (S.col(i) - x) * (Z_exp.col(i) - z_exp).transpose();
256  }
257 
258  // compute Kalman gain
259  const z_t inn = (z - z_exp); // innovation
260  const K_t K = computeKalmanGain(sc.x, inn, Pxz, Pzz);
261 
262  // check whether the inverse produced valid numbers
263  if (!K.array().isFinite().all())
264  {
265  ROS_ERROR("UKF: inverting of Pzz in correction update produced non-finite numbers!!! Fix your covariances (the measurement's is probably too low...)");
266  throw inverse_exception();
267  }
268 
269  // correct
270  statecov_t ret;
271  ret.x = x + K * inn;
272  ret.P = P - K * Pzz * K.transpose();
273  return ret;
274  }
275 
276  //}
277 
278 } // namespace mrs_lib
279 
280 #endif
mrs_lib::UKF::X_t
typename Eigen::Matrix< double, n, w > X_t
State sigma points matrix.
Definition: ukf.h:50
mrs_lib::UKF::setConstants
void setConstants(const double alpha, const double kappa, const double beta)
Changes the Unscented Transform parameters.
Definition: ukf.hpp:61
mrs_lib::UKF::UKF
UKF()
Convenience default constructor.
Definition: ukf.hpp:20
mrs_lib::KalmanFilter::R_t
Eigen::Matrix< double, p, p > R_t
Measurement noise covariance matrix type .
Definition: kalman_filter.h:39
mrs_lib::UKF::K_t
typename Eigen::Matrix< double, n, p > K_t
Kalman gain matrix.
Definition: ukf.h:53
mrs_lib::UKF::transition_model_t
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
mrs_lib::UKF::correct
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
mrs_lib::UKF::statecov_t
typename Base_class::statecov_t statecov_t
typedef of a helper struct for state and covariance
Definition: ukf.h:73
mrs_lib::UKF::W_t
typename Eigen::Matrix< double, w, 1 > W_t
weights vector (2n+1)*1 typedef
Definition: ukf.h:71
mrs_lib::KalmanFilter::u_t
Eigen::Matrix< double, m, 1 > u_t
Input vector type .
Definition: kalman_filter.h:35
mrs_lib::UKF::observation_model_t
typename std::function< z_t(const x_t &)> observation_model_t
function of the observation model typedef
Definition: ukf.h:77
mrs_lib::UKF::setObservationModel
void setObservationModel(const observation_model_t &observation_model)
Changes the observation model function.
Definition: ukf.hpp:85
mrs_lib::KalmanFilter::x_t
Eigen::Matrix< double, n, 1 > x_t
State vector type .
Definition: kalman_filter.h:34
mrs_lib
All mrs_lib functions, classes, variables and definitions are contained in this namespace.
Definition: attitude_converter.h:29
mrs_lib::UKF::setTransitionModel
void setTransitionModel(const transition_model_t &transition_model)
Changes the transition model function.
Definition: ukf.hpp:75
mrs_lib::UKF::Pzz_t
typename Eigen::Matrix< double, p, p > Pzz_t
Pzz helper matrix.
Definition: ukf.h:52
mrs_lib::UKF::P_t
typename Base_class::P_t P_t
state covariance n*n typedef
Definition: ukf.h:65
mrs_lib::UKF::predict
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
mrs_lib::KalmanFilter::Q_t
Eigen::Matrix< double, n, n > Q_t
Process noise covariance matrix type .
Definition: kalman_filter.h:40
mrs_lib::KalmanFilter::z_t
Eigen::Matrix< double, p, 1 > z_t
Measurement vector type .
Definition: kalman_filter.h:36
mrs_lib::UKF::Z_t
typename Eigen::Matrix< double, p, w > Z_t
Measurement sigma points matrix.
Definition: ukf.h:51
ukf.h
Defines UKF - a class implementing the Unscented Kalman Filter .
mrs_lib::UKF
Implementation of the Unscented Kalman filter .
Definition: ukf.h:39
mrs_lib::UKF::inverse_exception
is thrown when taking the inverse of a matrix fails during kalman gain calculation
Definition: ukf.h:89