Point Cloud Library (PCL)  1.14.0-dev
ia_fpcs.hpp
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2014-, Open Perception, Inc.
6  * Copyright (C) 2008 Ben Gurion University of the Negev, Beer Sheva, Israel.
7  *
8  * All rights reserved
9  *
10  * Redistribution and use in source and binary forms, with or without
11  * modification, are permitted provided that the following conditions are met
12  *
13  * * Redistributions of source code must retain the above copyright
14  * notice, this list of conditions and the following disclaimer.
15  * * Redistributions in binary form must reproduce the above
16  * copyright notice, this list of conditions and the following
17  * disclaimer in the documentation and/or other materials provided
18  * with the distribution.
19  * * Neither the name of the copyright holder(s) nor the names of its
20  * contributors may be used to endorse or promote products derived
21  * from this software without specific prior written permission.
22  *
23  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
24  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
25  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
26  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
27  * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
28  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
29  * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
30  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
31  * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
32  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
33  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
34  * POSSIBILITY OF SUCH DAMAGE.
35  *
36  */
37 
38 #ifndef PCL_REGISTRATION_IMPL_IA_FPCS_H_
39 #define PCL_REGISTRATION_IMPL_IA_FPCS_H_
40 
41 #include <pcl/common/distances.h>
42 #include <pcl/common/time.h>
43 #include <pcl/common/utils.h>
44 #include <pcl/registration/ia_fpcs.h>
45 #include <pcl/registration/transformation_estimation_3point.h>
46 #include <pcl/sample_consensus/sac_model_plane.h>
47 
48 #include <limits>
49 
50 ///////////////////////////////////////////////////////////////////////////////////////////
51 template <typename PointT>
52 inline float
54  float max_dist,
55  int nr_threads)
56 {
57  const float max_dist_sqr = max_dist * max_dist;
58  const std::size_t s = cloud->size();
59 
61  tree.setInputCloud(cloud);
62 
63  float mean_dist = 0.f;
64  int num = 0;
65  pcl::Indices ids(2);
66  std::vector<float> dists_sqr(2);
67 
68  pcl::utils::ignore(nr_threads);
69 #pragma omp parallel for default(none) shared(tree, cloud) \
70  firstprivate(ids, dists_sqr) reduction(+ : mean_dist, num) \
71  firstprivate(s, max_dist_sqr) num_threads(nr_threads)
72  for (int i = 0; i < 1000; i++) {
73  tree.nearestKSearch((*cloud)[rand() % s], 2, ids, dists_sqr);
74  if (dists_sqr[1] < max_dist_sqr) {
75  mean_dist += std::sqrt(dists_sqr[1]);
76  num++;
77  }
78  }
79 
80  return (mean_dist / num);
81 };
82 
83 ///////////////////////////////////////////////////////////////////////////////////////////
84 template <typename PointT>
85 inline float
87  const pcl::Indices& indices,
88  float max_dist,
89  int nr_threads)
90 {
91  const float max_dist_sqr = max_dist * max_dist;
92  const std::size_t s = indices.size();
93 
95  tree.setInputCloud(cloud);
96 
97  float mean_dist = 0.f;
98  int num = 0;
99  pcl::Indices ids(2);
100  std::vector<float> dists_sqr(2);
101 
102  pcl::utils::ignore(nr_threads);
103 #if OPENMP_LEGACY_CONST_DATA_SHARING_RULE
104 #pragma omp parallel for default(none) shared(tree, cloud, indices) \
105  firstprivate(ids, dists_sqr) reduction(+ : mean_dist, num) num_threads(nr_threads)
106 #else
107 #pragma omp parallel for default(none) shared(tree, cloud, indices, s, max_dist_sqr) \
108  firstprivate(ids, dists_sqr) reduction(+ : mean_dist, num) num_threads(nr_threads)
109 #endif
110  for (int i = 0; i < 1000; i++) {
111  tree.nearestKSearch((*cloud)[indices[rand() % s]], 2, ids, dists_sqr);
112  if (dists_sqr[1] < max_dist_sqr) {
113  mean_dist += std::sqrt(dists_sqr[1]);
114  num++;
115  }
116  }
117 
118  return (mean_dist / num);
119 };
120 
121 ///////////////////////////////////////////////////////////////////////////////////////////
122 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
125 : source_normals_()
126 , target_normals_()
127 , score_threshold_(std::numeric_limits<float>::max())
128 , fitness_score_(std::numeric_limits<float>::max())
129 {
130  reg_name_ = "pcl::registration::FPCSInitialAlignment";
131  max_iterations_ = 0;
132  ransac_iterations_ = 1000;
135 }
136 
137 ///////////////////////////////////////////////////////////////////////////////////////////
138 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
139 void
141  computeTransformation(PointCloudSource& output, const Eigen::Matrix4f& guess)
142 {
143  if (!initCompute())
144  return;
145 
146  final_transformation_ = guess;
147  bool abort = false;
148  std::vector<MatchingCandidates> all_candidates(max_iterations_);
149  pcl::StopWatch timer;
150 
151 #pragma omp parallel default(none) shared(abort, all_candidates, timer) \
152  num_threads(nr_threads_)
153  {
154 #ifdef _OPENMP
155  const unsigned int seed =
156  static_cast<unsigned int>(std::time(nullptr)) ^ omp_get_thread_num();
157  std::srand(seed);
158  PCL_DEBUG("[%s::computeTransformation] Using seed=%u\n", reg_name_.c_str(), seed);
159 #pragma omp for schedule(dynamic)
160 #endif
161  for (int i = 0; i < max_iterations_; i++) {
162 #pragma omp flush(abort)
163 
164  MatchingCandidates candidates(1);
165  pcl::Indices base_indices(4);
166  all_candidates[i] = candidates;
167 
168  if (!abort) {
169  float ratio[2];
170  // select four coplanar point base
171  if (selectBase(base_indices, ratio) == 0) {
172  // calculate candidate pair correspondences using diagonal lengths of base
173  pcl::Correspondences pairs_a, pairs_b;
174  if (bruteForceCorrespondences(base_indices[0], base_indices[1], pairs_a) ==
175  0 &&
176  bruteForceCorrespondences(base_indices[2], base_indices[3], pairs_b) ==
177  0) {
178  // determine candidate matches by combining pair correspondences based on
179  // segment distances
180  std::vector<pcl::Indices> matches;
181  if (determineBaseMatches(base_indices, matches, pairs_a, pairs_b, ratio) ==
182  0) {
183  // check and evaluate candidate matches and store them
184  handleMatches(base_indices, matches, candidates);
185  if (!candidates.empty())
186  all_candidates[i] = candidates;
187  }
188  }
189  }
190 
191  // check terminate early (time or fitness_score threshold reached)
192  abort = (!candidates.empty() ? candidates[0].fitness_score < score_threshold_
193  : abort);
194  abort = (abort ? abort : timer.getTimeSeconds() > max_runtime_);
195 
196 #pragma omp flush(abort)
197  }
198  }
199  }
200 
201  // determine best match over all tries
202  finalCompute(all_candidates);
203 
204  // apply the final transformation
205  pcl::transformPointCloud(*input_, output, final_transformation_);
206 
207  deinitCompute();
208 }
209 
210 ///////////////////////////////////////////////////////////////////////////////////////////
211 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
212 bool
214  initCompute()
215 {
216  const unsigned int seed = std::time(nullptr);
217  std::srand(seed);
218  PCL_DEBUG("[%s::initCompute] Using seed=%u\n", reg_name_.c_str(), seed);
219 
220  // basic pcl initialization
222  return (false);
223 
224  // check if source and target are given
225  if (!input_ || !target_) {
226  PCL_ERROR("[%s::initCompute] Source or target dataset not given!\n",
227  reg_name_.c_str());
228  return (false);
229  }
230 
231  if (!target_indices_ || target_indices_->empty()) {
232  target_indices_.reset(new pcl::Indices(target_->size()));
233  int index = 0;
234  for (auto& target_index : *target_indices_)
235  target_index = index++;
236  target_cloud_updated_ = true;
237  }
238 
239  // if a sample size for the point clouds is given; preferably no sampling of target
240  // cloud
241  if (nr_samples_ != 0) {
242  const int ss = static_cast<int>(indices_->size());
243  const int sample_fraction_src = std::max(1, static_cast<int>(ss / nr_samples_));
244 
245  source_indices_ = pcl::IndicesPtr(new pcl::Indices);
246  for (int i = 0; i < ss; i++)
247  if (rand() % sample_fraction_src == 0)
248  source_indices_->push_back((*indices_)[i]);
249  }
250  else
251  source_indices_ = indices_;
252 
253  // check usage of normals
254  if (source_normals_ && target_normals_ && source_normals_->size() == input_->size() &&
255  target_normals_->size() == target_->size())
256  use_normals_ = true;
257 
258  // set up tree structures
259  if (target_cloud_updated_) {
260  tree_->setInputCloud(target_, target_indices_);
261  target_cloud_updated_ = false;
262  }
263 
264  // set predefined variables
265  constexpr int min_iterations = 4;
266  constexpr float diameter_fraction = 0.3f;
267 
268  // get diameter of input cloud (distance between farthest points)
269  Eigen::Vector4f pt_min, pt_max;
270  pcl::getMinMax3D(*target_, *target_indices_, pt_min, pt_max);
271  diameter_ = (pt_max - pt_min).norm();
272 
273  // derive the limits for the random base selection
274  float max_base_diameter = diameter_ * approx_overlap_ * 2.f;
275  max_base_diameter_sqr_ = max_base_diameter * max_base_diameter;
276 
277  // normalize the delta
278  if (normalize_delta_) {
279  float mean_dist = getMeanPointDensity<PointTarget>(
280  target_, *target_indices_, 0.05f * diameter_, nr_threads_);
281  delta_ *= mean_dist;
282  }
283 
284  // heuristic determination of number of trials to have high probability of finding a
285  // good solution
286  if (max_iterations_ == 0) {
287  float first_est = std::log(small_error_) /
288  std::log(1.0 - std::pow(static_cast<double>(approx_overlap_),
289  static_cast<double>(min_iterations)));
290  max_iterations_ =
291  static_cast<int>(first_est / (diameter_fraction * approx_overlap_ * 2.f));
292  }
293 
294  // set further parameter
295  if (score_threshold_ == std::numeric_limits<float>::max())
296  score_threshold_ = 1.f - approx_overlap_;
297 
298  if (max_iterations_ < 4)
299  max_iterations_ = 4;
300 
301  if (max_runtime_ < 1)
302  max_runtime_ = std::numeric_limits<int>::max();
303 
304  // calculate internal parameters based on the the estimated point density
305  max_pair_diff_ = delta_ * 2.f;
306  max_edge_diff_ = delta_ * 4.f;
307  coincidation_limit_ = delta_ * 2.f; // EDITED: originally std::sqrt (delta_ * 2.f)
308  max_mse_ = powf(delta_ * 2.f, 2.f);
309  max_inlier_dist_sqr_ = powf(delta_ * 2.f, 2.f);
310 
311  // reset fitness_score
312  fitness_score_ = std::numeric_limits<float>::max();
313 
314  return (true);
315 }
316 
317 ///////////////////////////////////////////////////////////////////////////////////////////
318 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
319 int
321  selectBase(pcl::Indices& base_indices, float (&ratio)[2])
322 {
323  const float too_close_sqr = max_base_diameter_sqr_ * 0.01;
324 
325  Eigen::VectorXf coefficients(4);
327  plane.setIndices(target_indices_);
328  Eigen::Vector4f centre_pt;
329  float nearest_to_plane = std::numeric_limits<float>::max();
330 
331  // repeat base search until valid quadruple was found or ransac_iterations_ number of
332  // tries were unsuccessful
333  for (int i = 0; i < ransac_iterations_; i++) {
334  // random select an appropriate point triple
335  if (selectBaseTriangle(base_indices) < 0)
336  continue;
337 
338  pcl::Indices base_triple(base_indices.begin(), base_indices.end() - 1);
339  plane.computeModelCoefficients(base_triple, coefficients);
340  pcl::compute3DCentroid(*target_, base_triple, centre_pt);
341 
342  // loop over all points in source cloud to find most suitable fourth point
343  const PointTarget* pt1 = &((*target_)[base_indices[0]]);
344  const PointTarget* pt2 = &((*target_)[base_indices[1]]);
345  const PointTarget* pt3 = &((*target_)[base_indices[2]]);
346 
347  for (const auto& target_index : *target_indices_) {
348  const PointTarget* pt4 = &((*target_)[target_index]);
349 
350  float d1 = pcl::squaredEuclideanDistance(*pt4, *pt1);
351  float d2 = pcl::squaredEuclideanDistance(*pt4, *pt2);
352  float d3 = pcl::squaredEuclideanDistance(*pt4, *pt3);
353  float d4 = (pt4->getVector3fMap() - centre_pt.head<3>()).squaredNorm();
354 
355  // check distance between points w.r.t minimum sampling distance; EDITED -> 4th
356  // point now also limited by max base line
357  if (d1 < too_close_sqr || d2 < too_close_sqr || d3 < too_close_sqr ||
358  d4 < too_close_sqr || d1 > max_base_diameter_sqr_ ||
359  d2 > max_base_diameter_sqr_ || d3 > max_base_diameter_sqr_)
360  continue;
361 
362  // check distance to plane to get point closest to plane
363  float dist_to_plane = pcl::pointToPlaneDistance(*pt4, coefficients);
364  if (dist_to_plane < nearest_to_plane) {
365  base_indices[3] = target_index;
366  nearest_to_plane = dist_to_plane;
367  }
368  }
369 
370  // check if at least one point fulfilled the conditions
371  if (nearest_to_plane != std::numeric_limits<float>::max()) {
372  // order points to build largest quadrangle and calculate intersection ratios of
373  // diagonals
374  setupBase(base_indices, ratio);
375  return (0);
376  }
377  }
378 
379  // return unsuccessful if no quadruple was selected
380  return (-1);
381 }
382 
383 ///////////////////////////////////////////////////////////////////////////////////////////
384 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
385 int
387  selectBaseTriangle(pcl::Indices& base_indices)
388 {
389  const auto nr_points = target_indices_->size();
390  float best_t = 0.f;
391 
392  // choose random first point
393  base_indices[0] = (*target_indices_)[rand() % nr_points];
394  auto* index1 = base_indices.data();
395 
396  // random search for 2 other points (as far away as overlap allows)
397  for (int i = 0; i < ransac_iterations_; i++) {
398  auto* index2 = &(*target_indices_)[rand() % nr_points];
399  auto* index3 = &(*target_indices_)[rand() % nr_points];
400 
401  Eigen::Vector3f u =
402  (*target_)[*index2].getVector3fMap() - (*target_)[*index1].getVector3fMap();
403  Eigen::Vector3f v =
404  (*target_)[*index3].getVector3fMap() - (*target_)[*index1].getVector3fMap();
405  float t =
406  u.cross(v).squaredNorm(); // triangle area (0.5 * sqrt(t)) should be maximal
407 
408  // check for most suitable point triple
409  if (t > best_t && u.squaredNorm() < max_base_diameter_sqr_ &&
410  v.squaredNorm() < max_base_diameter_sqr_) {
411  best_t = t;
412  base_indices[1] = *index2;
413  base_indices[2] = *index3;
414  }
415  }
416 
417  // return if a triplet could be selected
418  return (best_t == 0.f ? -1 : 0);
419 }
420 
421 ///////////////////////////////////////////////////////////////////////////////////////////
422 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
423 void
425  setupBase(pcl::Indices& base_indices, float (&ratio)[2])
426 {
427  float best_t = std::numeric_limits<float>::max();
428  const pcl::Indices copy(base_indices.begin(), base_indices.end());
429  pcl::Indices temp(base_indices.begin(), base_indices.end());
430 
431  // loop over all combinations of base points
432  for (auto i = copy.begin(), i_e = copy.end(); i != i_e; ++i)
433  for (auto j = copy.begin(), j_e = copy.end(); j != j_e; ++j) {
434  if (i == j)
435  continue;
436 
437  for (auto k = copy.begin(), k_e = copy.end(); k != k_e; ++k) {
438  if (k == j || k == i)
439  continue;
440 
441  auto l = copy.begin();
442  while (l == i || l == j || l == k)
443  ++l;
444 
445  temp[0] = *i;
446  temp[1] = *j;
447  temp[2] = *k;
448  temp[3] = *l;
449 
450  // calculate diagonal intersection ratios and check for suitable segment to
451  // segment distances
452  float ratio_temp[2];
453  float t = segmentToSegmentDist(temp, ratio_temp);
454  if (t < best_t) {
455  best_t = t;
456  ratio[0] = ratio_temp[0];
457  ratio[1] = ratio_temp[1];
458  base_indices = temp;
459  }
460  }
461  }
462 }
463 
464 ///////////////////////////////////////////////////////////////////////////////////////////
465 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
466 float
468  segmentToSegmentDist(const pcl::Indices& base_indices, float (&ratio)[2])
469 {
470  // get point vectors
471  Eigen::Vector3f u = (*target_)[base_indices[1]].getVector3fMap() -
472  (*target_)[base_indices[0]].getVector3fMap();
473  Eigen::Vector3f v = (*target_)[base_indices[3]].getVector3fMap() -
474  (*target_)[base_indices[2]].getVector3fMap();
475  Eigen::Vector3f w = (*target_)[base_indices[0]].getVector3fMap() -
476  (*target_)[base_indices[2]].getVector3fMap();
477 
478  // calculate segment distances
479  float a = u.dot(u);
480  float b = u.dot(v);
481  float c = v.dot(v);
482  float d = u.dot(w);
483  float e = v.dot(w);
484  float D = a * c - b * b;
485  float sN = 0.f, sD = D;
486  float tN = 0.f, tD = D;
487 
488  // check segments
489  if (D < small_error_) {
490  sN = 0.f;
491  sD = 1.f;
492  tN = e;
493  tD = c;
494  }
495  else {
496  sN = (b * e - c * d);
497  tN = (a * e - b * d);
498 
499  if (sN < 0.f) {
500  sN = 0.f;
501  tN = e;
502  tD = c;
503  }
504  else if (sN > sD) {
505  sN = sD;
506  tN = e + b;
507  tD = c;
508  }
509  }
510 
511  if (tN < 0.f) {
512  tN = 0.f;
513 
514  if (-d < 0.f)
515  sN = 0.f;
516 
517  else if (-d > a)
518  sN = sD;
519 
520  else {
521  sN = -d;
522  sD = a;
523  }
524  }
525 
526  else if (tN > tD) {
527  tN = tD;
528 
529  if ((-d + b) < 0.f)
530  sN = 0.f;
531 
532  else if ((-d + b) > a)
533  sN = sD;
534 
535  else {
536  sN = (-d + b);
537  sD = a;
538  }
539  }
540 
541  // set intersection ratios
542  ratio[0] = (std::abs(sN) < small_error_) ? 0.f : sN / sD;
543  ratio[1] = (std::abs(tN) < small_error_) ? 0.f : tN / tD;
544 
545  Eigen::Vector3f x = w + (ratio[0] * u) - (ratio[1] * v);
546  return (x.norm());
547 }
548 
549 ///////////////////////////////////////////////////////////////////////////////////////////
550 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
551 int
553  bruteForceCorrespondences(int idx1, int idx2, pcl::Correspondences& pairs)
554 {
555  const float max_norm_diff = 0.5f * max_norm_diff_ * M_PI / 180.f;
556 
557  // calculate reference segment distance and normal angle
558  float ref_dist = pcl::euclideanDistance((*target_)[idx1], (*target_)[idx2]);
559  float ref_norm_angle =
560  (use_normals_ ? ((*target_normals_)[idx1].getNormalVector3fMap() -
561  (*target_normals_)[idx2].getNormalVector3fMap())
562  .norm()
563  : 0.f);
564 
565  // loop over all pairs of points in source point cloud
566  auto it_out = source_indices_->begin(), it_out_e = source_indices_->end() - 1;
567  auto it_in_e = source_indices_->end();
568  for (; it_out != it_out_e; it_out++) {
569  auto it_in = it_out + 1;
570  const PointSource* pt1 = &(*input_)[*it_out];
571  for (; it_in != it_in_e; it_in++) {
572  const PointSource* pt2 = &(*input_)[*it_in];
573 
574  // check point distance compared to reference dist (from base)
575  float dist = pcl::euclideanDistance(*pt1, *pt2);
576  if (std::abs(dist - ref_dist) < max_pair_diff_) {
577  // add here normal evaluation if normals are given
578  if (use_normals_) {
579  const NormalT* pt1_n = &((*source_normals_)[*it_out]);
580  const NormalT* pt2_n = &((*source_normals_)[*it_in]);
581 
582  float norm_angle_1 =
583  (pt1_n->getNormalVector3fMap() - pt2_n->getNormalVector3fMap()).norm();
584  float norm_angle_2 =
585  (pt1_n->getNormalVector3fMap() + pt2_n->getNormalVector3fMap()).norm();
586 
587  float norm_diff = std::min<float>(std::abs(norm_angle_1 - ref_norm_angle),
588  std::abs(norm_angle_2 - ref_norm_angle));
589  if (norm_diff > max_norm_diff)
590  continue;
591  }
592 
593  pairs.push_back(pcl::Correspondence(*it_in, *it_out, dist));
594  pairs.push_back(pcl::Correspondence(*it_out, *it_in, dist));
595  }
596  }
597  }
598 
599  // return success if at least one correspondence was found
600  return (pairs.empty() ? -1 : 0);
601 }
602 
603 ///////////////////////////////////////////////////////////////////////////////////////////
604 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
605 int
607  determineBaseMatches(const pcl::Indices& base_indices,
608  std::vector<pcl::Indices>& matches,
609  const pcl::Correspondences& pairs_a,
610  const pcl::Correspondences& pairs_b,
611  const float (&ratio)[2])
612 {
613  // calculate edge lengths of base
614  float dist_base[4];
615  dist_base[0] =
616  pcl::euclideanDistance((*target_)[base_indices[0]], (*target_)[base_indices[2]]);
617  dist_base[1] =
618  pcl::euclideanDistance((*target_)[base_indices[0]], (*target_)[base_indices[3]]);
619  dist_base[2] =
620  pcl::euclideanDistance((*target_)[base_indices[1]], (*target_)[base_indices[2]]);
621  dist_base[3] =
622  pcl::euclideanDistance((*target_)[base_indices[1]], (*target_)[base_indices[3]]);
623 
624  // loop over first point pair correspondences and store intermediate points 'e' in new
625  // point cloud
626  PointCloudSourcePtr cloud_e(new PointCloudSource);
627  cloud_e->resize(pairs_a.size() * 2);
628  auto it_pt = cloud_e->begin();
629  for (const auto& pair : pairs_a) {
630  const PointSource* pt1 = &((*input_)[pair.index_match]);
631  const PointSource* pt2 = &((*input_)[pair.index_query]);
632 
633  // calculate intermediate points using both ratios from base (r1,r2)
634  for (int i = 0; i < 2; i++, it_pt++) {
635  it_pt->x = pt1->x + ratio[i] * (pt2->x - pt1->x);
636  it_pt->y = pt1->y + ratio[i] * (pt2->y - pt1->y);
637  it_pt->z = pt1->z + ratio[i] * (pt2->z - pt1->z);
638  }
639  }
640 
641  // initialize new kd tree of intermediate points from first point pair correspondences
643  tree_e->setInputCloud(cloud_e);
644 
645  pcl::Indices ids;
646  std::vector<float> dists_sqr;
647 
648  // loop over second point pair correspondences
649  for (const auto& pair : pairs_b) {
650  const PointTarget* pt1 = &((*input_)[pair.index_match]);
651  const PointTarget* pt2 = &((*input_)[pair.index_query]);
652 
653  // calculate intermediate points using both ratios from base (r1,r2)
654  for (const float& r : ratio) {
655  PointTarget pt_e;
656  pt_e.x = pt1->x + r * (pt2->x - pt1->x);
657  pt_e.y = pt1->y + r * (pt2->y - pt1->y);
658  pt_e.z = pt1->z + r * (pt2->z - pt1->z);
659 
660  // search for corresponding intermediate points
661  tree_e->radiusSearch(pt_e, coincidation_limit_, ids, dists_sqr);
662  for (const auto& id : ids) {
663  pcl::Indices match_indices(4);
664 
665  match_indices[0] =
666  pairs_a[static_cast<int>(std::floor((id / 2.f)))].index_match;
667  match_indices[1] =
668  pairs_a[static_cast<int>(std::floor((id / 2.f)))].index_query;
669  match_indices[2] = pair.index_match;
670  match_indices[3] = pair.index_query;
671 
672  // EDITED: added coarse check of match based on edge length (due to rigid-body )
673  if (checkBaseMatch(match_indices, dist_base) < 0)
674  continue;
675 
676  matches.push_back(match_indices);
677  }
678  }
679  }
680 
681  // return unsuccessful if no match was found
682  return (!matches.empty() ? 0 : -1);
683 }
684 
685 ///////////////////////////////////////////////////////////////////////////////////////////
686 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
687 int
689  checkBaseMatch(const pcl::Indices& match_indices, const float (&dist_ref)[4])
690 {
691  float d0 =
692  pcl::euclideanDistance((*input_)[match_indices[0]], (*input_)[match_indices[2]]);
693  float d1 =
694  pcl::euclideanDistance((*input_)[match_indices[0]], (*input_)[match_indices[3]]);
695  float d2 =
696  pcl::euclideanDistance((*input_)[match_indices[1]], (*input_)[match_indices[2]]);
697  float d3 =
698  pcl::euclideanDistance((*input_)[match_indices[1]], (*input_)[match_indices[3]]);
699 
700  // check edge distances of match w.r.t the base
701  return (std::abs(d0 - dist_ref[0]) < max_edge_diff_ &&
702  std::abs(d1 - dist_ref[1]) < max_edge_diff_ &&
703  std::abs(d2 - dist_ref[2]) < max_edge_diff_ &&
704  std::abs(d3 - dist_ref[3]) < max_edge_diff_)
705  ? 0
706  : -1;
707 }
708 
709 ///////////////////////////////////////////////////////////////////////////////////////////
710 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
711 void
713  handleMatches(const pcl::Indices& base_indices,
714  std::vector<pcl::Indices>& matches,
715  MatchingCandidates& candidates)
716 {
717  candidates.resize(1);
718  float fitness_score = std::numeric_limits<float>::max();
719 
720  // loop over all Candidate matches
721  for (auto& match : matches) {
722  Eigen::Matrix4f transformation_temp;
723  pcl::Correspondences correspondences_temp;
724 
725  // determine correspondences between base and match according to their distance to
726  // centroid
727  linkMatchWithBase(base_indices, match, correspondences_temp);
728 
729  // check match based on residuals of the corresponding points after
730  if (validateMatch(base_indices, match, correspondences_temp, transformation_temp) <
731  0)
732  continue;
733 
734  // check resulting using a sub sample of the source point cloud and compare to
735  // previous matches
736  if (validateTransformation(transformation_temp, fitness_score) < 0)
737  continue;
738 
739  // store best match as well as associated fitness_score and transformation
740  candidates[0].fitness_score = fitness_score;
741  candidates[0].transformation = transformation_temp;
742  correspondences_temp.erase(correspondences_temp.end() - 1);
743  candidates[0].correspondences = correspondences_temp;
744  }
745 }
746 
747 ///////////////////////////////////////////////////////////////////////////////////////////
748 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
749 void
751  linkMatchWithBase(const pcl::Indices& base_indices,
752  pcl::Indices& match_indices,
753  pcl::Correspondences& correspondences)
754 {
755  // calculate centroid of base and target
756  Eigen::Vector4f centre_base{0, 0, 0, 0}, centre_match{0, 0, 0, 0};
757  pcl::compute3DCentroid(*target_, base_indices, centre_base);
758  pcl::compute3DCentroid(*input_, match_indices, centre_match);
759 
760  PointTarget centre_pt_base;
761  centre_pt_base.x = centre_base[0];
762  centre_pt_base.y = centre_base[1];
763  centre_pt_base.z = centre_base[2];
764 
765  PointSource centre_pt_match;
766  centre_pt_match.x = centre_match[0];
767  centre_pt_match.y = centre_match[1];
768  centre_pt_match.z = centre_match[2];
769 
770  // find corresponding points according to their distance to the centroid
771  pcl::Indices copy = match_indices;
772 
773  auto it_match_orig = match_indices.begin();
774  for (auto it_base = base_indices.cbegin(), it_base_e = base_indices.cend();
775  it_base != it_base_e;
776  it_base++, it_match_orig++) {
777  float dist_sqr_1 =
778  pcl::squaredEuclideanDistance((*target_)[*it_base], centre_pt_base);
779  float best_diff_sqr = std::numeric_limits<float>::max();
780  int best_index = -1;
781 
782  for (const auto& match_index : copy) {
783  // calculate difference of distances to centre point
784  float dist_sqr_2 =
785  pcl::squaredEuclideanDistance((*input_)[match_index], centre_pt_match);
786  float diff_sqr = std::abs(dist_sqr_1 - dist_sqr_2);
787 
788  if (diff_sqr < best_diff_sqr) {
789  best_diff_sqr = diff_sqr;
790  best_index = match_index;
791  }
792  }
793 
794  // assign new correspondence and update indices of matched targets
795  correspondences.push_back(pcl::Correspondence(best_index, *it_base, best_diff_sqr));
796  *it_match_orig = best_index;
797  }
798 }
799 
800 ///////////////////////////////////////////////////////////////////////////////////////////
801 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
802 int
804  validateMatch(const pcl::Indices& base_indices,
805  const pcl::Indices& match_indices,
806  const pcl::Correspondences& correspondences,
807  Eigen::Matrix4f& transformation)
808 {
809  // only use triplet of points to simplify process (possible due to planar case)
810  pcl::Correspondences correspondences_temp = correspondences;
811  correspondences_temp.erase(correspondences_temp.end() - 1);
812 
813  // estimate transformation between correspondence set
814  transformation_estimation_->estimateRigidTransformation(
815  *input_, *target_, correspondences_temp, transformation);
816 
817  // transform base points
818  PointCloudSource match_transformed;
819  pcl::transformPointCloud(*input_, match_indices, match_transformed, transformation);
820 
821  // calculate residuals of transformation and check against maximum threshold
822  std::size_t nr_points = correspondences_temp.size();
823  float mse = 0.f;
824  for (std::size_t i = 0; i < nr_points; i++)
825  mse += pcl::squaredEuclideanDistance(match_transformed.points[i],
826  target_->points[base_indices[i]]);
827 
828  mse /= nr_points;
829  return (mse < max_mse_ ? 0 : -1);
830 }
831 
832 ///////////////////////////////////////////////////////////////////////////////////////////
833 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
834 int
836  validateTransformation(Eigen::Matrix4f& transformation, float& fitness_score)
837 {
838  // transform source point cloud
839  PointCloudSource source_transformed;
841  *input_, *source_indices_, source_transformed, transformation);
842 
843  std::size_t nr_points = source_transformed.size();
844  std::size_t terminate_value =
845  fitness_score > 1 ? 0
846  : static_cast<std::size_t>((1.f - fitness_score) * nr_points);
847 
848  float inlier_score_temp = 0;
849  pcl::Indices ids;
850  std::vector<float> dists_sqr;
851  auto it = source_transformed.begin();
852 
853  for (std::size_t i = 0; i < nr_points; it++, i++) {
854  // search for nearest point using kd tree search
855  tree_->nearestKSearch(*it, 1, ids, dists_sqr);
856  inlier_score_temp += (dists_sqr[0] < max_inlier_dist_sqr_ ? 1 : 0);
857 
858  // early terminating
859  if (nr_points - i + inlier_score_temp < terminate_value)
860  break;
861  }
862 
863  // check current costs and return unsuccessful if larger than previous ones
864  inlier_score_temp /= static_cast<float>(nr_points);
865  float fitness_score_temp = 1.f - inlier_score_temp;
866 
867  if (fitness_score_temp > fitness_score)
868  return (-1);
869 
870  fitness_score = fitness_score_temp;
871  return (0);
872 }
873 
874 ///////////////////////////////////////////////////////////////////////////////////////////
875 template <typename PointSource, typename PointTarget, typename NormalT, typename Scalar>
876 void
878  finalCompute(const std::vector<MatchingCandidates>& candidates)
879 {
880  // get best fitness_score over all tries
881  int nr_candidates = static_cast<int>(candidates.size());
882  int best_index = -1;
883  float best_score = std::numeric_limits<float>::max();
884  for (int i = 0; i < nr_candidates; i++) {
885  const float& fitness_score = candidates[i][0].fitness_score;
886  if (fitness_score < best_score) {
887  best_score = fitness_score;
888  best_index = i;
889  }
890  }
891 
892  // check if a valid candidate was available
893  if (!(best_index < 0)) {
894  fitness_score_ = candidates[best_index][0].fitness_score;
895  final_transformation_ = candidates[best_index][0].transformation;
896  *correspondences_ = candidates[best_index][0].correspondences;
897 
898  // here we define convergence if resulting fitness_score is below 1-threshold
899  converged_ = fitness_score_ < score_threshold_;
900  }
901 }
902 
903 ///////////////////////////////////////////////////////////////////////////////////////////
904 
905 #endif // PCL_REGISTRATION_IMPL_IA_4PCS_H_
PCL base class.
Definition: pcl_base.h:70
std::size_t size() const
Definition: point_cloud.h:443
shared_ptr< const PointCloud< PointT > > ConstPtr
Definition: point_cloud.h:414
typename KdTreeReciprocal::Ptr KdTreeReciprocalPtr
Definition: registration.h:74
std::string reg_name_
The registration method name.
Definition: registration.h:548
typename PointCloudSource::Ptr PointCloudSourcePtr
Definition: registration.h:77
int ransac_iterations_
The number of iterations RANSAC should run for.
Definition: registration.h:566
TransformationEstimationPtr transformation_estimation_
A TransformationEstimation object, used to calculate the 4x4 rigid transformation.
Definition: registration.h:625
int max_iterations_
The maximum number of iterations the internal optimization should run for.
Definition: registration.h:563
void setIndices(const IndicesPtr &indices)
Provide a pointer to the vector of indices that represents the input data.
Definition: sac_model.h:324
SampleConsensusModelPlane defines a model for 3D plane segmentation.
bool computeModelCoefficients(const Indices &samples, Eigen::VectorXf &model_coefficients) const override
Check whether the given index samples can form a valid plane model, compute the model coefficients fr...
Simple stopwatch.
Definition: time.h:59
double getTimeSeconds() const
Retrieve the time in seconds spent since the last call to reset().
Definition: time.h:74
virtual void finalCompute(const std::vector< MatchingCandidates > &candidates)
Final computation of best match out of vector of best matches.
Definition: ia_fpcs.hpp:878
void setupBase(pcl::Indices &base_indices, float(&ratio)[2])
Setup the base (four coplanar points) by ordering the points and computing intersection ratios and se...
Definition: ia_fpcs.hpp:425
int selectBaseTriangle(pcl::Indices &base_indices)
Select randomly a triplet of points with large point-to-point distances.
Definition: ia_fpcs.hpp:387
virtual int bruteForceCorrespondences(int idx1, int idx2, pcl::Correspondences &pairs)
Search for corresponding point pairs given the distance between two base points.
Definition: ia_fpcs.hpp:553
int selectBase(pcl::Indices &base_indices, float(&ratio)[2])
Select an approximately coplanar set of four points from the source cloud.
Definition: ia_fpcs.hpp:321
virtual int validateTransformation(Eigen::Matrix4f &transformation, float &fitness_score)
Validate the transformation by calculating the number of inliers after transforming the source cloud.
Definition: ia_fpcs.hpp:836
virtual int determineBaseMatches(const pcl::Indices &base_indices, std::vector< pcl::Indices > &matches, const pcl::Correspondences &pairs_a, const pcl::Correspondences &pairs_b, const float(&ratio)[2])
Determine base matches by combining the point pair candidate and search for coinciding intersection p...
Definition: ia_fpcs.hpp:607
virtual void linkMatchWithBase(const pcl::Indices &base_indices, pcl::Indices &match_indices, pcl::Correspondences &correspondences)
Sets the correspondences between the base B and the match M by using the distance of each point to th...
Definition: ia_fpcs.hpp:751
virtual void handleMatches(const pcl::Indices &base_indices, std::vector< pcl::Indices > &matches, MatchingCandidates &candidates)
Method to handle current candidate matches.
Definition: ia_fpcs.hpp:713
int checkBaseMatch(const pcl::Indices &match_indices, const float(&ds)[4])
Check if outer rectangle distance of matched points fit with the base rectangle.
Definition: ia_fpcs.hpp:689
float segmentToSegmentDist(const pcl::Indices &base_indices, float(&ratio)[2])
Calculate intersection ratios and segment to segment distances of base diagonals.
Definition: ia_fpcs.hpp:468
virtual int validateMatch(const pcl::Indices &base_indices, const pcl::Indices &match_indices, const pcl::Correspondences &correspondences, Eigen::Matrix4f &transformation)
Validate the matching by computing the transformation between the source and target based on the four...
Definition: ia_fpcs.hpp:804
void computeTransformation(PointCloudSource &output, const Eigen::Matrix4f &guess) override
Rigid transformation computation method.
Definition: ia_fpcs.hpp:141
virtual bool initCompute()
Internal computation initialization.
Definition: ia_fpcs.hpp:214
TransformationEstimation3Points represents the class for transformation estimation based on:
search::KdTree is a wrapper class which inherits the pcl::KdTree class for performing search function...
Definition: kdtree.h:62
int nearestKSearch(const PointT &point, int k, Indices &k_indices, std::vector< float > &k_sqr_distances) const override
Search for the k-nearest neighbors for the given query point.
Definition: kdtree.hpp:88
bool setInputCloud(const PointCloudConstPtr &cloud, const IndicesConstPtr &indices=IndicesConstPtr()) override
Provide a pointer to the input dataset.
Definition: kdtree.hpp:76
Define standard C methods to do distance calculations.
Define methods for measuring time spent in code blocks.
void getMinMax3D(const pcl::PointCloud< PointT > &cloud, PointT &min_pt, PointT &max_pt)
Get the minimum and maximum values on each of the 3 (x-y-z) dimensions in a given pointcloud.
Definition: common.hpp:295
void transformPointCloud(const pcl::PointCloud< PointT > &cloud_in, pcl::PointCloud< PointT > &cloud_out, const Eigen::Matrix< Scalar, 4, 4 > &transform, bool copy_all_fields)
Apply a rigid transform defined by a 4x4 matrix.
Definition: transforms.hpp:221
unsigned int compute3DCentroid(ConstCloudIterator< PointT > &cloud_iterator, Eigen::Matrix< Scalar, 4, 1 > &centroid)
Compute the 3D (X-Y-Z) centroid of a set of points and return it as a 3D vector.
Definition: centroid.hpp:57
double pointToPlaneDistance(const Point &p, double a, double b, double c, double d)
Get the distance from a point to a plane (unsigned) defined by ax+by+cz+d=0.
std::vector< MatchingCandidate, Eigen::aligned_allocator< MatchingCandidate > > MatchingCandidates
void ignore(const T &...)
Utility function to eliminate unused variable warnings.
Definition: utils.h:62
float squaredEuclideanDistance(const PointType1 &p1, const PointType2 &p2)
Calculate the squared euclidean distance between the two given points.
Definition: distances.h:182
float getMeanPointDensity(const typename pcl::PointCloud< PointT >::ConstPtr &cloud, float max_dist, int nr_threads=1)
Compute the mean point density of a given point cloud.
Definition: ia_fpcs.hpp:53
std::vector< pcl::Correspondence, Eigen::aligned_allocator< pcl::Correspondence > > Correspondences
float euclideanDistance(const PointType1 &p1, const PointType2 &p2)
Calculate the euclidean distance between the two given points.
Definition: distances.h:204
IndicesAllocator<> Indices
Type used for indices in PCL.
Definition: types.h:133
shared_ptr< Indices > IndicesPtr
Definition: pcl_base.h:58
#define M_PI
Definition: pcl_macros.h:201
Correspondence represents a match between two entities (e.g., points, descriptors,...
A point structure representing normal coordinates and the surface curvature estimate.