MRPT  1.9.9
CGridMapAligner.cpp
Go to the documentation of this file.
1 /* +------------------------------------------------------------------------+
2  | Mobile Robot Programming Toolkit (MRPT) |
3  | https://www.mrpt.org/ |
4  | |
5  | Copyright (c) 2005-2019, Individual contributors, see AUTHORS file |
6  | See: https://www.mrpt.org/Authors - All rights reserved. |
7  | Released under BSD License. See: https://www.mrpt.org/License |
8  +------------------------------------------------------------------------+ */
9 
10 #include "slam-precomp.h" // Precompiled headers
11 
17 #include <mrpt/math/TPose2D.h>
19 #include <mrpt/math/geometry.h>
24 #include <mrpt/random.h>
26 #include <mrpt/slam/CICP.h>
27 #include <mrpt/system/filesystem.h>
28 #include <mrpt/tfest/se2.h>
29 #include <Eigen/Dense>
30 
31 using namespace mrpt::math;
32 using namespace mrpt::slam;
33 using namespace mrpt::maps;
34 using namespace mrpt::poses;
35 using namespace mrpt::random;
36 using namespace mrpt::img;
37 using namespace mrpt::system;
38 using namespace mrpt::vision;
39 using namespace std;
40 
41 CPosePDF::Ptr CGridMapAligner::AlignPDF(
42  const mrpt::maps::CMetricMap* mm1, const mrpt::maps::CMetricMap* mm2,
43  const CPosePDFGaussian& initialEstimationPDF, float* runningTime,
44  void* info)
45 {
47 
48  switch (options.methodSelection)
49  {
50  case CGridMapAligner::amCorrelation:
51  return AlignPDF_correlation(
52  mm1, mm2, initialEstimationPDF, runningTime, info);
53 
54  case CGridMapAligner::amModifiedRANSAC:
55  case CGridMapAligner::amRobustMatch:
56  // The same function has an internal switch for the specific method:
57  return AlignPDF_robustMatch(
58  mm1, mm2, initialEstimationPDF, runningTime, info);
59 
60  default:
61  THROW_EXCEPTION("Wrong value found in 'options.methodSelection'!!");
62  }
63 
64  MRPT_END
65 }
66 
67 bool myVectorOrder(const pair<size_t, float>& o1, const pair<size_t, float>& o2)
68 {
69  return o1.second < o2.second;
70 }
71 
72 /*---------------------------------------------------------------
73  AlignPDF_robustMatch
74 ---------------------------------------------------------------*/
75 CPosePDF::Ptr CGridMapAligner::AlignPDF_robustMatch(
76  const mrpt::maps::CMetricMap* mm1, const mrpt::maps::CMetricMap* mm2,
77  const CPosePDFGaussian& initialEstimationPDF, float* runningTime,
78  void* info)
79 {
81 
82  ASSERT_(
83  options.methodSelection == CGridMapAligner::amRobustMatch ||
84  options.methodSelection == CGridMapAligner::amModifiedRANSAC);
85 
86  TReturnInfo outInfo;
87  mrpt::tfest::TMatchingPairList& correspondences =
88  outInfo.correspondences; // Use directly this placeholder to save 1
89  // variable & 1 copy.
90  mrpt::tfest::TMatchingPairList largestConsensusCorrs;
91 
92  CTicTac* tictac = nullptr;
93 
94  CPose2D grossEst = initialEstimationPDF.mean;
95 
98 
99  std::vector<size_t> idxs1, idxs2;
100  std::vector<size_t>::iterator it1, it2;
101 
102  // Asserts:
103  // -----------------
104  const CMultiMetricMap* multimap1 = nullptr;
105  const CMultiMetricMap* multimap2 = nullptr;
106  const COccupancyGridMap2D* m1 = nullptr;
107  const COccupancyGridMap2D* m2 = nullptr;
108 
109  if (IS_CLASS(*mm1, CMultiMetricMap) && IS_CLASS(*mm2, CMultiMetricMap))
110  {
111  multimap1 = static_cast<const CMultiMetricMap*>(mm1);
112  multimap2 = static_cast<const CMultiMetricMap*>(mm2);
113 
116 
117  m1 = multimap1->mapByClass<COccupancyGridMap2D>().get();
118  m2 = multimap2->mapByClass<COccupancyGridMap2D>().get();
119  }
120  else if (
121  IS_CLASS(*mm1, COccupancyGridMap2D) &&
123  {
124  m1 = static_cast<const COccupancyGridMap2D*>(mm1);
125  m2 = static_cast<const COccupancyGridMap2D*>(mm2);
126  }
127  else
129  "Metric maps must be of classes COccupancyGridMap2D or "
130  "CMultiMetricMap");
131 
132  ASSERT_(m1);
133  ASSERT_(m2);
134 
135  ASSERTMSG_(
136  m1->getResolution() == m2->getResolution(),
137  mrpt::format(
138  "Different resolutions for m1, m2:\n"
139  "\tres(m1) = %f\n\tres(m2) = %f\n",
140  m1->getResolution(), m2->getResolution()));
141 
142  // The algorithm output auxiliar info:
143  // -------------------------------------------------
144  outInfo.goodness = 1.0f;
145 
146  outInfo.landmarks_map1 = lm1;
147  outInfo.landmarks_map2 = lm2;
148 
149  // The PDF to estimate:
150  // ------------------------------------------------------
151  CPosePDFSOG::Ptr pdf_SOG = std::make_shared<CPosePDFSOG>();
152 
153  // Extract features from grid-maps:
154  // ------------------------------------------------------
155  const size_t N1 =
156  std::max(40, mrpt::round(m1->getArea() * options.featsPerSquareMeter));
157  const size_t N2 =
158  std::max(40, mrpt::round(m2->getArea() * options.featsPerSquareMeter));
159 
160  m_grid_feat_extr.extractFeatures(
161  *m1, *lm1, N1, options.feature_descriptor,
162  options.feature_detector_options);
163  m_grid_feat_extr.extractFeatures(
164  *m2, *lm2, N2, options.feature_descriptor,
165  options.feature_detector_options);
166 
167  if (runningTime)
168  {
169  tictac = new CTicTac();
170  tictac->Tic();
171  }
172 
173  const size_t nLM1 = lm1->size();
174  const size_t nLM2 = lm2->size();
175 
176  // At least two landmarks at each map!
177  // ------------------------------------------------------
178  if (nLM1 < 2 || nLM2 < 2)
179  {
180  outInfo.goodness = 0;
181  }
182  else
183  {
184  //#define METHOD_FFT
185  //#define DEBUG_SHOW_CORRELATIONS
186 
187  // Compute correlation between landmarks:
188  // ---------------------------------------------
189  CMatrixFloat CORR(lm1->size(), lm2->size()), auxCorr;
190  CImage im1, im2; // Grayscale
191  CVectorFloat corr;
192  unsigned int corrsCount = 0;
193  std::vector<bool> hasCorr(nLM1, false);
194 
195  const double thres_max = options.threshold_max;
196  const double thres_delta = options.threshold_delta;
197 
198  // CDisplayWindowPlots::Ptr auxWin;
199  if (options.debug_show_corrs)
200  {
201  // auxWin = CDisplayWindowPlots::Ptr( new
202  // CDisplayWindowPlots("Individual corr.") );
203  std::cerr << "Warning: options.debug_show_corrs has no effect "
204  "since MRPT 0.9.1\n";
205  }
206 
207  for (size_t idx1 = 0; idx1 < nLM1; idx1++)
208  {
209  // CVectorFloat corrs_indiv;
210  vector<pair<size_t, float>> corrs_indiv; // (index, distance);
211  // Index is used to
212  // recover the original
213  // index after sorting.
214  vector<float> corrs_indiv_only;
215  corrs_indiv.reserve(nLM2);
216  corrs_indiv_only.reserve(nLM2);
217 
218  for (size_t idx2 = 0; idx2 < nLM2; idx2++)
219  {
220  float minDist;
221  minDist =
222  lm1->landmarks.get(idx1)->features[0].descriptorDistanceTo(
223  lm2->landmarks.get(idx2)->features[0]);
224 
225  corrs_indiv.emplace_back(idx2, minDist);
226  corrs_indiv_only.push_back(minDist);
227  } // end for idx2
228 
229  // double corr_mean,corr_std;
230  // mrpt::math::meanAndStd(corrs_indiv_only,corr_mean,corr_std);
231  const double corr_best = mrpt::math::minimum(corrs_indiv_only);
232  // cout << "M: " << corr_mean << " std: " << corr_std << " best: "
233  // << corr_best << endl;
234 
235  // Sort the list and keep the N best features:
236  std::sort(corrs_indiv.begin(), corrs_indiv.end(), myVectorOrder);
237 
238  // const size_t nBestToKeep = std::min( (size_t)30,
239  // corrs_indiv.size() );
240  const size_t nBestToKeep = corrs_indiv.size();
241 
242  for (size_t w = 0; w < nBestToKeep; w++)
243  {
244  if (corrs_indiv[w].second <= thres_max &&
245  corrs_indiv[w].second <= (corr_best + thres_delta))
246  {
247  idxs1.push_back(idx1);
248  idxs2.push_back(corrs_indiv[w].first);
249  outInfo.correspondences_dists_maha.emplace_back(
250  idx1, corrs_indiv[w].first, corrs_indiv[w].second);
251  }
252  }
253  } // end for idx1
254 
255  // Save an image with each feature and its matches:
256  if (options.save_feat_coors)
257  {
259  mrpt::system::createDirectory("grid_feats");
260 
261  std::map<size_t, std::set<size_t>> corrs_by_idx;
262  for (size_t l = 0; l < idxs1.size(); l++)
263  corrs_by_idx[idxs1[l]].insert(idxs2[l]);
264 
265  for (auto& it : corrs_by_idx)
266  {
267  CMatrixFloat descriptor1;
268  lm1->landmarks.get(it.first)
269  ->features[0]
270  .getFirstDescriptorAsMatrix(descriptor1);
271 
272  im1.setFromMatrix(descriptor1, true /*normalized*/);
273 
274  const size_t FEAT_W = im1.getWidth();
275  const size_t FEAT_H = im1.getHeight();
276  const size_t nF = it.second.size();
277 
278  CImage img_compose(FEAT_W * 2 + 15, 10 + (5 + FEAT_H) * nF);
279  img_compose.filledRectangle(
280  0, 0, img_compose.getWidth() - 1,
281  img_compose.getHeight() - 1, TColor::black());
282 
283  img_compose.drawImage(5, 5, im1);
284 
285  size_t j;
286  std::set<size_t>::iterator it_j;
287  string fil =
288  format("grid_feats/_FEAT_MATCH_%03i", (int)it.first);
289 
290  for (j = 0, it_j = it.second.begin(); j < nF; ++j, ++it_j)
291  {
292  fil += format("_%u", static_cast<unsigned int>(*it_j));
293 
294  CMatrixFloat descriptor2;
295  lm2->landmarks.get(*it_j)
296  ->features[0]
297  .getFirstDescriptorAsMatrix(descriptor2);
298  im2.setFromMatrix(descriptor2, true);
299  img_compose.drawImage(
300  10 + FEAT_W, 5 + j * (FEAT_H + 5), im2);
301  }
302  fil += ".png";
303  img_compose.saveToFile(fil);
304  } // end for map
305  }
306 
307  // ------------------------------------------------------------------
308  // Create the list of correspondences from the lists: idxs1 & idxs2
309  // ------------------------------------------------------------------
310  correspondences.clear();
311  for (it1 = idxs1.begin(), it2 = idxs2.begin(); it1 != idxs1.end();
312  ++it1, ++it2)
313  {
315  mp.this_idx = *it1;
316  mp.this_x = lm1->landmarks.get(*it1)->pose_mean.x;
317  mp.this_y = lm1->landmarks.get(*it1)->pose_mean.y;
318  mp.this_z = lm1->landmarks.get(*it1)->pose_mean.z;
319 
320  mp.other_idx = *it2;
321  mp.other_x = lm2->landmarks.get(*it2)->pose_mean.x;
322  mp.other_y = lm2->landmarks.get(*it2)->pose_mean.y;
323  mp.other_z = lm2->landmarks.get(*it2)->pose_mean.z;
324  correspondences.push_back(mp);
325 
326  if (!hasCorr[*it1])
327  {
328  hasCorr[*it1] = true;
329  corrsCount++;
330  }
331  } // end for corres.
332 
333  outInfo.goodness = (2.0f * corrsCount) / (nLM1 + nLM2);
334 
335  // Compute the estimation using ALL the correspondences (NO ROBUST):
336  // ----------------------------------------------------------------------
337  mrpt::math::TPose2D noRobustEst;
338  if (!mrpt::tfest::se2_l2(correspondences, noRobustEst))
339  {
340  // There's no way to match the maps! e.g. no correspondences
341  outInfo.goodness = 0;
342  pdf_SOG->clear();
343  outInfo.sog1 = pdf_SOG;
344  outInfo.sog2 = pdf_SOG;
345  outInfo.sog3 = pdf_SOG;
346  }
347  else
348  {
349  outInfo.noRobustEstimation = mrpt::poses::CPose2D(noRobustEst);
351  "[CGridMapAligner] Overall estimation(%u corrs, total: "
352  "%u): (%.03f,%.03f,%.03fdeg)\n",
353  corrsCount, (unsigned)correspondences.size(),
354  outInfo.noRobustEstimation.x(), outInfo.noRobustEstimation.y(),
355  RAD2DEG(outInfo.noRobustEstimation.phi())));
356 
357  // The list of SOG modes & their corresponding sub-sets of
358  // matchings:
359  using TMapMatchingsToPoseMode = std::map<
361  TMapMatchingsToPoseMode sog_modes;
362 
363  // ---------------------------------------------------------------
364  // Now, we have to choose between the methods:
365  // - CGridMapAligner::amRobustMatch ("normal" RANSAC)
366  // - CGridMapAligner::amModifiedRANSAC
367  // ---------------------------------------------------------------
368  if (options.methodSelection == CGridMapAligner::amRobustMatch)
369  {
370  // ====================================================
371  // METHOD: "Normal" RANSAC
372  // ====================================================
373 
374  // RANSAC over the found correspondences:
375  // -------------------------------------------------
376  const unsigned int min_inliers =
377  options.ransac_minSetSizeRatio * (nLM1 + nLM2) / 2;
378 
379  mrpt::tfest::TSE2RobustParams tfest_params;
380  tfest_params.ransac_minSetSize = min_inliers;
381  tfest_params.ransac_maxSetSize = nLM1 + nLM2;
383  options.ransac_mahalanobisDistanceThreshold;
384  tfest_params.ransac_nSimulations = 0; // 0=auto
385  tfest_params.ransac_fuseByCorrsMatch = true;
386  tfest_params.ransac_fuseMaxDiffXY = 0.01;
387  tfest_params.ransac_fuseMaxDiffPhi = DEG2RAD(0.1);
388  tfest_params.ransac_algorithmForLandmarks = true;
389  tfest_params.probability_find_good_model =
390  options.ransac_prob_good_inliers;
391  tfest_params.verbose = false;
392 
393  mrpt::tfest::TSE2RobustResult tfest_result;
395  correspondences, options.ransac_SOG_sigma_m, tfest_params,
396  tfest_result);
397 
398  ASSERT_(pdf_SOG);
399  *pdf_SOG = tfest_result.transformation;
400  largestConsensusCorrs = tfest_result.largestSubSet;
401 
402  // Simplify the SOG by merging close modes:
403  // -------------------------------------------------
404  size_t nB = pdf_SOG->size();
405  outInfo.sog1 = pdf_SOG;
406 
407  // Keep only the most-likely Gaussian mode:
408  CPosePDFSOG::TGaussianMode best_mode;
409  best_mode.log_w = -std::numeric_limits<double>::max();
410 
411  for (size_t n = 0; n < pdf_SOG->size(); n++)
412  {
413  const CPosePDFSOG::TGaussianMode& m = (*pdf_SOG)[n];
414 
415  if (m.log_w > best_mode.log_w) best_mode = m;
416  }
417 
418  pdf_SOG->clear();
419  pdf_SOG->push_back(best_mode);
420  outInfo.sog2 = pdf_SOG;
421 
423  "[CGridMapAligner] amRobustMatch: "
424  << nB << " SOG modes reduced to " << pdf_SOG->size()
425  << " (most-likely) (min.inliers=" << min_inliers << ")\n");
426 
427  } // end of amRobustMatch
428  else
429  {
430  // ====================================================
431  // METHOD: "Modified" RANSAC
432  // ====================================================
433  mrpt::tfest::TMatchingPairList all_corrs = correspondences;
434 
435  const size_t nCorrs = all_corrs.size();
436  ASSERT_(nCorrs >= 2);
437 
438  pdf_SOG->clear(); // Start with 0 Gaussian modes
439 
440  // Build a points-map for exploiting its KD-tree:
441  // ----------------------------------------------------
442  CSimplePointsMap lm1_pnts, lm2_pnts;
443  lm1_pnts.reserve(nLM1);
444  for (size_t i = 0; i < nLM1; i++)
445  lm1_pnts.insertPoint(lm1->landmarks.get(i)->pose_mean);
446  lm2_pnts.reserve(nLM2);
447  for (size_t i = 0; i < nLM2; i++)
448  lm2_pnts.insertPoint(lm2->landmarks.get(i)->pose_mean);
449 
450  // RANSAC loop
451  // ---------------------
452  const size_t minInliersTOaccept =
453  round(options.ransac_minSetSizeRatio * 0.5 * (nLM1 + nLM2));
454  // Set an initial # of iterations:
455  const unsigned int ransac_min_nSimulations =
456  2 * (nLM1 + nLM2); // 1000;
457  unsigned int ransac_nSimulations =
458  10; // It doesn't matter actually, since will be changed in
459  // the first loop
460  const double probability_find_good_model = 0.9999;
461 
462  const double chi2_thres_dim1 =
463  mrpt::math::chi2inv(options.ransac_chi2_quantile, 1);
464  const double chi2_thres_dim2 =
465  mrpt::math::chi2inv(options.ransac_chi2_quantile, 2);
466 
467  // Generic 2x2 covariance matrix for all features in their local
468  // coords:
469  CMatrixDouble22 COV_pnt;
470  COV_pnt(0, 0) = COV_pnt(1, 1) =
471  square(options.ransac_SOG_sigma_m);
472 
473  // The absolute number of trials.
474  // in practice it's only important for a reduced number of
475  // correspondences, to avoid a dead-lock:
476  // It's the binomial coefficient:
477  // / n \ n! n * (n-1)
478  // | | = ----------- = -----------
479  // \ 2 / 2! (n-2)! 2
480  //
481  const unsigned int max_trials =
482  (nCorrs * (nCorrs - 1) / 2) *
483  5; // "*5" is just for safety...
484 
485  unsigned int iter = 0; // Valid iterations (those passing the
486  // first mahalanobis test)
487  unsigned int trials = 0; // counter of all iterations,
488  // including "iter" + failing ones.
489  while (iter < ransac_nSimulations &&
490  trials <
491  max_trials) // ransac_nSimulations can be dynamic
492  {
493  trials++;
494 
495  mrpt::tfest::TMatchingPairList tentativeSubSet;
496 
497  // Pick 2 random correspondences:
498  uint32_t idx1, idx2;
499  idx1 = getRandomGenerator().drawUniform32bit() % nCorrs;
500  do
501  {
502  idx2 = getRandomGenerator().drawUniform32bit() % nCorrs;
503  } while (idx1 == idx2); // Avoid a degenerated case!
504 
505  // Uniqueness of features:
506  if (all_corrs[idx1].this_idx == all_corrs[idx2].this_idx ||
507  all_corrs[idx1].this_idx == all_corrs[idx2].other_idx)
508  continue;
509  if (all_corrs[idx1].other_idx == all_corrs[idx2].this_idx ||
510  all_corrs[idx1].other_idx == all_corrs[idx2].other_idx)
511  continue;
512 
513  // Check the feasibility of this pair "idx1"-"idx2":
514  // The distance between the pair of points in MAP1 must be
515  // very close
516  // to that of their correspondences in MAP2:
517  const double corrs_dist1 =
519  all_corrs[idx1].this_x, all_corrs[idx1].this_y,
520  all_corrs[idx1].this_z, all_corrs[idx2].this_x,
521  all_corrs[idx2].this_y, all_corrs[idx2].this_z);
522 
523  const double corrs_dist2 =
525  all_corrs[idx1].other_x, all_corrs[idx1].other_y,
526  all_corrs[idx1].other_z, all_corrs[idx2].other_x,
527  all_corrs[idx2].other_y, all_corrs[idx2].other_z);
528 
529  // Is is a consistent possibility?
530  // We use a chi2 test (see paper for the derivation)
531  const double corrs_dist_chi2 =
532  square(square(corrs_dist1) - square(corrs_dist2)) /
533  (8.0 * square(options.ransac_SOG_sigma_m) *
534  (square(corrs_dist1) + square(corrs_dist2)));
535 
536  if (corrs_dist_chi2 > chi2_thres_dim1) continue; // Nope
537 
538  iter++; // Do not count iterations if they fail the test
539  // above.
540 
541  // before proceeding with this hypothesis, is it an old one?
542  bool is_new_hyp = true;
543  for (auto& sog_mode : sog_modes)
544  {
545  if (sog_mode.first.contains(all_corrs[idx1]) &&
546  sog_mode.first.contains(all_corrs[idx2]))
547  {
548  // Increment weight:
549  sog_mode.second.log_w =
550  std::log(std::exp(sog_mode.second.log_w) + 1.0);
551  is_new_hyp = false;
552  break;
553  }
554  }
555  if (!is_new_hyp) continue;
556 
557  // Ok, it's a new hypothesis:
558  tentativeSubSet.push_back(all_corrs[idx1]);
559  tentativeSubSet.push_back(all_corrs[idx2]);
560 
561  // Maintain a list of already used landmarks IDs in both
562  // maps to avoid repetitions:
563  std::vector<bool> used_landmarks1(nLM1, false);
564  std::vector<bool> used_landmarks2(nLM2, false);
565 
566  used_landmarks1[all_corrs[idx1].this_idx] = true;
567  used_landmarks1[all_corrs[idx2].this_idx] = true;
568  used_landmarks2[all_corrs[idx1].other_idx] = true;
569  used_landmarks2[all_corrs[idx2].other_idx] = true;
570 
571  // Build the transformation for these temptative
572  // correspondences:
573  bool keep_incorporating = true;
574  CPosePDFGaussian temptPose;
575  do // Incremently incorporate inliers:
576  {
577  if (!mrpt::tfest::se2_l2(tentativeSubSet, temptPose))
578  continue; // Invalid matching...
579 
580  // The computed cov is "normalized", i.e. must be
581  // multiplied by std^2_xy
582  temptPose.cov *= square(options.ransac_SOG_sigma_m);
583 
584  // cout << "q: " << temptPose << endl;
585 
586  // Find the landmark in MAP2 with the best (maximum)
587  // product-integral:
588  // (i^* , j^*) = arg max_(i,j) \int p_i()p_j()
589  //----------------------------------------------------------------------
590  const double ccos = cos(temptPose.mean.phi());
591  const double ssin = sin(temptPose.mean.phi());
592 
593  CMatrixDouble22 Hc; // Jacobian wrt point_j
594  Hc(1, 1) = ccos;
595  Hc(0, 0) = ccos;
596  Hc(1, 0) = ssin;
597  Hc(0, 1) = -ssin;
598 
600  Hq; // Jacobian wrt transformation q
601  Hq(0, 0) = 1;
602  Hq(1, 1) = 1;
603 
604  TPoint2D p2_j_local;
605  vector<float> matches_dist;
606  std::vector<size_t> matches_idx;
607 
608  CPoint2DPDFGaussian pdf_M2_j;
609  CPoint2DPDFGaussian pdf_M1_i;
610 
611 // Use integral-of-product vs. mahalanobis distances to match:
612 #define GRIDMAP_USE_PROD_INTEGRAL
613 
614 #ifdef GRIDMAP_USE_PROD_INTEGRAL
615  double best_pair_value = 0;
616 #else
617  double best_pair_value =
618  std::numeric_limits<double>::max();
619 #endif
620  double best_pair_d2 =
621  std::numeric_limits<double>::max();
622  pair<size_t, size_t> best_pair_ij;
623 
624  //#define SHOW_CORRS
625 
626 #ifdef SHOW_CORRS
627  CDisplayWindowPlots win("Matches");
628 #endif
629  for (size_t j = 0; j < nLM2; j++)
630  {
631  if (used_landmarks2[j]) continue;
632 
633  lm2_pnts.getPoint(
634  j, p2_j_local); // In local coords.
635  pdf_M2_j.mean = mrpt::poses::CPoint2D(
636  temptPose.mean +
637  p2_j_local); // In (temptative) global coords:
638  pdf_M2_j.cov(0, 0) = pdf_M2_j.cov(1, 1) =
639  square(options.ransac_SOG_sigma_m);
640 
641 #ifdef SHOW_CORRS
642  win.plotEllipse(
643  pdf_M2_j.mean.x(), pdf_M2_j.mean.y(),
644  pdf_M2_j.cov, 2, "b-",
645  format("M2_%u", (unsigned)j), true);
646 #endif
647 
648  static const unsigned int N_KDTREE_SEARCHED = 3;
649 
650  // Look for a few close features which may be
651  // potential matches:
652  lm1_pnts.kdTreeNClosestPoint2DIdx(
653  pdf_M2_j.mean.x(), pdf_M2_j.mean.y(),
654  N_KDTREE_SEARCHED, matches_idx, matches_dist);
655 
656  // And for each one, compute the product-integral:
657  for (unsigned long u : matches_idx)
658  {
659  if (used_landmarks1[u]) continue;
660 
661  // Jacobian wrt transformation q
662  Hq(0, 2) =
663  -p2_j_local.x * ssin - p2_j_local.y * ccos;
664  Hq(1, 2) =
665  p2_j_local.x * ccos - p2_j_local.y * ssin;
666 
667  // COV_j = Hq \Sigma_q Hq^t + Hc Cj Hc^t
668  pdf_M1_i.cov =
669  mrpt::math::multiply_HCHt(Hc, COV_pnt);
670  pdf_M1_i.cov += mrpt::math::multiply_HCHt(
671  Hq, temptPose.cov);
672 
673  float px, py;
674  lm1_pnts.getPoint(u, px, py);
675  pdf_M1_i.mean.x(px);
676  pdf_M1_i.mean.y(py);
677 
678 #ifdef SHOW_CORRS
679  win.plotEllipse(
680  pdf_M1_i.mean.x(), pdf_M1_i.mean.y(),
681  pdf_M1_i.cov, 2, "r-",
682  format("M1_%u", (unsigned)matches_idx[u]),
683  true);
684 #endif
685 
686 // And now compute the product integral:
687 #ifdef GRIDMAP_USE_PROD_INTEGRAL
688  const double prod_ij =
689  pdf_M1_i.productIntegralWith(pdf_M2_j);
690  // const double prod_ij_d2 = square(
691  // pdf_M1_i.mahalanobisDistanceTo(pdf_M2_j) );
692 
693  if (prod_ij > best_pair_value)
694 #else
695  const double prod_ij =
696  pdf_M1_i.mean.distanceTo(pdf_M2_j.mean);
697  if (prod_ij < best_pair_value)
698 #endif
699  {
700  // const double prodint_ij =
701  // pdf_M1_i.productIntegralWith2D(pdf_M2_j);
702 
703  best_pair_value = prod_ij;
704  best_pair_ij.first = u;
705  best_pair_ij.second = j;
706 
707  best_pair_d2 =
708  square(pdf_M1_i.mahalanobisDistanceTo(
709  pdf_M2_j));
710 
711  // cout << "P1: " << pdf_M1_i.mean << " C= "
712  // << pdf_M1_i.cov.inMatlabFormat() << endl;
713  // cout << "P2: " << pdf_M2_j.mean << " C= "
714  // << pdf_M2_j.cov.inMatlabFormat() << endl;
715  // cout << " -> " << format("%e",prod_ij)
716  // << " d2: " << best_pair_d2 << endl <<
717  // endl;
718  }
719  } // end for u (closest matches of LM2 in MAP 1)
720 
721 #ifdef SHOW_CORRS
722  win.axis_fit(true);
723  win.waitForKey();
724  win.clear();
725 #endif
726 
727  } // end for each LM2
728 
729  // Stop when the best choice has a bad mahal. dist.
730  keep_incorporating = false;
731 
732  // For the best (i,j), gate by the mahalanobis distance:
733  if (best_pair_d2 < chi2_thres_dim2)
734  {
735  // AND, also, check if the descriptors have some
736  // resemblance!!
737  // ----------------------------------------------------------------
738  // double feat_dist =
739  // lm1->landmarks.get(best_pair_ij.first)->features[0]->descriptorDistanceTo(*lm1->landmarks.get(best_pair_ij.second)->features[0]);
740  // if (feat_dist< options.threshold_max)
741  {
742  float p1_i_localx, p1_i_localy;
743  float p2_j_localx, p2_j_localy;
744  lm1_pnts.getPoint(
745  best_pair_ij.first, p1_i_localx,
746  p1_i_localy);
747  lm2_pnts.getPoint(
748  best_pair_ij.second, p2_j_localx,
749  p2_j_localy); // In local coords.
750 
751  used_landmarks1[best_pair_ij.first] = true;
752  used_landmarks2[best_pair_ij.second] = true;
753 
754  tentativeSubSet.push_back(
756  best_pair_ij.first, best_pair_ij.second,
757  p1_i_localx, p1_i_localy, 0, // MAP1
758  p2_j_localx, p2_j_localy, 0 // MAP2
759  ));
760 
761  keep_incorporating = true;
762  }
763  }
764 
765  } while (keep_incorporating);
766 
767  // Consider this pairing?
768  const size_t ninliers = tentativeSubSet.size();
769  if (ninliers > minInliersTOaccept)
770  {
772  newGauss.log_w = 0; // log(1); //
773  // std::log(static_cast<double>(nCoincidences));
774  newGauss.mean = temptPose.mean;
775  newGauss.cov = temptPose.cov;
776 
777  sog_modes[tentativeSubSet] = newGauss;
778 
779  // cout << "ITER: " << iter << " #inliers: " << ninliers
780  // << " q: " << temptPose.mean << endl;
781  }
782 
783  // Keep the largest consensus & dynamic # of steps:
784  if (largestConsensusCorrs.size() < ninliers)
785  {
786  largestConsensusCorrs = tentativeSubSet;
787 
788  // Update estimate of N, the number of trials to ensure
789  // we pick,
790  // with probability p, a data set with no outliers.
791  const double fracinliers =
792  ninliers /
793  static_cast<double>(std::min(nLM1, nLM2));
794  double pNoOutliers =
795  1 - pow(fracinliers,
796  static_cast<double>(
797  2.0)); // minimumSizeSamplesToFit
798 
799  pNoOutliers = std::max(
800  std::numeric_limits<double>::epsilon(),
801  pNoOutliers); // Avoid division by -Inf
802  pNoOutliers = std::min(
803  1.0 - std::numeric_limits<double>::epsilon(),
804  pNoOutliers); // Avoid division by 0.
805  // Number of
806  ransac_nSimulations =
807  log(1 - probability_find_good_model) /
808  log(pNoOutliers);
809 
810  if (ransac_nSimulations < ransac_min_nSimulations)
811  ransac_nSimulations = ransac_min_nSimulations;
812 
813  // if (verbose)
814  // cout << "[Align] Iter #" << iter << " Est. #iters: "
815  //<< ransac_nSimulations << " pNoOutliers=" <<
816  // pNoOutliers << " #inliers: " << ninliers << endl;
817  }
818 
819  } // end of RANSAC loop
820 
821  // Move SOG modes into pdf_SOG:
822  pdf_SOG->clear();
823  for (auto s = sog_modes.begin(); s != sog_modes.end(); ++s)
824  {
825  cout << "SOG mode: " << s->second.mean
826  << " inliers: " << s->first.size() << endl;
827  pdf_SOG->push_back(s->second);
828  }
829 
830  // Normalize:
831  if (pdf_SOG->size() > 0) pdf_SOG->normalizeWeights();
832 
833  // Simplify the SOG by merging close modes:
834  // -------------------------------------------------
835  size_t nB = pdf_SOG->size();
836  outInfo.sog1 = pdf_SOG;
837 
838  CTicTac merge_clock;
839  pdf_SOG->mergeModes(options.maxKLd_for_merge, false);
840  const double merge_clock_tim = merge_clock.Tac();
841 
842  outInfo.sog2 = pdf_SOG;
843  size_t nA = pdf_SOG->size();
845  "[CGridMapAligner] amModifiedRANSAC: %u SOG modes "
846  "merged to %u in %.03fsec\n",
847  (unsigned int)nB, (unsigned int)nA, merge_clock_tim));
848 
849  } // end of amModifiedRANSAC
850 
851  // Save best corrs:
852  if (options.debug_save_map_pairs)
853  {
854  static unsigned int NN = 0;
855  static const COccupancyGridMap2D* lastM1 = nullptr;
856  if (lastM1 != m1)
857  {
858  lastM1 = m1;
859  NN = 0;
860  }
861  printf(
862  " Largest consensus: %u\n",
863  static_cast<unsigned>(largestConsensusCorrs.size()));
864  CEnhancedMetaFile::LINUX_IMG_WIDTH(
865  m1->getSizeX() + m2->getSizeX() + 50);
866  CEnhancedMetaFile::LINUX_IMG_HEIGHT(
867  max(m1->getSizeY(), m2->getSizeY()) + 50);
868 
869  for (auto s = sog_modes.begin(); s != sog_modes.end(); ++s)
870  {
871  COccupancyGridMap2D::saveAsEMFTwoMapsWithCorrespondences(
872  format("__debug_corrsGrid_%05u.emf", NN), m1, m2,
873  s->first);
874  COccupancyGridMap2D::saveAsBitmapTwoMapsWithCorrespondences(
875  format("__debug_corrsGrid_%05u.png", NN), m1, m2,
876  s->first);
877  ++NN;
878  }
879  }
880 
881  // --------------------------------------------------------------------
882  // At this point:
883  // - "pdf_SOG": has the resulting PDF with the SOG (from whatever
884  // method)
885  // - "largestConsensusCorrs": The 'best' set of correspondences
886  //
887  // Now: If we had a multi-metric map, use the points map to improve
888  // the estimation with ICP.
889  // --------------------------------------------------------------------
890  if (multimap1 && multimap2 &&
891  multimap1->countMapsByClass<CSimplePointsMap>() != 0 &&
892  multimap2->countMapsByClass<CSimplePointsMap>() != 0)
893  {
894  auto pnts1 = multimap1->mapByClass<CSimplePointsMap>();
895  auto pnts2 = multimap2->mapByClass<CSimplePointsMap>();
896 
897  CICP icp;
898  CICP::TReturnInfo icpInfo;
899 
900  icp.options.maxIterations = 20;
901  icp.options.smallestThresholdDist = 0.05f;
902  icp.options.thresholdDist = 0.75f;
903 
904  // Invoke ICP once for each mode in the SOG:
905  size_t cnt = 0;
906  outInfo.icp_goodness_all_sog_modes.clear();
907  for (auto i = pdf_SOG->begin(); i != pdf_SOG->end(); ++cnt)
908  {
909  CPosePDF::Ptr icp_est = icp.Align(
910  pnts1.get(), pnts2.get(), (i)->mean, nullptr, &icpInfo);
911 
912  //(i)->cov(0,0) += square( 0.05 );
913  //(i)->cov(1,1) += square( 0.05 );
914  //(i)->cov(2,2) += square( DEG2RAD(0.05) );
915 
916  CPosePDFGaussian i_gauss(i->mean, i->cov);
917  CPosePDFGaussian icp_gauss(icp_est->getMeanVal(), i->cov);
918 
919  const double icp_maha_dist =
920  i_gauss.mahalanobisDistanceTo(icp_gauss);
921 
922  cout << "ICP " << cnt << " log-w: " << i->log_w
923  << " Goodness: " << 100 * icpInfo.goodness
924  << " D_M= " << icp_maha_dist << endl;
925  cout << " final pos: " << icp_est->getMeanVal() << endl;
926  cout << " org pos: " << i->mean << endl;
927 
928  outInfo.icp_goodness_all_sog_modes.push_back(
929  icpInfo.goodness);
930 
931  // Discard?
932  if (icpInfo.goodness > options.min_ICP_goodness &&
933  icp_maha_dist < options.max_ICP_mahadist)
934  {
935  icp_est->getMean((i)->mean);
936  ++i;
937  }
938  else
939  {
940  // Delete:
941  i = pdf_SOG->erase(i);
942  }
943 
944  } // end for i
945 
946  // Merge:
947  outInfo.sog3 = pdf_SOG;
948 
949  pdf_SOG->mergeModes(options.maxKLd_for_merge, false);
951  "[CGridMapAligner] " << pdf_SOG->size()
952  << " SOG modes merged after ICP.");
953 
954  } // end multimapX
955 
956  } // end of, yes, we have enough correspondences
957 
958  } // end of: yes, there are landmarks in the grid maps!
959 
960  // Copy the output info if requested:
961  // -------------------------------------------------
962  MRPT_TODO(
963  "Refactor `info` so it is polymorphic and can use dynamic_cast<> here");
964  if (info)
965  {
966  auto* info_ = static_cast<TReturnInfo*>(info);
967  *info_ = outInfo;
968  }
969 
970  if (runningTime)
971  {
972  *runningTime = tictac->Tac();
973  delete tictac;
974  }
975 
976  return pdf_SOG;
977 
978  MRPT_END
979 }
980 
981 /*---------------------------------------------------------------
982  AlignPDF_correlation
983 ---------------------------------------------------------------*/
984 CPosePDF::Ptr CGridMapAligner::AlignPDF_correlation(
985  const mrpt::maps::CMetricMap* mm1, const mrpt::maps::CMetricMap* mm2,
986  const CPosePDFGaussian& initialEstimationPDF, float* runningTime,
987  void* info)
988 {
989  MRPT_UNUSED_PARAM(initialEstimationPDF);
990  MRPT_UNUSED_PARAM(info);
991 
992  MRPT_START
993 
994  //#define CORRELATION_SHOW_DEBUG
995 
996  CTicTac* tictac = nullptr;
997 
998  // Asserts:
999  // -----------------
1002  const auto* m1 = static_cast<const COccupancyGridMap2D*>(mm1);
1003  const auto* m2 = static_cast<const COccupancyGridMap2D*>(mm2);
1004 
1005  ASSERT_(m1->getResolution() == m2->getResolution());
1006 
1007  if (runningTime)
1008  {
1009  tictac = new CTicTac();
1010  tictac->Tic();
1011  }
1012 
1013  // The PDF to estimate:
1014  // ------------------------------------------------------
1015  CPosePDFGaussian::Ptr PDF = std::make_shared<CPosePDFGaussian>();
1016 
1017  // Determine the extension to compute the correlation into:
1018  // ----------------------------------------------------------
1019  float phiResolution = DEG2RAD(0.2f);
1020  float phiMin = -M_PIf + 0.5f * phiResolution;
1021  float phiMax = M_PIf;
1022 
1023  // Compute the difference between maps for each (u,v) pair!
1024  // --------------------------------------------------------------
1025  float phi;
1026  float pivotPt_x = 0.5f * (m1->getXMax() + m1->getXMin());
1027  float pivotPt_y = 0.5f * (m1->getYMax() + m1->getYMin());
1028  COccupancyGridMap2D map2_mod;
1029  CImage map1_img, map2_img;
1030  float currentMaxCorr = -1e6f;
1031 
1032  m1->getAsImage(map1_img);
1033 
1034  map2_mod.setSize(
1035  m1->getXMin(), m1->getXMax(), m1->getYMin(), m1->getYMax(),
1036  m1->getResolution());
1037  size_t map2_lx = map2_mod.getSizeX();
1038  size_t map2_ly = map2_mod.getSizeY();
1039  CMatrixF outCrossCorr, bestCrossCorr;
1040  float map2width_2 = 0.5f * (map2_mod.getXMax() - map2_mod.getXMin());
1041 
1042 #ifdef CORRELATION_SHOW_DEBUG
1043  CDisplayWindow* win = new CDisplayWindow("corr");
1044  CDisplayWindow* win2 = new CDisplayWindow("map2_mod");
1045 #endif
1046 
1047  // --------------------------------------------------------
1048  // Use FFT-based correlation for each orientation:
1049  // --------------------------------------------------------
1050  for (phi = phiMin; phi < phiMax; phi += phiResolution)
1051  {
1052  // Create the displaced map2 grid, for the size of map1:
1053  // --------------------------------------------------------
1054  CPose2D v2(
1055  pivotPt_x - cos(phi) * map2width_2,
1056  pivotPt_y - sin(phi) * map2width_2,
1057  phi); // Rotation point: the centre of img1:
1058  CPoint2D v1, v3;
1059  v2 = CPose2D(0, 0, 0) - v2; // Inverse
1060 
1061  for (unsigned int cy2 = 0; cy2 < map2_ly; cy2++)
1062  {
1063  COccupancyGridMap2D::cellType* row = map2_mod.getRow(cy2);
1064  for (unsigned int cx2 = 0; cx2 < map2_lx; cx2++)
1065  {
1066  v3 = v2 + CPoint2D(map2_mod.idx2x(cx2), map2_mod.idx2y(cy2));
1067  *row++ = m2->p2l(m2->getPos(v3.x(), v3.y()));
1068  }
1069  }
1070 
1071  map2_mod.getAsImage(map2_img);
1072  // map2_img.saveToBMP("__debug_map2mod.bmp");
1073 
1074  // Compute the correlation:
1075  map1_img.cross_correlation_FFT(
1076  map2_img, outCrossCorr, -1, -1, -1, -1,
1077  127, // Bias to be substracted
1078  127 // Bias to be substracted
1079  );
1080 
1081  float corrPeak = outCrossCorr.maxCoeff();
1082  printf("phi = %fdeg \tmax corr=%f\n", RAD2DEG(phi), corrPeak);
1083 
1084  // Keep the maximum:
1085  if (corrPeak > currentMaxCorr)
1086  {
1087  currentMaxCorr = corrPeak;
1088  bestCrossCorr = outCrossCorr;
1089  PDF->mean.phi(phi);
1090  }
1091 
1092 #ifdef CORRELATION_SHOW_DEBUG
1093  outCrossCorr.adjustRange(0, 1);
1094  CImage aux(outCrossCorr, true);
1095  win->showImage(aux);
1096  win2->showImage(map2_img);
1097  std::this_thread::sleep_for(5ms);
1098 #endif
1099 
1100  } // end for phi
1101 
1102  if (runningTime)
1103  {
1104  *runningTime = tictac->Tac();
1105  delete tictac;
1106  }
1107 
1108  CImage aux;
1109  aux.setFromMatrix(bestCrossCorr, false /* do normalization [0,1]*/);
1110  aux.saveToFile("_debug_best_corr.png");
1111 
1112 #ifdef CORRELATION_SHOW_DEBUG
1113  delete win;
1114  delete win2;
1115 #endif
1116 
1117  // Transform the best corr matrix peak into coordinates:
1118  std::size_t uMax, vMax;
1119  currentMaxCorr = bestCrossCorr.maxCoeff(uMax, vMax);
1120 
1121  PDF->mean.x(m1->idx2x(uMax));
1122  PDF->mean.y(m1->idx2y(vMax));
1123 
1124  return PDF;
1125 
1126  // Done!
1127  MRPT_END
1128 }
1129 
1130 /*---------------------------------------------------------------
1131  TConfigParams
1132  ---------------------------------------------------------------*/
1133 CGridMapAligner::TConfigParams::TConfigParams() : feature_detector_options() {}
1134 /*---------------------------------------------------------------
1135  dumpToTextStream
1136  ---------------------------------------------------------------*/
1138 {
1139  out << mrpt::format(
1140  "\n----------- [CGridMapAligner::TConfigParams] ------------ \n\n");
1141 
1142  LOADABLEOPTS_DUMP_VAR(methodSelection, int)
1143  LOADABLEOPTS_DUMP_VAR(featsPerSquareMeter, float)
1144  LOADABLEOPTS_DUMP_VAR(threshold_max, float)
1145  LOADABLEOPTS_DUMP_VAR(threshold_delta, float)
1146  LOADABLEOPTS_DUMP_VAR(min_ICP_goodness, float)
1147  LOADABLEOPTS_DUMP_VAR(max_ICP_mahadist, double)
1148  LOADABLEOPTS_DUMP_VAR(maxKLd_for_merge, float)
1149  LOADABLEOPTS_DUMP_VAR(ransac_minSetSizeRatio, float)
1151  LOADABLEOPTS_DUMP_VAR(ransac_chi2_quantile, double)
1152  LOADABLEOPTS_DUMP_VAR(ransac_prob_good_inliers, double)
1153  LOADABLEOPTS_DUMP_VAR(ransac_SOG_sigma_m, float)
1154  LOADABLEOPTS_DUMP_VAR(save_feat_coors, bool)
1155  LOADABLEOPTS_DUMP_VAR(debug_show_corrs, bool)
1156  LOADABLEOPTS_DUMP_VAR(debug_save_map_pairs, bool)
1157 
1158  LOADABLEOPTS_DUMP_VAR(feature_descriptor, int)
1159 
1160  feature_detector_options.dumpToTextStream(out);
1161 
1162  out << mrpt::format("\n");
1163 }
1164 
1165 /*---------------------------------------------------------------
1166  loadFromConfigFile
1167  ---------------------------------------------------------------*/
1169  const mrpt::config::CConfigFileBase& iniFile, const std::string& section)
1170 {
1171  methodSelection =
1172  iniFile.read_enum(section, "methodSelection", methodSelection);
1173 
1175  featsPerSquareMeter, float, iniFile, section)
1176  MRPT_LOAD_CONFIG_VAR(ransac_SOG_sigma_m, float, iniFile, section)
1177 
1178  MRPT_LOAD_CONFIG_VAR_NO_DEFAULT(threshold_max, float, iniFile, section)
1179  MRPT_LOAD_CONFIG_VAR_NO_DEFAULT(threshold_delta, float, iniFile, section)
1180 
1181  MRPT_LOAD_CONFIG_VAR_NO_DEFAULT(min_ICP_goodness, float, iniFile, section)
1182  MRPT_LOAD_CONFIG_VAR_NO_DEFAULT(max_ICP_mahadist, double, iniFile, section)
1183 
1184  MRPT_LOAD_CONFIG_VAR_NO_DEFAULT(maxKLd_for_merge, float, iniFile, section)
1186  ransac_minSetSizeRatio, float, iniFile, section)
1190  ransac_chi2_quantile, double, iniFile, section)
1192  ransac_prob_good_inliers, double, iniFile, section)
1193 
1194  MRPT_LOAD_CONFIG_VAR(save_feat_coors, bool, iniFile, section)
1195  MRPT_LOAD_CONFIG_VAR(debug_show_corrs, bool, iniFile, section)
1196  MRPT_LOAD_CONFIG_VAR(debug_save_map_pairs, bool, iniFile, section)
1197 
1198  feature_descriptor = iniFile.read_enum(
1199  section, "feature_descriptor", feature_descriptor, true);
1200  feature_detector_options.loadFromConfigFile(iniFile, section);
1201 }
1202 
1204  const mrpt::maps::CMetricMap* m1, const mrpt::maps::CMetricMap* m2,
1205  const CPose3DPDFGaussian& initialEstimationPDF, float* runningTime,
1206  void* info)
1207 {
1208  MRPT_UNUSED_PARAM(m1);
1209  MRPT_UNUSED_PARAM(m2);
1210  MRPT_UNUSED_PARAM(initialEstimationPDF);
1211  MRPT_UNUSED_PARAM(runningTime);
1212  MRPT_UNUSED_PARAM(info);
1213  THROW_EXCEPTION("Align3D method not applicable to gridmap-aligner");
1214 }
Scalar maxCoeff() const
Maximum value in the matrix/vector.
A namespace of pseudo-random numbers generators of diferent distributions.
uint32_t drawUniform32bit()
Generate a uniformly distributed pseudo-random number using the MT19937 algorithm, in the whole range of 32-bit integers.
float getXMax() const
Returns the "x" coordinate of right side of grid map.
std::size_t countMapsByClass() const
Count how many maps exist of the given class (or derived class)
double Tac() noexcept
Stops the stopwatch.
Definition: CTicTac.cpp:86
double x
X,Y coordinates.
Definition: TPoint2D.h:23
void getPoint(size_t index, float &x, float &y, float &z) const
Access to a given point from map, as a 2D point.
Definition: CPointsMap.cpp:200
double getArea() const
Returns the area of the gridmap, in square meters.
A compile-time fixed-size numeric matrix container.
Definition: CMatrixFixed.h:33
bool createDirectory(const std::string &dirName)
Creates a directory.
Definition: filesystem.cpp:161
#define MRPT_START
Definition: exceptions.h:241
double productIntegralWith(const CPoint2DPDFGaussian &p) const
Computes the "correspondence likelihood" of this PDF with another one: This is implemented as the int...
CPose2D mean
The mean value.
std::vector< TPairPlusDistance > correspondences_dists_maha
Mahalanobis distance for each potential correspondence.
#define min(a, b)
double RAD2DEG(const double x)
Radians to degrees.
The struct for each mode:
Definition: CPosePDFSOG.h:41
float getResolution() const
Returns the resolution of the grid map.
Template for column vectors of dynamic size, compatible with Eigen.
#define THROW_EXCEPTION(msg)
Definition: exceptions.h:67
void reserve(size_t newLength) override
Reserves memory for a given number of points: the size of the map does not change, it only reserves the memory.
A gaussian distribution for 2D points.
bool se2_l2(const mrpt::tfest::TMatchingPairList &in_correspondences, mrpt::math::TPose2D &out_transformation, mrpt::math::CMatrixDouble33 *out_estimateCovariance=nullptr)
Least-squares (L2 norm) solution to finding the optimal SE(2) (x,y,yaw) between two reference frames...
Definition: se2_l2.cpp:46
Several implementations of ICP (Iterative closest point) algorithms for aligning two point maps or a ...
Definition: CICP.h:64
GLint * first
Definition: glext.h:3833
Parameters for se2_l2_robust().
Definition: se2.h:73
double DEG2RAD(const double x)
Degrees to radians.
void dumpToTextStream(std::ostream &out) const override
This method should clearly display all the contents of the structure in textual form, sending it to a std::ostream.
GLenum GLsizei n
Definition: glext.h:5136
This file implements several operations that operate element-wise on individual or pairs of container...
bool ransac_fuseByCorrsMatch
(Default = true) If true, the weight of Gaussian modes will be increased when an exact match in the s...
Definition: se2.h:92
unsigned int ransac_nSimulations
(Default = 0) If set to 0, an adaptive algorithm is used to determine the number of iterations...
Definition: se2.h:85
void drawImage(int x, int y, const mrpt::img::CImage &img) override
Draws an image as a bitmap at a given position.
Definition: CImage.cpp:1160
A high-performance stopwatch, with typical resolution of nanoseconds.
unsigned int ransac_minSetSize
(Default=3)
Definition: se2.h:76
mrpt::containers::deepcopy_poly_ptr< mrpt::poses::CPosePDFSOG::Ptr > sog3
const float ransac_mahalanobisDistanceThreshold
bool myVectorOrder(const pair< size_t, float > &o1, const pair< size_t, float > &o2)
size_t getHeight() const override
Returns the height of the image in pixels.
Definition: CImage.cpp:878
A cloud of points in 2D or 3D, which can be built from a sequence of laser scans. ...
void insertPoint(float x, float y, float z=0)
Provides a way to insert (append) individual points into the map: the missing fields of child classes...
Definition: CPointsMap.h:625
STL namespace.
double mahalanobisDistanceTo(const CPosePDFGaussian &theOther)
Computes the Mahalanobis distance between the centers of two Gaussians.
void kdTreeNClosestPoint2DIdx(float x0, float y0, size_t knn, std::vector< size_t > &out_idx, std::vector< float > &out_dist_sqr) const
KD Tree-based search for the N closest point to some given 2D coordinates and returns their indexes...
GLdouble s
Definition: glext.h:3682
The ICP algorithm return information.
mrpt::math::CMatrixDouble33 cov
The 3x3 covariance matrix.
GLubyte GLubyte GLubyte GLubyte w
Definition: glext.h:4199
double distanceTo(const CPoseOrPoint< OTHERCLASS, DIM2 > &b) const
Returns the Euclidean distance to another pose/point:
Definition: CPoseOrPoint.h:214
mrpt::system::CTicTac CTicTac
Definition: utils/CTicTac.h:5
mrpt::poses::CPose3DPDF::Ptr Align3DPDF(const mrpt::maps::CMetricMap *m1, const mrpt::maps::CMetricMap *m2, const mrpt::poses::CPose3DPDFGaussian &initialEstimationPDF, float *runningTime=nullptr, void *info=nullptr) override
Not applicable in this class, will launch an exception.
bool verbose
(Default=false)
Definition: se2.h:111
T square(const T x)
Inline function for the square of a number.
#define ASSERT_(f)
Defines an assertion mechanism.
Definition: exceptions.h:120
float goodness
A goodness measure for the alignment, it is a [0,1] range indicator of percentage of correspondences...
Definition: CICP.h:197
mrpt::tfest::TMatchingPairList correspondences
All the found correspondences (not consistent)
This class allows loading and storing values and vectors of different types from a configuration text...
This base provides a set of functions for maths stuff.
size_t getWidth() const override
Returns the width of the image in pixels.
Definition: CImage.cpp:847
void loadFromConfigFile(const mrpt::config::CConfigFileBase &source, const std::string &section) override
This method load the options from a ".ini"-like file or memory-stored string list.
bool se2_l2_robust(const mrpt::tfest::TMatchingPairList &in_correspondences, const double in_normalizationStd, const TSE2RobustParams &in_ransac_params, TSE2RobustResult &out_results)
Robust least-squares (L2 norm) solution to finding the optimal SE(2) (x,y,yaw) between two reference ...
#define CLASS_ID(T)
Access to runtime class ID for a defined class name.
Definition: CObject.h:89
unsigned int ransac_maxSetSize
(Default = 20)
Definition: se2.h:78
A class for storing a map of 3D probabilistic landmarks.
Definition: CLandmarksMap.h:74
Declares a class that represents a Probability Density function (PDF) of a 2D pose ...
float getXMin() const
Returns the "x" coordinate of left side of grid map.
virtual void filledRectangle(int x0, int y0, int x1, int y1, const mrpt::img::TColor color)
Draws a filled rectangle.
Definition: CCanvas.cpp:213
float idx2x(const size_t cx) const
Transform a cell index into a coordinate value.
A list of TMatchingPair.
Definition: TMatchingPair.h:70
bool ransac_algorithmForLandmarks
(Default = true) Use Mahalanobis distance (true) or Euclidean dist (false)
Definition: se2.h:99
mrpt::math::CMatrixDouble33 cov
Definition: CPosePDFSOG.h:45
GLfloat GLfloat GLfloat GLfloat v3
Definition: glext.h:4125
#define M_PIf
Definition: common.h:61
string iniFile(myDataDir+string("benchmark-options.ini"))
unsigned int getSizeX() const
Returns the horizontal size of grid map in cells count.
void setFromMatrix(const MAT &m, bool matrix_is_normalized=true)
Set the image from a matrix, interpreted as grayscale intensity values, in the range [0...
Definition: img/CImage.h:827
Classes for computer vision, detectors, features, etc.
Definition: CDifodo.h:17
#define IS_CLASS(obj, class_name)
True if the given reference to object (derived from mrpt::rtti::CObject) is of the given class...
Definition: CObject.h:133
#define MRPT_LOAD_CONFIG_VAR_NO_DEFAULT( variableName, variableType, configFileObject, sectionNameStr)
#define ASSERTMSG_(f, __ERROR_MSG)
Defines an assertion mechanism.
Definition: exceptions.h:108
double x() const
Common members of all points & poses classes.
Definition: CPoseOrPoint.h:143
#define MRPT_LOG_DEBUG_STREAM(__CONTENTS)
Use: MRPT_LOG_DEBUG_STREAM("Var=" << value << " foo=" << foo_var);
mrpt::gui::CDisplayWindow3D::Ptr win
GLsizei const GLchar ** string
Definition: glext.h:4116
double ransac_mahalanobisDistanceThreshold
(Default = 3.0)
Definition: se2.h:80
A class used to store a 2D point.
Definition: CPoint2D.h:32
double ransac_fuseMaxDiffPhi
(Default=0.1degree) (In radians)
Definition: se2.h:96
Classes for 2D/3D geometry representation, both of single values and probability density distribution...
A structure for holding correspondences between two sets of points or points-like entities in 2D or 3...
Definition: TMatchingPair.h:31
float goodness
A goodness measure for the alignment, it is a [0,1] range indicator of percentage of correspondences...
void cross_correlation_FFT(const CImage &in_img, math::CMatrixFloat &out_corr, int u_search_ini=-1, int v_search_ini=-1, int u_search_size=-1, int v_search_size=-1, float biasThisImg=0, float biasInImg=0) const
Computes the correlation matrix between this image and another one.
Definition: CImage.cpp:1375
#define LOADABLEOPTS_DUMP_VAR(variableName, variableType)
Macro for dumping a variable to a stream, within the method "dumpToTextStream(out)" (Variable types a...
#define MRPT_LOAD_CONFIG_VAR( variableName, variableType, configFileObject, sectionNameStr)
An useful macro for loading variables stored in a INI-like file under a key with the same name that t...
#define MRPT_LOG_INFO_STREAM(__CONTENTS)
CONTAINER::Scalar minimum(const CONTAINER &v)
mrpt::maps::CLandmarksMap::Ptr landmarks_map1
The landmarks of each map (the indices of these landmarks correspond to those in "correspondences") ...
A class for storing an occupancy grid map.
void getAsImage(mrpt::img::CImage &img, bool verticalFlip=false, bool forceRGB=false, bool tricolor=false) const
Returns the grid as a 8-bit graylevel image, where each pixel is a cell (output image is RGB only if ...
std::vector< double > icp_goodness_all_sog_modes
The ICP goodness of all potential SOG modes at the stage "sog2", thus before the removing of "bad" IC...
This class is a "CSerializable" wrapper for "CMatrixFloat".
Definition: CMatrixF.h:22
virtual const mrpt::rtti::TRuntimeClassId * GetRuntimeClass() const override
Returns information about the class of an object in runtime.
#define MRPT_TODO(x)
Definition: common.h:129
double chi2inv(double P, unsigned int dim=1)
The "quantile" of the Chi-Square distribution, for dimension "dim" and probability 0<P<1 (the inverse...
Definition: math.cpp:43
GLfloat GLfloat v1
Definition: glext.h:4121
Declares a virtual base class for all metric maps storage classes.
Definition: CMetricMap.h:52
A class used to store a 2D pose, including the 2D coordinate point and a heading (phi) angle...
Definition: CPose2D.h:39
cellType * getRow(int cy)
Access to a "row": mainly used for drawing grid as a bitmap efficiently, do not use it normally...
T::Ptr mapByClass(const size_t &ith=0) const
Returns the i&#39;th map of a given class (or of a derived class), or empty smart pointer if there is no ...
std::string format(const char *fmt,...) MRPT_printf_format_check(1
A std::string version of C sprintf.
Definition: format.cpp:16
const double & phi() const
Get the phi angle of the 2D pose (in radians)
Definition: CPose2D.h:86
GLenum GLenum GLvoid * row
Definition: glext.h:3580
#define MRPT_END
Definition: exceptions.h:245
The ICP algorithm return information.
Definition: CICP.h:190
Lightweight 2D pose.
Definition: TPose2D.h:22
double mean(const CONTAINER &v)
Computes the mean value of a vector.
float idx2y(const size_t cy) const
mrpt::poses::CPosePDFSOG transformation
The output as a set of transformations (sum of Gaussians)
Definition: se2.h:133
Declares a class that represents a Probability Density function (PDF) of a 3D pose ...
OccGridCellTraits::cellType cellType
The type of the map cells:
bool deleteFilesInDirectory(const std::string &s, bool deleteDirectoryAsWell=false)
Delete all the files in a given directory (nothing done if directory does not exists, or path is a file).
Definition: filesystem.cpp:218
mrpt::maps::CLandmarksMap::Ptr landmarks_map2
T distanceBetweenPoints(const T x1, const T y1, const T x2, const T y2)
Returns the distance between 2 points in 2D.
Definition: geometry.h:942
GLfloat GLfloat GLfloat v2
Definition: glext.h:4123
bool saveToFile(const std::string &fileName, int jpeg_quality=95) const
Save the image to a file, whose format is determined from the extension (internally uses OpenCV)...
Definition: CImage.cpp:336
void setSize(float x_min, float x_max, float y_min, float y_max, float resolution, float default_value=0.5f)
Change the size of gridmap, erasing all its previous contents.
void Tic() noexcept
Starts the stopwatch.
Definition: CTicTac.cpp:75
mrpt::containers::deepcopy_poly_ptr< mrpt::poses::CPosePDFSOG::Ptr > sog2
mrpt::tfest::TMatchingPairList largestSubSet
the largest consensus sub-set
Definition: se2.h:135
This class stores any customizable set of metric maps.
unsigned int getSizeY() const
Returns the vertical size of grid map in cells count.
unsigned __int32 uint32_t
Definition: rptypes.h:50
CRandomGenerator & getRandomGenerator()
A static instance of a CRandomGenerator class, for use in single-thread applications.
Lightweight 2D point.
Definition: TPoint2D.h:31
mrpt::poses::CPose2D noRobustEstimation
The "brute" estimation from using all the available correspondences, provided just for comparison pur...
This template class provides the basic functionality for a general 2D any-size, resizable container o...
mrpt::containers::deepcopy_poly_ptr< mrpt::poses::CPosePDFSOG::Ptr > sog1
The different SOG densities at different steps of the algorithm:
mrpt::math::CMatrixDouble22 cov
The 2x2 covariance matrix.
Output placeholder for se2_l2_robust()
Definition: se2.h:130
double ransac_fuseMaxDiffXY
(Default = 0.01)
Definition: se2.h:94
A class for storing images as grayscale or RGB bitmaps.
Definition: img/CImage.h:147
double probability_find_good_model
(Default = 0.999) See parameter ransac_nSimulations.
Definition: se2.h:103
double mahalanobisDistanceTo(const CPoint2DPDFGaussian &other) const
Returns the Mahalanobis distance from this PDF to another PDF, that is, it&#39;s evaluation at (0...
#define MRPT_LOG_INFO(_STRING)
#define MRPT_UNUSED_PARAM(a)
Determines whether this is an X86 or AMD64 platform.
Definition: common.h:186
void multiply_HCHt(const MAT_H &H, const MAT_C &C, MAT_R &R, bool accumResultInOutput=false)
R = H * C * H^t.
Definition: ops_matrices.h:28
int round(const T value)
Returns the closer integer (int) to x.
Definition: round.h:23



Page generated by Doxygen 1.8.14 for MRPT 1.9.9 Git: 8fe78517f Sun Jul 14 19:43:28 2019 +0200 at lun oct 28 02:10:00 CET 2019