mrs_lib
Various reusable classes, functions and utilities for use in MRS projects
nckf.h
1 // clang: MatousFormat
6 #ifndef NCKFSYSTEMMODELS_H
7 #define NCKFSYSTEMMODELS_H
8 
9 #include <mrs_lib/lkf.h>
10 #include <mrs_lib/ukf.h>
11 
12 namespace mrs_lib
13 {
14 
15  /* class NCLKF //{ */
27  template <int n_states, int n_inputs, int n_measurements>
28  class NCLKF : public LKF<n_states, n_inputs, n_measurements>
29  {
30  public:
31  /* NCLKF definitions (typedefs, constants etc) //{ */
32  static const int n = n_states;
33  static const int m = n_inputs;
34  static const int p = n_measurements;
37  using x_t = typename Base_class::x_t;
38  using u_t = typename Base_class::u_t;
39  using z_t = typename Base_class::z_t;
40  using P_t = typename Base_class::P_t;
41  using R_t = typename Base_class::R_t;
42  using Q_t = typename Base_class::Q_t;
45  using A_t = typename Base_class::A_t;
46  using B_t = typename Base_class::B_t;
47  using H_t = typename Base_class::H_t;
48  using K_t = typename Base_class::K_t;
49  //}
50 
51  public:
58  NCLKF(){};
59 
68  NCLKF(const A_t& A, const B_t& B, const H_t& H, const double l) : Base_class(A, B, H), l(l) {};
69 
70  protected:
71  double l;
72 
73  protected:
74  /* computeKalmanGain() method //{ */
75  virtual K_t computeKalmanGain(const statecov_t& sc, const z_t& z, const R_t& R, const H_t& H) const override
76  {
77  const R_t W = H * sc.P * H.transpose() + R;
78  const R_t W_inv = Base_class::invert_W(W);
79  const K_t K_orig = sc.P * H.transpose() * W_inv;
80 
81  const z_t inn = z - (H * sc.x); // innovation
82  const x_t x = sc.x + K_orig * inn;
83  const double inn_scale = inn.transpose() * W_inv * inn;
84 
85  const double x_norm = x.norm();
86  const K_t K = K_orig + (l/x_norm - 1.0) * x * (inn.transpose() * W_inv) / inn_scale;
87 
88  return K;
89  }
90  //}
91 
92  };
93  //}
94 
95  /* class NCLKF_partial //{ */
96 
111  template <int n_states, int n_inputs, int n_measurements, int n_norm_constrained_states>
112  class NCLKF_partial : public NCLKF<n_states, n_inputs, n_measurements>
113  {
114  public:
115  /* NCLKF_partial definitions (typedefs, constants etc) //{ */
116  static const int n = n_states;
117  static const int m = n_inputs;
118  static const int p = n_measurements;
119  static const int nq = n_norm_constrained_states;
122  using x_t = typename Base_class::x_t;
123  using u_t = typename Base_class::u_t;
124  using z_t = typename Base_class::z_t;
125  using P_t = typename Base_class::P_t;
126  using R_t = typename Base_class::R_t;
127  using Q_t = typename Base_class::Q_t;
130  using A_t = typename Base_class::A_t;
131  using B_t = typename Base_class::B_t;
132  using H_t = typename Base_class::H_t;
133  using K_t = typename Base_class::K_t;
135  using indices_t = std::array<int, nq>;
136  using xq_t = Eigen::Matrix<double, nq, 1>;
137  using Hq_t = Eigen::Matrix<double, p, nq>;
138  using Kq_t = Eigen::Matrix<double, nq, p>;
139  //}
140 
141  public:
149 
159  NCLKF_partial(const A_t& A, const B_t& B, const H_t& H, const double l, const indices_t& norm_constrained_indices)
160  : Base_class(A, B, H, l),
161  norm_indices(norm_constrained_indices)
162  {
163  if (n >= 0)
164  {
165  for (auto& idx : norm_indices)
166  {
167  if (idx > n)
168  {
169  ROS_ERROR("[NCLKF_partial]: Index of a norm-constrained state (%d) cannot be higher than the number of states (%d)! Setting it to zero.", idx, n);
170  idx = 0;
171  }
172  if (idx < 0)
173  {
174  ROS_ERROR("[NCLKF_partial]: Index of a norm-constrained state (%d) cannot be less than zero! Setting it to zero.", idx);
175  idx = 0;
176  }
177  }
178  }
179  };
180 
181  private:
182  indices_t norm_indices;
183 
184  /* helper methods for Eigen matrix subscripting //{ */
185 
186  template <typename T, int rows, int cols, size_t out_rows, size_t out_cols>
187  Eigen::Matrix<T, out_rows, out_cols> select_indices(const Eigen::Matrix<T, rows, cols>& mat, const std::array<int, out_rows>& row_indices, const std::array<int, out_cols>& col_indices) const
188  {
189  Eigen::Matrix<T, out_rows, cols> tmp;
190  for (int it = 0; it < (int)out_rows; it++)
191  {
192  const auto row = row_indices[it];
193  assert(rows < 0 || row < rows);
194  tmp.row(it) = mat.row(row);
195  }
196 
197  Eigen::Matrix<T, out_rows, out_cols> ret;
198  for (int it = 0; it < (int)out_cols; it++)
199  {
200  const auto col = col_indices[it];
201  assert(cols < 0 || col < cols);
202  tmp.col(it) = tmp.col(col);
203  }
204 
205  return ret;
206  }
207 
208  template <typename T, int rows, int cols, size_t out_rows>
209  Eigen::Matrix<T, out_rows, cols> select_rows(const Eigen::Matrix<T, rows, cols>& mat, const std::array<int, out_rows>& row_indices) const
210  {
211  Eigen::Matrix<T, out_rows, cols> ret;
212  for (int it = 0; it < (int)out_rows; it++)
213  {
214  const auto row = row_indices[it];
215  assert(rows < 0 || row < rows);
216  ret.row(it) = mat.row(row);
217  }
218  return ret;
219  }
220 
221  template <typename T, int rows, int cols, size_t out_cols>
222  Eigen::Matrix<T, rows, out_cols> select_cols(const Eigen::Matrix<T, rows, cols>& mat, const std::array<int, out_cols>& col_indices) const
223  {
224  Eigen::Matrix<T, rows, out_cols> ret;
225  for (int it = 0; it < (int)out_cols; it++)
226  {
227  const auto col = col_indices[it];
228  assert(cols < 0 || col < cols);
229  ret.col(it) = mat.col(col);
230  }
231  return ret;
232  }
233 
234  template <typename T, int rows, size_t out_rows>
235  Eigen::Matrix<T, out_rows, 1> select_indices(const Eigen::Matrix<T, rows, 1>& vec, const std::array<int, out_rows>& row_indices) const
236  {
237  return select_rows(vec, row_indices);
238  }
239 
240  template <typename T, int rows, int cols, size_t n_rows>
241  void set_rows(const Eigen::Matrix<T, (size_t)n_rows, cols>& from_mat, Eigen::Matrix<T, rows, cols>& to_mat, const std::array<int, n_rows>& row_indices) const
242  {
243  for (int rit = 0; rit < (int)n_rows; rit++)
244  {
245  const auto ridx = row_indices.at(rit);
246  assert(rows < 0 || ridx < rows);
247  to_mat.row(ridx) = from_mat.row(rit);
248  }
249  }
250 
251  //}
252 
253  protected:
254  /* computeKalmanGain() method //{ */
255  virtual K_t computeKalmanGain(const statecov_t& sc, const z_t& z, const R_t& R, const H_t& H) const override
256  {
257  const R_t W = H * sc.P * H.transpose() + R;
258  const R_t W_inv = Base_class::invert_W(W);
259  K_t K = sc.P * H.transpose() * W_inv;
260 
261  // calculate the kalman gain for the norm-constrained states
262  {
263  const Kq_t K_orig = select_rows(K, norm_indices);
264  const xq_t xq = select_indices(sc.x, norm_indices);
265  const Hq_t Hq = select_cols(H, norm_indices);
266  const z_t inn = z - (Hq * xq); // innovation
267  const xq_t x = xq + K_orig * inn;
268  const double inn_scale = inn.transpose() * W_inv * inn;
269 
270  const double x_norm = x.norm();
271  const Kq_t Kq = K_orig + (Base_class::l/x_norm - 1.0) * x * (inn.transpose() * W_inv) / inn_scale;
272  set_rows(Kq, K, norm_indices);
273  }
274 
275  return K;
276  }
277  //}
278 
279  };
280  //}
281 
282  /* class NCUKF //{ */
283  template <int n_states, int n_inputs, int n_measurements>
284  class NCUKF : public UKF<n_states, n_inputs, n_measurements>
285  {
286  public:
287  /* LKF definitions (typedefs, constants etc) //{ */
288  static const int n = n_states;
289  static const int m = n_inputs;
290  static const int p = n_measurements;
291  using Base_class = UKF<n, m, p>;
292 
293  using x_t = typename Base_class::x_t; // state vector n*1
294  using u_t = typename Base_class::u_t; // input vector m*1
295  using z_t = typename Base_class::z_t; // measurement vector p*1
296  using P_t = typename Base_class::P_t; // state covariance n*n
297  using R_t = typename Base_class::R_t; // measurement covariance p*p
298  using Q_t = typename Base_class::Q_t; // measurement covariance p*p
299  using statecov_t = typename Base_class::statecov_t; // helper struct for state and covariance
300 
301  using transition_model_t = typename Base_class::transition_model_t;
302  using observation_model_t = typename Base_class::observation_model_t;
303 
304  using X_t = typename Base_class::X_t; // state sigma points matrix n*w
305  using Z_t = typename Base_class::Z_t; // measurement sigma points matrix p*w
306  using Pzz_t = typename Base_class::Pzz_t; // Pzz helper matrix p*n
307  using K_t = typename Base_class::K_t; // kalman gain n*p
308  //}
309 
310  public:
311  NCUKF(){};
312 
313  NCUKF(const transition_model_t& transition_model, const observation_model_t& observation_model, const double l, const double alpha = 1e-3, const double kappa = 1, const double beta = 2)
314  : Base_class(transition_model, observation_model, alpha, kappa, beta), l(l) {};
315 
316  public:
317  /* correct() method //{ */
318  virtual statecov_t correct(const statecov_t& sc, const z_t& z, const R_t& R) const override
319  {
320  const auto& x = sc.x;
321  const auto& P = sc.P;
322  const X_t S = Base_class::computeSigmas(x, P);
323 
324  // propagate sigmas through the observation model
325  Z_t Z_exp;
326  for (int i = 0; i < Base_class::w; i++)
327  {
328  Z_exp.col(i) = Base_class::m_observation_model(S.col(i));
329  }
330 
331  // compute expected measurement
332  z_t z_exp = z_t::Zero();
333  for (int i = 0; i < Base_class::w; i++)
334  {
335  z_exp += Base_class::m_Wm(i) * Z_exp.col(i);
336  }
337 
338  // compute the covariance of measurement
339  Pzz_t Pzz = Pzz_t::Zero();
340  for (int i = 0; i < Base_class::w; i++)
341  {
342  Pzz += Base_class::m_Wc(i) * (Z_exp.col(i) - z_exp) * (Z_exp.col(i) - z_exp).transpose();
343  }
344  Pzz += R;
345 
346  // compute cross covariance
347  K_t Pxz = K_t::Zero();
348  for (int i = 0; i < Base_class::w; i++)
349  {
350  Pxz += Base_class::m_Wc(i) * (S.col(i) - x) * (Z_exp.col(i) - z_exp).transpose();
351  }
352 
353  // compute Kalman gain
354  const z_t inn = (z - z_exp); // innovation
355  const K_t K = computeKalmanGain(sc.x, inn, Pxz, Pzz);
356 
357  // check whether the inverse produced valid numbers
358  if (!K.array().isFinite().all())
359  {
360  ROS_ERROR("UKF: inverting of Pzz in correction update produced non-finite numbers!!! Fix your covariances (the measurement's is probably too low...)");
361  throw typename Base_class::inverse_exception();
362  }
363 
364  // correct
365  statecov_t ret;
366  ret.x = x + K * inn;
367  ret.P = P - K * Pxz.transpose() - Pxz * K.transpose() + K * Pzz * K.transpose();
368  return ret;
369  }
370  //}
371 
372  protected:
373  double l;
374 
375  protected:
376  /* computeKalmanGain() method //{ */
377  virtual K_t computeKalmanGain(const x_t& x, const z_t& inn, const K_t& Pxz, const Pzz_t& Pzz) const override
378  {
379  const Pzz_t Pzz_inv = Base_class::computeInverse(Pzz);
380  const K_t K_orig = Pxz * Pzz_inv;
381 
382  const x_t x_pred = x + K_orig * inn;
383  const double inn_scale = inn.transpose() * Pzz_inv * inn;
384 
385  const double x_norm = x_pred.norm();
386  const K_t K = K_orig + (l/x_norm - 1.0) * x_pred * (inn.transpose() * Pzz_inv) / inn_scale;
387 
388  return K;
389  }
390  //}
391 
392  };
393  //}
394 
395 } // namespace mrs_lib
396 
397 #endif // NCKFSYSTEMMODELS_H
398 
mrs_lib::UKF::X_t
typename Eigen::Matrix< double, n, w > X_t
State sigma points matrix.
Definition: ukf.h:50
mrs_lib::UKF::Q_t
typename Base_class::Q_t Q_t
process covariance n*n typedef
Definition: ukf.h:69
mrs_lib::LKF::Q_t
typename Base_class::Q_t Q_t
Process noise covariance matrix type .
Definition: lkf.h:52
mrs_lib::NCLKF::B_t
typename Base_class::B_t B_t
Input to state mapping matrix type .
Definition: nckf.h:46
mrs_lib::LKF::H_t
Eigen::Matrix< double, p, n > H_t
State to measurement mapping matrix type .
Definition: lkf.h:57
mrs_lib::LKF::P_t
typename Base_class::P_t P_t
State uncertainty covariance matrix type .
Definition: lkf.h:50
mrs_lib::NCLKF::m
static const int m
Length of the input vector of the system.
Definition: nckf.h:33
mrs_lib::NCLKF::u_t
typename Base_class::u_t u_t
Input vector type .
Definition: nckf.h:38
mrs_lib::LKF::H
H_t H
The state to measurement mapping matrix .
Definition: lkf.h:150
mrs_lib::NCLKF_partial::m
static const int m
Length of the input vector of the system.
Definition: nckf.h:117
mrs_lib::UKF::z_t
typename Base_class::z_t z_t
measurement vector p*1 typedef
Definition: ukf.h:63
mrs_lib::NCLKF::p
static const int p
Length of the measurement vector of the system.
Definition: nckf.h:34
mrs_lib::NCLKF_partial::NCLKF_partial
NCLKF_partial(const A_t &A, const B_t &B, const H_t &H, const double l, const indices_t &norm_constrained_indices)
The main constructor.
Definition: nckf.h:159
mrs_lib::LKF::statecov_t
typename Base_class::statecov_t statecov_t
Helper struct for passing around the state and its covariance in one variable.
Definition: lkf.h:53
mrs_lib::NCUKF::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: nckf.h:318
mrs_lib::LKF::R_t
typename Base_class::R_t R_t
Measurement noise covariance matrix type .
Definition: lkf.h:51
mrs_lib::NCLKF::NCLKF
NCLKF()
Convenience default constructor.
Definition: nckf.h:58
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::NCLKF::z_t
typename Base_class::z_t z_t
Measurement vector type .
Definition: nckf.h:39
mrs_lib::NCLKF::Q_t
typename Base_class::Q_t Q_t
Process noise covariance matrix type .
Definition: nckf.h:42
mrs_lib::NCLKF::x_t
typename Base_class::x_t x_t
State vector type .
Definition: nckf.h:37
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::NCLKF_partial::n
static const int n
Length of the state vector of the system.
Definition: nckf.h:116
mrs_lib::NCLKF_partial::Hq_t
Eigen::Matrix< double, p, nq > Hq_t
Norm-constrained measurement mapping type .
Definition: nckf.h:137
mrs_lib::KalmanFilter::P_t
Eigen::Matrix< double, n, n > P_t
State uncertainty covariance matrix type .
Definition: kalman_filter.h:38
mrs_lib::NCLKF_partial
This class implements the partially norm-constrained linear Kalman filter .
Definition: nckf.h:112
mrs_lib::UKF::u_t
typename Base_class::u_t u_t
input vector m*1 typedef
Definition: ukf.h:61
mrs_lib::NCLKF::P_t
typename Base_class::P_t P_t
State uncertainty covariance matrix type .
Definition: nckf.h:40
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::NCLKF_partial::nq
static const int nq
Number of states to which the norm constraint applies.
Definition: nckf.h:119
mrs_lib::NCLKF::A_t
typename Base_class::A_t A_t
System transition matrix type .
Definition: nckf.h:45
mrs_lib::KalmanFilter::u_t
Eigen::Matrix< double, m, 1 > u_t
Input vector type .
Definition: kalman_filter.h:35
mrs_lib::UKF::w
static constexpr int w
Number of sigma points/weights.
Definition: ukf.h:46
mrs_lib::LKF::B_t
Eigen::Matrix< double, n, m > B_t
Input to state mapping matrix type .
Definition: lkf.h:56
mrs_lib::LKF::x_t
typename Base_class::x_t x_t
State vector type .
Definition: lkf.h:47
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::KalmanFilter::x_t
Eigen::Matrix< double, n, 1 > x_t
State vector type .
Definition: kalman_filter.h:34
mrs_lib::NCLKF_partial::xq_t
Eigen::Matrix< double, nq, 1 > xq_t
Norm-constrained states vector type .
Definition: nckf.h:136
mrs_lib::NCLKF_partial::NCLKF_partial
NCLKF_partial()
Convenience default constructor.
Definition: nckf.h:148
mrs_lib
All mrs_lib functions, classes, variables and definitions are contained in this namespace.
Definition: attitude_converter.h:29
mrs_lib::UKF::R_t
typename Base_class::R_t R_t
measurement covariance p*p typedef
Definition: ukf.h:67
mrs_lib::NCLKF::R_t
typename Base_class::R_t R_t
Measurement noise covariance matrix type .
Definition: nckf.h:41
mrs_lib::LKF::A_t
Eigen::Matrix< double, n, n > A_t
System transition matrix type .
Definition: lkf.h:55
mrs_lib::LKF::B
B_t B
The input to state mapping matrix .
Definition: lkf.h:149
mrs_lib::UKF::Pzz_t
typename Eigen::Matrix< double, p, p > Pzz_t
Pzz helper matrix.
Definition: ukf.h:52
mrs_lib::NCLKF_partial::indices_t
std::array< int, nq > indices_t
Indices of the norm-constrained states type.
Definition: nckf.h:135
mrs_lib::UKF::P_t
typename Base_class::P_t P_t
state covariance n*n typedef
Definition: ukf.h:65
mrs_lib::NCLKF::NCLKF
NCLKF(const A_t &A, const B_t &B, const H_t &H, const double l)
The main constructor.
Definition: nckf.h:68
mrs_lib::NCLKF::n
static const int n
Length of the state vector of the system.
Definition: nckf.h:32
mrs_lib::LKF
Implementation of the Linear Kalman filter .
Definition: lkf.h:38
mrs_lib::LKF::K_t
Eigen::Matrix< double, n, p > K_t
Kalman gain matrix type .
Definition: lkf.h:58
mrs_lib::NCLKF::K_t
typename Base_class::K_t K_t
Kalman gain matrix type .
Definition: nckf.h:48
mrs_lib::NCLKF::statecov_t
typename Base_class::statecov_t statecov_t
Helper struct for passing around the state and its covariance in one variable.
Definition: nckf.h:43
mrs_lib::NCLKF_partial::p
static const int p
Length of the measurement vector of the system.
Definition: nckf.h:118
mrs_lib::NCLKF::H_t
typename Base_class::H_t H_t
State to measurement mapping matrix type .
Definition: nckf.h:47
mrs_lib::UKF::x_t
typename Base_class::x_t x_t
state vector n*1 typedef
Definition: ukf.h:59
mrs_lib::KalmanFilter::Q_t
Eigen::Matrix< double, n, n > Q_t
Process noise covariance matrix type .
Definition: kalman_filter.h:40
mrs_lib::LKF::z_t
typename Base_class::z_t z_t
Measurement vector type .
Definition: lkf.h:49
mrs_lib::KalmanFilter::z_t
Eigen::Matrix< double, p, 1 > z_t
Measurement vector type .
Definition: kalman_filter.h:36
mrs_lib::NCUKF
Definition: nckf.h:284
mrs_lib::NCLKF_partial::K_t
typename Base_class::K_t K_t
Kalman gain matrix type .
Definition: nckf.h:133
mrs_lib::UKF::Z_t
typename Eigen::Matrix< double, p, w > Z_t
Measurement sigma points matrix.
Definition: ukf.h:51
mrs_lib::LKF::A
A_t A
The system transition matrix .
Definition: lkf.h:144
ukf.h
Defines UKF - a class implementing the Unscented Kalman Filter .
mrs_lib::LKF::u_t
typename Base_class::u_t u_t
Input vector type .
Definition: lkf.h:48
mrs_lib::NCLKF
This class implements the norm-constrained linear Kalman filter .
Definition: nckf.h:28
mrs_lib::NCLKF_partial::Kq_t
Eigen::Matrix< double, nq, p > Kq_t
Norm-constrained kalman gain type .
Definition: nckf.h:138
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
mrs_lib::KalmanFilter
This abstract class defines common interfaces and types for a generic Kalman filter.
Definition: kalman_filter.h:26