TransformBtwRobotsUnaryFactorEM.h
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4  * Atlanta, Georgia 30332-0415
5  * All Rights Reserved
6  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8  * See LICENSE for the license information
9 
10  * -------------------------------------------------------------------------- */
11 
17 #pragma once
18 
24 #include <gtsam/base/Testable.h>
25 #include <gtsam/base/Lie.h>
26 
27 #include <ostream>
28 
29 namespace gtsam {
30 
36  template<class VALUE>
38 
39  public:
40 
41  typedef VALUE T;
42 
43  private:
44 
47 
49 
50  VALUE measured_;
52  Values valA_; // given values for robot A map\trajectory
53  Values valB_; // given values for robot B map\trajectory
54  Key keyA_; // key of robot A to which the measurement refers
55  Key keyB_; // key of robot B to which the measurement refers
56 
57  // TODO: create function to update valA_ and valB_
58 
61 
62  double prior_inlier_;
64 
66  mutable bool start_with_M_step_;
67 
71 
72  public:
73 
74  // shorthand for a smart pointer to a factor
76 
79 
82  const Values& valA, const Values& valB,
83  const SharedGaussian& model_inlier, const SharedGaussian& model_outlier,
84  const double prior_inlier, const double prior_outlier,
85  const bool flag_bump_up_near_zero_probs = false,
86  const bool start_with_M_step = false) :
87  Base(cref_list_of<1>(key)), key_(key), measured_(measured), keyA_(keyA), keyB_(keyB),
88  model_inlier_(model_inlier), model_outlier_(model_outlier),
89  prior_inlier_(prior_inlier), prior_outlier_(prior_outlier), flag_bump_up_near_zero_probs_(flag_bump_up_near_zero_probs),
90  start_with_M_step_(false){
91 
92  setValAValB(valA, valB);
93 
94  }
95 
97 
98 
100  NonlinearFactor::shared_ptr clone() const override { return boost::make_shared<This>(*this); }
101 
102 
106  void print(const std::string& s, const KeyFormatter& keyFormatter = DefaultKeyFormatter) const override {
107  std::cout << s << "TransformBtwRobotsUnaryFactorEM("
108  << keyFormatter(key_) << ")\n";
109  std::cout << "MR between factor keys: "
110  << keyFormatter(keyA_) << ","
111  << keyFormatter(keyB_) << "\n";
112  measured_.print(" measured: ");
113  model_inlier_->print(" noise model inlier: ");
114  model_outlier_->print(" noise model outlier: ");
115  std::cout << "(prior_inlier, prior_outlier_) = ("
116  << prior_inlier_ << ","
117  << prior_outlier_ << ")\n";
118  // Base::print(s, keyFormatter);
119  }
120 
122  bool equals(const NonlinearFactor& f, double tol=1e-9) const override {
123  const This *t = dynamic_cast<const This*> (&f);
124 
125  if(t && Base::equals(f))
126  return key_ == t->key_ && measured_.equals(t->measured_) &&
127  // model_inlier_->equals(t->model_inlier_ ) && // TODO: fix here
128  // model_outlier_->equals(t->model_outlier_ ) &&
129  prior_outlier_ == t->prior_outlier_ && prior_inlier_ == t->prior_inlier_;
130  else
131  return false;
132  }
133 
136  /* ************************************************************************* */
137  void setValAValB(const Values& valA, const Values& valB){
138  if ( (!valA.exists(keyA_)) && (!valB.exists(keyA_)) && (!valA.exists(keyB_)) && (!valB.exists(keyB_)) )
139  throw("something is wrong!");
140 
141  // TODO: make sure the two keys belong to different robots
142 
143  if (valA.exists(keyA_)){
144  valA_ = valA;
145  valB_ = valB;
146  }
147  else {
148  valA_ = valB;
149  valB_ = valA;
150  }
151  }
152 
153  /* ************************************************************************* */
154  double error(const Values& x) const override {
155  return whitenedError(x).squaredNorm();
156  }
157 
158  /* ************************************************************************* */
164  /* This version of linearize recalculates the noise model each time */
165  boost::shared_ptr<GaussianFactor> linearize(const Values& x) const override {
166  // Only linearize if the factor is active
167  if (!this->active(x))
168  return boost::shared_ptr<JacobianFactor>();
169 
170  //std::cout<<"About to linearize"<<std::endl;
171  Matrix A1;
172  std::vector<Matrix> A(this->size());
173  Vector b = -whitenedError(x, A);
174  A1 = A[0];
175 
177  new JacobianFactor(key_, A1, b, noiseModel::Unit::Create(b.size())));
178  }
179 
180 
181  /* ************************************************************************* */
183  boost::optional<std::vector<Matrix>&> H = boost::none) const {
184 
185  bool debug = true;
186 
187  Matrix H_compose, H_between1, H_dummy;
188 
189  T orgA_T_currA = valA_.at<T>(keyA_);
190  T orgB_T_currB = valB_.at<T>(keyB_);
191 
192  T orgA_T_orgB = x.at<T>(key_);
193 
194  T orgA_T_currB = orgA_T_orgB.compose(orgB_T_currB, H_compose, H_dummy);
195 
196  T currA_T_currB_pred = orgA_T_currA.between(orgA_T_currB, H_dummy, H_between1);
197 
198  T currA_T_currB_msr = measured_;
199 
200  Vector err = currA_T_currB_msr.localCoordinates(currA_T_currB_pred);
201 
202  // Calculate indicator probabilities (inlier and outlier)
203  Vector p_inlier_outlier = calcIndicatorProb(x, err);
204  double p_inlier = p_inlier_outlier[0];
205  double p_outlier = p_inlier_outlier[1];
206 
207  if (start_with_M_step_){
208  start_with_M_step_ = false;
209 
210  p_inlier = 0.5;
211  p_outlier = 0.5;
212  }
213 
214  Vector err_wh_inlier = model_inlier_->whiten(err);
215  Vector err_wh_outlier = model_outlier_->whiten(err);
216 
217  Matrix invCov_inlier = model_inlier_->R().transpose() * model_inlier_->R();
218  Matrix invCov_outlier = model_outlier_->R().transpose() * model_outlier_->R();
219 
220  Vector err_wh_eq;
221  err_wh_eq.resize(err_wh_inlier.rows()*2);
222  err_wh_eq << sqrt(p_inlier) * err_wh_inlier.array() , sqrt(p_outlier) * err_wh_outlier.array();
223 
224  Matrix H_unwh = H_compose * H_between1;
225 
226  if (H){
227 
228  Matrix H_inlier = sqrt(p_inlier)*model_inlier_->Whiten(H_unwh);
229  Matrix H_outlier = sqrt(p_outlier)*model_outlier_->Whiten(H_unwh);
230  Matrix H_aug = stack(2, &H_inlier, &H_outlier);
231 
232  (*H)[0].resize(H_aug.rows(),H_aug.cols());
233  (*H)[0] = H_aug;
234  }
235 
236  if (debug){
237  // std::cout<<"H_compose - rows, cols, : "<<H_compose.rows()<<", "<< H_compose.cols()<<std::endl;
238  // std::cout<<"H_between1 - rows, cols, : "<<H_between1.rows()<<", "<< H_between1.cols()<<std::endl;
239  // std::cout<<"H_unwh - rows, cols, : "<<H_unwh.rows()<<", "<< H_unwh.cols()<<std::endl;
240  // std::cout<<"H_unwh: "<<std:endl<<H_unwh[0]
241 
242  }
243 
244 
245  return err_wh_eq;
246  }
247 
248  /* ************************************************************************* */
250 
251  Vector err = unwhitenedError(x);
252 
253  return this->calcIndicatorProb(x, err);
254  }
255 
256  /* ************************************************************************* */
257  Vector calcIndicatorProb(const Values& x, const Vector& err) const {
258 
259  // Calculate indicator probabilities (inlier and outlier)
260  Vector err_wh_inlier = model_inlier_->whiten(err);
261  Vector err_wh_outlier = model_outlier_->whiten(err);
262 
263  Matrix invCov_inlier = model_inlier_->R().transpose() * model_inlier_->R();
264  Matrix invCov_outlier = model_outlier_->R().transpose() * model_outlier_->R();
265 
266  double p_inlier = prior_inlier_ * sqrt(invCov_inlier.norm()) * exp( -0.5 * err_wh_inlier.dot(err_wh_inlier) );
267  double p_outlier = prior_outlier_ * sqrt(invCov_outlier.norm()) * exp( -0.5 * err_wh_outlier.dot(err_wh_outlier) );
268 
269  double sumP = p_inlier + p_outlier;
270  p_inlier /= sumP;
271  p_outlier /= sumP;
272 
273  if (flag_bump_up_near_zero_probs_){
274  // Bump up near-zero probabilities (as in linerFlow.h)
275  double minP = 0.05; // == 0.1 / 2 indicator variables
276  if (p_inlier < minP || p_outlier < minP){
277  if (p_inlier < minP)
278  p_inlier = minP;
279  if (p_outlier < minP)
280  p_outlier = minP;
281  sumP = p_inlier + p_outlier;
282  p_inlier /= sumP;
283  p_outlier /= sumP;
284  }
285  }
286 
287  return (Vector(2) << p_inlier, p_outlier).finished();
288  }
289 
290  /* ************************************************************************* */
291  Vector unwhitenedError(const Values& x) const {
292 
293  T orgA_T_currA = valA_.at<T>(keyA_);
294  T orgB_T_currB = valB_.at<T>(keyB_);
295 
296  T orgA_T_orgB = x.at<T>(key_);
297 
298  T orgA_T_currB = orgA_T_orgB.compose(orgB_T_currB);
299 
300  T currA_T_currB_pred = orgA_T_currA.between(orgA_T_currB);
301 
302  T currA_T_currB_msr = measured_;
303 
304  return currA_T_currB_msr.localCoordinates(currA_T_currB_pred);
305  }
306 
307  /* ************************************************************************* */
309  return model_inlier_;
310  }
311 
312  /* ************************************************************************* */
314  return model_outlier_;
315  }
316 
317  /* ************************************************************************* */
319  return (model_inlier_->R().transpose()*model_inlier_->R()).inverse();
320  }
321 
322  /* ************************************************************************* */
324  return (model_outlier_->R().transpose()*model_outlier_->R()).inverse();
325  }
326 
327  /* ************************************************************************* */
329  /* given marginals version, don't need to marginal multiple times if update a lot */
330 
331  KeyVector Keys;
332  Keys.push_back(keyA_);
333  Keys.push_back(keyB_);
334  JointMarginal joint_marginal12 = marginals.jointMarginalCovariance(Keys);
335  Matrix cov1 = joint_marginal12(keyA_, keyA_);
336  Matrix cov2 = joint_marginal12(keyB_, keyB_);
337  Matrix cov12 = joint_marginal12(keyA_, keyB_);
338 
339  updateNoiseModels_givenCovs(values, cov1, cov2, cov12);
340  }
341 
342  /* ************************************************************************* */
344  /* Update model_inlier_ and model_outlier_ to account for uncertainty in robot trajectories
345  * (note these are given in the E step, where indicator probabilities are calculated).
346  *
347  * Principle: R += [H1 H2] * joint_cov12 * [H1 H2]', where H1, H2 are Jacobians of the
348  * unwhitened error w.r.t. states, and R is the measurement covariance (inlier or outlier modes).
349  *
350  * TODO: improve efficiency (info form)
351  */
352 
353  // get joint covariance of the involved states
354 
355  Marginals marginals(graph, values, Marginals::QR);
356 
357  this->updateNoiseModels(values, marginals);
358  }
359 
360  /* ************************************************************************* */
361  void updateNoiseModels_givenCovs(const Values& values, const Matrix& cov1, const Matrix& cov2, const Matrix& cov12){
362  /* Update model_inlier_ and model_outlier_ to account for uncertainty in robot trajectories
363  * (note these are given in the E step, where indicator probabilities are calculated).
364  *
365  * Principle: R += [H1 H2] * joint_cov12 * [H1 H2]', where H1, H2 are Jacobians of the
366  * unwhitened error w.r.t. states, and R is the measurement covariance (inlier or outlier modes).
367  *
368  * TODO: improve efficiency (info form)
369  */
370 
371  const T& p1 = values.at<T>(keyA_);
372  const T& p2 = values.at<T>(keyB_);
373 
374  Matrix H1, H2;
375  p1.between(p2, H1, H2); // h(x)
376 
377  Matrix H;
378  H.resize(H1.rows(), H1.rows()+H2.rows());
379  H << H1, H2; // H = [H1 H2]
380 
381  Matrix joint_cov;
382  joint_cov.resize(cov1.rows()+cov2.rows(), cov1.cols()+cov2.cols());
383  joint_cov << cov1, cov12,
384  cov12.transpose(), cov2;
385 
386  Matrix cov_state = H*joint_cov*H.transpose();
387 
388  // model_inlier_->print("before:");
389 
390  // update inlier and outlier noise models
391  Matrix covRinlier = (model_inlier_->R().transpose()*model_inlier_->R()).inverse();
392  model_inlier_ = noiseModel::Gaussian::Covariance(covRinlier + cov_state);
393 
394  Matrix covRoutlier = (model_outlier_->R().transpose()*model_outlier_->R()).inverse();
395  model_outlier_ = noiseModel::Gaussian::Covariance(covRoutlier + cov_state);
396 
397  // model_inlier_->print("after:");
398  // std::cout<<"covRinlier + cov_state: "<<covRinlier + cov_state<<std::endl;
399  }
400 
401 
402  /* ************************************************************************* */
403 
404  size_t dim() const override {
405  return model_inlier_->R().rows() + model_inlier_->R().cols();
406  }
407 
408  private:
409 
412  template<class ARCHIVE>
413  void serialize(ARCHIVE & ar, const unsigned int /*version*/) {
414  ar & boost::serialization::make_nvp("NonlinearFactor",
415  boost::serialization::base_object<Base>(*this));
416  //ar & BOOST_SERIALIZATION_NVP(measured_);
417  }
418  }; // \class TransformBtwRobotsUnaryFactorEM
419 
421  template<class VALUE>
423  public Testable<TransformBtwRobotsUnaryFactorEM<VALUE> > {
424  };
425 
426 }
bool equals(const This &other, double tol=1e-9) const
check equality
Definition: Factor.cpp:42
TransformBtwRobotsUnaryFactorEM< VALUE > This
bool exists(Key j) const
Definition: Values.cpp:104
static shared_ptr Covariance(const Matrix &covariance, bool smart=true)
Definition: NoiseModel.cpp:116
Concept check for values that can be used in unit tests.
Vector3f p1
Factor Graph consisting of non-linear factors.
EIGEN_DEVICE_FUNC const ExpReturnType exp() const
double error(const Values &x) const override
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:43
A factor with a quadratic error function - a Gaussian.
leaf::MyValues values
Vector whitenedError(const Values &x, boost::optional< std::vector< Matrix > & > H=boost::none) const
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
JointMarginal jointMarginalCovariance(const KeyVector &variables) const
Definition: Marginals.cpp:134
void print(const std::string &s, const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override
boost::shared_ptr< GaussianFactor > linearize(const Values &x) const override
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate set rrange[*:*] noreverse nowriteback set trange[*:*] noreverse nowriteback set urange[*:*] noreverse nowriteback set vrange[*:*] noreverse nowriteback set xlabel matrix size set x2label set timefmt d m y n H
#define GTSAM_CONCEPT_TESTABLE_TYPE(T)
Definition: Testable.h:175
NonlinearFactorGraph graph
static const KeyFormatter DefaultKeyFormatter
Definition: Key.h:43
virtual bool active(const Values &) const
static shared_ptr Create(size_t dim)
Definition: NoiseModel.h:608
void updateNoiseModels(const Values &values, const Marginals &marginals)
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:86
Eigen::VectorXd Vector
Definition: Vector.h:38
const ValueType at(Key j) const
Definition: Values-inl.h:342
NonlinearFactor::shared_ptr clone() const override
Vector calcIndicatorProb(const Values &x, const Vector &err) const
void updateNoiseModels_givenCovs(const Values &values, const Matrix &cov1, const Matrix &cov2, const Matrix &cov12)
boost::shared_ptr< This > shared_ptr
std::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:35
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
static bool debug
Array< double, 1, 3 > e(1./3., 0.5, 2.)
RealScalar s
bool equals(const NonlinearFactor &f, double tol=1e-9) const override
virtual bool equals(const NonlinearFactor &f, double tol=1e-9) const
const G & b
Definition: Group.h:83
boost::shared_ptr< This > shared_ptr
shared_ptr to this class
Base class and basic functions for Lie types.
traits
Definition: chartTesting.h:28
#define GTSAM_CONCEPT_LIE_TYPE(T)
Definition: Lie.h:356
void setValAValB(const Values &valA, const Values &valB)
boost::shared_ptr< Factor > shared_ptr
A shared_ptr to this class.
Definition: Factor.h:60
void serialize(ARCHIVE &ar, const unsigned int)
Non-linear factor base classes.
size_t size() const
Definition: Factor.h:135
Matrix stack(size_t nrMatrices,...)
Definition: Matrix.cpp:396
static Point3 p2
A class for computing marginals in a NonlinearFactorGraph.
const G double tol
Definition: Group.h:83
TransformBtwRobotsUnaryFactorEM(Key key, const VALUE &measured, Key keyA, Key keyB, const Values &valA, const Values &valB, const SharedGaussian &model_inlier, const SharedGaussian &model_outlier, const double prior_inlier, const double prior_outlier, const bool flag_bump_up_near_zero_probs=false, const bool start_with_M_step=false)
Point3 measured
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
EIGEN_DEVICE_FUNC const InverseReturnType inverse() const
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:61
Marginals marginals(graph, result)
noiseModel::Gaussian::shared_ptr SharedGaussian
Definition: NoiseModel.h:735
Point2 t(10, 10)
void updateNoiseModels(const Values &values, const NonlinearFactorGraph &graph)


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:51:06