Point Cloud Library (PCL)  1.14.1-dev
gicp.hpp
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2010, Willow Garage, Inc.
6  * Copyright (c) 2012-, Open Perception, Inc.
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
12  * are met:
13  *
14  * * Redistributions of source code must retain the above copyright
15  * notice, this list of conditions and the following disclaimer.
16  * * Redistributions in binary form must reproduce the above
17  * copyright notice, this list of conditions and the following
18  * disclaimer in the documentation and/or other materials provided
19  * with the distribution.
20  * * Neither the name of the copyright holder(s) nor the names of its
21  * contributors may be used to endorse or promote products derived
22  * from this software without specific prior written permission.
23  *
24  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
28  * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
29  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
30  * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
31  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
32  * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
33  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
34  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
35  * POSSIBILITY OF SUCH DAMAGE.
36  *
37  * $Id$
38  *
39  */
40 
41 #ifndef PCL_REGISTRATION_IMPL_GICP_HPP_
42 #define PCL_REGISTRATION_IMPL_GICP_HPP_
43 
44 #include <pcl/registration/exceptions.h>
45 
46 namespace pcl {
47 
48 template <typename PointSource, typename PointTarget, typename Scalar>
49 void
51  unsigned int nr_threads)
52 {
53 #ifdef _OPENMP
54  if (nr_threads == 0)
55  threads_ = omp_get_num_procs();
56  else
57  threads_ = nr_threads;
58  PCL_DEBUG("[pcl::GeneralizedIterativeClosestPoint::setNumberOfThreads] Setting "
59  "number of threads to %u.\n",
60  threads_);
61 #else
62  threads_ = 1;
63  if (nr_threads != 1)
64  PCL_WARN("[pcl::GeneralizedIterativeClosestPoint::setNumberOfThreads] "
65  "Parallelization is requested, but OpenMP is not available! Continuing "
66  "without parallelization.\n");
67 #endif // _OPENMP
68 }
69 
70 template <typename PointSource, typename PointTarget, typename Scalar>
71 template <typename PointT>
72 void
74  typename pcl::PointCloud<PointT>::ConstPtr cloud,
75  const typename pcl::search::KdTree<PointT>::Ptr kdtree,
76  MatricesVector& cloud_covariances)
77 {
78  if (k_correspondences_ > static_cast<int>(cloud->size())) {
79  PCL_ERROR("[pcl::GeneralizedIterativeClosestPoint::computeCovariances] Number or "
80  "points in cloud (%lu) is less than k_correspondences_ (%lu)!\n",
81  cloud->size(),
82  k_correspondences_);
83  return;
84  }
85 
86  Eigen::Vector3d mean;
87  Eigen::Matrix3d cov;
88  pcl::Indices nn_indices(k_correspondences_);
89  std::vector<float> nn_dist_sq(k_correspondences_);
90 
91  // We should never get there but who knows
92  if (cloud_covariances.size() < cloud->size())
93  cloud_covariances.resize(cloud->size());
94 
95 #pragma omp parallel for default(none) num_threads(threads_) schedule(dynamic, 32) \
96  shared(cloud, cloud_covariances, kdtree) \
97  firstprivate(mean, cov, nn_indices, nn_dist_sq)
98  for (std::ptrdiff_t i = 0; i < static_cast<std::ptrdiff_t>(cloud->size()); ++i) {
99  const PointT& query_point = (*cloud)[i];
100  // Zero out the cov and mean
101  cov.setZero();
102  mean.setZero();
103 
104  // Search for the K nearest neighbours
105  kdtree->nearestKSearch(query_point, k_correspondences_, nn_indices, nn_dist_sq);
106 
107  // Find the covariance matrix
108  for (int j = 0; j < k_correspondences_; j++) {
109  // de-mean neighbourhood to avoid inaccuracies when far away from origin
110  const double ptx = (*cloud)[nn_indices[j]].x - query_point.x,
111  pty = (*cloud)[nn_indices[j]].y - query_point.y,
112  ptz = (*cloud)[nn_indices[j]].z - query_point.z;
113 
114  mean[0] += ptx;
115  mean[1] += pty;
116  mean[2] += ptz;
117 
118  cov(0, 0) += ptx * ptx;
119 
120  cov(1, 0) += pty * ptx;
121  cov(1, 1) += pty * pty;
122 
123  cov(2, 0) += ptz * ptx;
124  cov(2, 1) += ptz * pty;
125  cov(2, 2) += ptz * ptz;
126  }
127 
128  mean /= static_cast<double>(k_correspondences_);
129  // Get the actual covariance
130  for (int k = 0; k < 3; k++)
131  for (int l = 0; l <= k; l++) {
132  cov(k, l) /= static_cast<double>(k_correspondences_);
133  cov(k, l) -= mean[k] * mean[l];
134  cov(l, k) = cov(k, l);
135  }
136 
137  // Compute the SVD (covariance matrix is symmetric so U = V')
138  Eigen::JacobiSVD<Eigen::Matrix3d> svd(cov, Eigen::ComputeFullU);
139  cov.setZero();
140  Eigen::Matrix3d U = svd.matrixU();
141  // Reconstitute the covariance matrix with modified singular values using the column
142  // // vectors in V.
143  for (int k = 0; k < 3; k++) {
144  Eigen::Vector3d col = U.col(k);
145  double v = 1.; // biggest 2 singular values replaced by 1
146  if (k == 2) // smallest singular value replaced by gicp_epsilon
147  v = gicp_epsilon_;
148  cov += v * col * col.transpose();
149  }
150  cloud_covariances[i] = cov;
151  }
152 }
153 
154 template <typename PointSource, typename PointTarget, typename Scalar>
155 void
157  double phi,
158  double theta,
159  double psi,
160  Eigen::Matrix3d& dR_dPhi,
161  Eigen::Matrix3d& dR_dTheta,
162  Eigen::Matrix3d& dR_dPsi) const
163 {
164  const double cphi = std::cos(phi), sphi = std::sin(phi);
165  const double ctheta = std::cos(theta), stheta = std::sin(theta);
166  const double cpsi = std::cos(psi), spsi = std::sin(psi);
167  dR_dPhi(0, 0) = 0.;
168  dR_dPhi(1, 0) = 0.;
169  dR_dPhi(2, 0) = 0.;
170 
171  dR_dPhi(0, 1) = sphi * spsi + cphi * cpsi * stheta;
172  dR_dPhi(1, 1) = -cpsi * sphi + cphi * spsi * stheta;
173  dR_dPhi(2, 1) = cphi * ctheta;
174 
175  dR_dPhi(0, 2) = cphi * spsi - cpsi * sphi * stheta;
176  dR_dPhi(1, 2) = -cphi * cpsi - sphi * spsi * stheta;
177  dR_dPhi(2, 2) = -ctheta * sphi;
178 
179  dR_dTheta(0, 0) = -cpsi * stheta;
180  dR_dTheta(1, 0) = -spsi * stheta;
181  dR_dTheta(2, 0) = -ctheta;
182 
183  dR_dTheta(0, 1) = cpsi * ctheta * sphi;
184  dR_dTheta(1, 1) = ctheta * sphi * spsi;
185  dR_dTheta(2, 1) = -sphi * stheta;
186 
187  dR_dTheta(0, 2) = cphi * cpsi * ctheta;
188  dR_dTheta(1, 2) = cphi * ctheta * spsi;
189  dR_dTheta(2, 2) = -cphi * stheta;
190 
191  dR_dPsi(0, 0) = -ctheta * spsi;
192  dR_dPsi(1, 0) = cpsi * ctheta;
193  dR_dPsi(2, 0) = 0.;
194 
195  dR_dPsi(0, 1) = -cphi * cpsi - sphi * spsi * stheta;
196  dR_dPsi(1, 1) = -cphi * spsi + cpsi * sphi * stheta;
197  dR_dPsi(2, 1) = 0.;
198 
199  dR_dPsi(0, 2) = cpsi * sphi - cphi * spsi * stheta;
200  dR_dPsi(1, 2) = sphi * spsi + cphi * cpsi * stheta;
201  dR_dPsi(2, 2) = 0.;
202 }
203 
204 template <typename PointSource, typename PointTarget, typename Scalar>
205 void
207  const Vector6d& x, const Eigen::Matrix3d& dCost_dR_T, Vector6d& g) const
208 {
209  Eigen::Matrix3d dR_dPhi;
210  Eigen::Matrix3d dR_dTheta;
211  Eigen::Matrix3d dR_dPsi;
212  getRDerivatives(x[3], x[4], x[5], dR_dPhi, dR_dTheta, dR_dPsi);
213 
214  g[3] = (dR_dPhi * dCost_dR_T).trace();
215  g[4] = (dR_dTheta * dCost_dR_T).trace();
216  g[5] = (dR_dPsi * dCost_dR_T).trace();
217 }
218 
219 template <typename PointSource, typename PointTarget, typename Scalar>
220 void
222  double phi,
223  double theta,
224  double psi,
225  Eigen::Matrix3d& ddR_dPhi_dPhi,
226  Eigen::Matrix3d& ddR_dPhi_dTheta,
227  Eigen::Matrix3d& ddR_dPhi_dPsi,
228  Eigen::Matrix3d& ddR_dTheta_dTheta,
229  Eigen::Matrix3d& ddR_dTheta_dPsi,
230  Eigen::Matrix3d& ddR_dPsi_dPsi) const
231 {
232  const double sphi = std::sin(phi), stheta = std::sin(theta), spsi = std::sin(psi);
233  const double cphi = std::cos(phi), ctheta = std::cos(theta), cpsi = std::cos(psi);
234  ddR_dPhi_dPhi(0, 0) = 0.0;
235  ddR_dPhi_dPhi(1, 0) = 0.0;
236  ddR_dPhi_dPhi(2, 0) = 0.0;
237  ddR_dPhi_dPhi(0, 1) = -cpsi * stheta * sphi + spsi * cphi;
238  ddR_dPhi_dPhi(1, 1) = -cpsi * cphi - spsi * stheta * sphi;
239  ddR_dPhi_dPhi(2, 1) = -ctheta * sphi;
240  ddR_dPhi_dPhi(0, 2) = -spsi * sphi - cpsi * stheta * cphi;
241  ddR_dPhi_dPhi(1, 2) = -spsi * stheta * cphi + cpsi * sphi;
242  ddR_dPhi_dPhi(2, 2) = -ctheta * cphi;
243 
244  ddR_dPhi_dTheta(0, 0) = 0.0;
245  ddR_dPhi_dTheta(1, 0) = 0.0;
246  ddR_dPhi_dTheta(2, 0) = 0.0;
247  ddR_dPhi_dTheta(0, 1) = cpsi * ctheta * cphi;
248  ddR_dPhi_dTheta(1, 1) = spsi * ctheta * cphi;
249  ddR_dPhi_dTheta(2, 1) = -stheta * cphi;
250  ddR_dPhi_dTheta(0, 2) = -cpsi * ctheta * sphi;
251  ddR_dPhi_dTheta(1, 2) = -spsi * ctheta * sphi;
252  ddR_dPhi_dTheta(2, 2) = stheta * sphi;
253 
254  ddR_dPhi_dPsi(0, 0) = 0.0;
255  ddR_dPhi_dPsi(1, 0) = 0.0;
256  ddR_dPhi_dPsi(2, 0) = 0.0;
257  ddR_dPhi_dPsi(0, 1) = -spsi * stheta * cphi + cpsi * sphi;
258  ddR_dPhi_dPsi(1, 1) = spsi * sphi + cpsi * stheta * cphi;
259  ddR_dPhi_dPsi(2, 1) = 0.0;
260  ddR_dPhi_dPsi(0, 2) = cpsi * cphi + spsi * stheta * sphi;
261  ddR_dPhi_dPsi(1, 2) = -cpsi * stheta * sphi + spsi * cphi;
262  ddR_dPhi_dPsi(2, 2) = 0.0;
263 
264  ddR_dTheta_dTheta(0, 0) = -cpsi * ctheta;
265  ddR_dTheta_dTheta(1, 0) = -spsi * ctheta;
266  ddR_dTheta_dTheta(2, 0) = stheta;
267  ddR_dTheta_dTheta(0, 1) = -cpsi * stheta * sphi;
268  ddR_dTheta_dTheta(1, 1) = -spsi * stheta * sphi;
269  ddR_dTheta_dTheta(2, 1) = -ctheta * sphi;
270  ddR_dTheta_dTheta(0, 2) = -cpsi * stheta * cphi;
271  ddR_dTheta_dTheta(1, 2) = -spsi * stheta * cphi;
272  ddR_dTheta_dTheta(2, 2) = -ctheta * cphi;
273 
274  ddR_dTheta_dPsi(0, 0) = spsi * stheta;
275  ddR_dTheta_dPsi(1, 0) = -cpsi * stheta;
276  ddR_dTheta_dPsi(2, 0) = 0.0;
277  ddR_dTheta_dPsi(0, 1) = -spsi * ctheta * sphi;
278  ddR_dTheta_dPsi(1, 1) = cpsi * ctheta * sphi;
279  ddR_dTheta_dPsi(2, 1) = 0.0;
280  ddR_dTheta_dPsi(0, 2) = -spsi * ctheta * cphi;
281  ddR_dTheta_dPsi(1, 2) = cpsi * ctheta * cphi;
282  ddR_dTheta_dPsi(2, 2) = 0.0;
283 
284  ddR_dPsi_dPsi(0, 0) = -cpsi * ctheta;
285  ddR_dPsi_dPsi(1, 0) = -spsi * ctheta;
286  ddR_dPsi_dPsi(2, 0) = 0.0;
287  ddR_dPsi_dPsi(0, 1) = -cpsi * stheta * sphi + spsi * cphi;
288  ddR_dPsi_dPsi(1, 1) = -cpsi * cphi - spsi * stheta * sphi;
289  ddR_dPsi_dPsi(2, 1) = 0.0;
290  ddR_dPsi_dPsi(0, 2) = -spsi * sphi - cpsi * stheta * cphi;
291  ddR_dPsi_dPsi(1, 2) = -spsi * stheta * cphi + cpsi * sphi;
292  ddR_dPsi_dPsi(2, 2) = 0.0;
293 }
294 
295 template <typename PointSource, typename PointTarget, typename Scalar>
296 void
299  const pcl::Indices& indices_src,
300  const PointCloudTarget& cloud_tgt,
301  const pcl::Indices& indices_tgt,
302  Matrix4& transformation_matrix)
303 {
304  // need at least min_number_correspondences_ samples
305  if (indices_src.size() < min_number_correspondences_) {
306  PCL_THROW_EXCEPTION(
308  "[pcl::GeneralizedIterativeClosestPoint::estimateRigidTransformationBFGS] Need "
309  "at least "
310  << min_number_correspondences_
311  << " points to estimate a transform! "
312  "Source and target have "
313  << indices_src.size() << " points!");
314  return;
315  }
316  // Set the initial solution
317  Vector6d x = Vector6d::Zero();
318  // translation part
319  x[0] = transformation_matrix(0, 3);
320  x[1] = transformation_matrix(1, 3);
321  x[2] = transformation_matrix(2, 3);
322  // rotation part (Z Y X euler angles convention)
323  // see: https://en.wikipedia.org/wiki/Rotation_matrix#General_rotations
324  x[3] = std::atan2(transformation_matrix(2, 1), transformation_matrix(2, 2));
325  x[4] = asin(-transformation_matrix(2, 0));
326  x[5] = std::atan2(transformation_matrix(1, 0), transformation_matrix(0, 0));
327 
328  // Set temporary pointers
329  tmp_src_ = &cloud_src;
330  tmp_tgt_ = &cloud_tgt;
331  tmp_idx_src_ = &indices_src;
332  tmp_idx_tgt_ = &indices_tgt;
333 
334  // Optimize using BFGS
335  OptimizationFunctorWithIndices functor(this);
337  bfgs.parameters.sigma = 0.01;
338  bfgs.parameters.rho = 0.01;
339  bfgs.parameters.tau1 = 9;
340  bfgs.parameters.tau2 = 0.05;
341  bfgs.parameters.tau3 = 0.5;
342  bfgs.parameters.order = 3;
343 
344  int inner_iterations_ = 0;
345  int result = bfgs.minimizeInit(x);
346  result = BFGSSpace::Running;
347  do {
348  inner_iterations_++;
349  result = bfgs.minimizeOneStep(x);
350  if (result) {
351  break;
352  }
353  result = bfgs.testGradient();
354  } while (result == BFGSSpace::Running && inner_iterations_ < max_inner_iterations_);
355  if (result == BFGSSpace::NoProgress || result == BFGSSpace::Success ||
356  inner_iterations_ == max_inner_iterations_) {
357  PCL_DEBUG("[pcl::registration::TransformationEstimationBFGS::"
358  "estimateRigidTransformation]");
359  PCL_DEBUG("BFGS solver finished with exit code %i \n", result);
360  transformation_matrix.setIdentity();
361  applyState(transformation_matrix, x);
362  }
363  else
364  PCL_THROW_EXCEPTION(
366  "[pcl::" << getClassName()
367  << "::TransformationEstimationBFGS::estimateRigidTransformation] BFGS "
368  "solver didn't converge!");
369 }
370 
371 template <typename PointSource, typename PointTarget, typename Scalar>
372 void
375  const pcl::Indices& indices_src,
376  const PointCloudTarget& cloud_tgt,
377  const pcl::Indices& indices_tgt,
378  Matrix4& transformation_matrix)
379 {
380  // need at least min_number_correspondences_ samples
381  if (indices_src.size() < min_number_correspondences_) {
382  PCL_THROW_EXCEPTION(NotEnoughPointsException,
383  "[pcl::GeneralizedIterativeClosestPoint::"
384  "estimateRigidTransformationNewton] Need "
385  "at least "
386  << min_number_correspondences_
387  << " points to estimate a transform! "
388  "Source and target have "
389  << indices_src.size() << " points!");
390  return;
391  }
392  // Set the initial solution
393  Vector6d x = Vector6d::Zero();
394  // translation part
395  x[0] = transformation_matrix(0, 3);
396  x[1] = transformation_matrix(1, 3);
397  x[2] = transformation_matrix(2, 3);
398  // rotation part (Z Y X euler angles convention)
399  // see: https://en.wikipedia.org/wiki/Rotation_matrix#General_rotations
400  x[3] = std::atan2(transformation_matrix(2, 1), transformation_matrix(2, 2));
401  x[4] = std::asin(
402  std::min<double>(1.0, std::max<double>(-1.0, -transformation_matrix(2, 0))));
403  x[5] = std::atan2(transformation_matrix(1, 0), transformation_matrix(0, 0));
404 
405  // Set temporary pointers
406  tmp_src_ = &cloud_src;
407  tmp_tgt_ = &cloud_tgt;
408  tmp_idx_src_ = &indices_src;
409  tmp_idx_tgt_ = &indices_tgt;
410 
411  // Optimize using Newton
412  OptimizationFunctorWithIndices functor(this);
413  Eigen::Matrix<double, 6, 6> hessian;
414  Eigen::Matrix<double, 6, 1> gradient;
415  double current_x_value = functor(x);
416  functor.dfddf(x, gradient, hessian);
417  Eigen::Matrix<double, 6, 1> delta;
418  int inner_iterations_ = 0;
419  do {
420  ++inner_iterations_;
421  // compute descent direction from hessian and gradient. Take special measures if
422  // hessian is not positive-definite (positive Eigenvalues)
423  Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double, 6, 6>> eigensolver(hessian);
424  Eigen::Matrix<double, 6, 6> inverted_eigenvalues =
425  Eigen::Matrix<double, 6, 6>::Zero();
426  for (int i = 0; i < 6; ++i) {
427  const double ev = eigensolver.eigenvalues()[i];
428  if (ev < 0)
429  inverted_eigenvalues(i, i) = 1.0 / eigensolver.eigenvalues()[5];
430  else
431  inverted_eigenvalues(i, i) = 1.0 / ev;
432  }
433  delta = eigensolver.eigenvectors() * inverted_eigenvalues *
434  eigensolver.eigenvectors().transpose() * gradient;
435 
436  // simple line search to guarantee a decrease in the function value
437  double alpha = 1.0;
438  double candidate_x_value;
439  bool improvement_found = false;
440  for (int i = 0; i < 10; ++i, alpha /= 2) {
441  Vector6d candidate_x = x - alpha * delta;
442  candidate_x_value = functor(candidate_x);
443  if (candidate_x_value < current_x_value) {
444  PCL_DEBUG("[estimateRigidTransformationNewton] Using stepsize=%g, function "
445  "value previously: %g, now: %g, "
446  "improvement: %g\n",
447  alpha,
448  current_x_value,
449  candidate_x_value,
450  current_x_value - candidate_x_value);
451  x = candidate_x;
452  current_x_value = candidate_x_value;
453  improvement_found = true;
454  break;
455  }
456  }
457  if (!improvement_found) {
458  PCL_DEBUG("[estimateRigidTransformationNewton] finishing because no progress\n");
459  break;
460  }
461  functor.dfddf(x, gradient, hessian);
462  if (gradient.head<3>().norm() < translation_gradient_tolerance_ &&
463  gradient.tail<3>().norm() < rotation_gradient_tolerance_) {
464  PCL_DEBUG("[estimateRigidTransformationNewton] finishing because gradient below "
465  "threshold: translation: %g<%g, rotation: %g<%g\n",
466  gradient.head<3>().norm(),
467  translation_gradient_tolerance_,
468  gradient.tail<3>().norm(),
469  rotation_gradient_tolerance_);
470  break;
471  }
472  } while (inner_iterations_ < max_inner_iterations_);
473  PCL_DEBUG("[estimateRigidTransformationNewton] solver finished after %i iterations "
474  "(of max %i)\n",
475  inner_iterations_,
476  max_inner_iterations_);
477  transformation_matrix.setIdentity();
478  applyState(transformation_matrix, x);
479 }
480 
481 template <typename PointSource, typename PointTarget, typename Scalar>
482 inline double
485 {
486  Matrix4 transformation_matrix = gicp_->base_transformation_;
487  gicp_->applyState(transformation_matrix, x);
488  double f = 0;
489  int m = static_cast<int>(gicp_->tmp_idx_src_->size());
490  for (int i = 0; i < m; ++i) {
491  // The last coordinate, p_src[3] is guaranteed to be set to 1.0 in registration.hpp
492  Vector4fMapConst p_src =
493  (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
494  // The last coordinate, p_tgt[3] is guaranteed to be set to 1.0 in registration.hpp
495  Vector4fMapConst p_tgt =
496  (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
497  Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
498  // Estimate the distance (cost function)
499  // The last coordinate is still guaranteed to be set to 1.0
500  // The d here is the negative of the d in the paper
501  Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
502  p_trans_src[1] - p_tgt[1],
503  p_trans_src[2] - p_tgt[2]);
504  Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
505  // increment= d'*Md/num_matches = d'*M*d/num_matches (we postpone
506  // 1/num_matches after the loop closes)
507  f += static_cast<double>(d.transpose() * Md);
508  }
509  return f / m;
510 }
511 
512 template <typename PointSource, typename PointTarget, typename Scalar>
513 inline void
516 {
517  Matrix4 transformation_matrix = gicp_->base_transformation_;
518  gicp_->applyState(transformation_matrix, x);
519  // Zero out g
520  g.setZero();
521  // Eigen::Vector3d g_t = g.head<3> ();
522  // the transpose of the derivative of the cost function w.r.t rotation matrix
523  Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
524  int m = static_cast<int>(gicp_->tmp_idx_src_->size());
525  for (int i = 0; i < m; ++i) {
526  // The last coordinate, p_src[3] is guaranteed to be set to 1.0 in registration.hpp
527  Vector4fMapConst p_src =
528  (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
529  // The last coordinate, p_tgt[3] is guaranteed to be set to 1.0 in registration.hpp
530  Vector4fMapConst p_tgt =
531  (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
532 
533  Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
534  // The last coordinate is still guaranteed to be set to 1.0
535  // The d here is the negative of the d in the paper
536  Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
537  p_trans_src[1] - p_tgt[1],
538  p_trans_src[2] - p_tgt[2]);
539  // Md = M*d
540  Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
541  // Increment translation gradient
542  // g.head<3> ()+= 2*M*d/num_matches (we postpone 2/num_matches after the loop
543  // closes)
544  g.head<3>() += Md;
545  // Increment rotation gradient
546  p_trans_src = gicp_->base_transformation_.template cast<float>() * p_src;
547  Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
548  dCost_dR_T += p_base_src * Md.transpose();
549  }
550  g.head<3>() *= 2.0 / m;
551  dCost_dR_T *= 2.0 / m;
552  gicp_->computeRDerivative(x, dCost_dR_T, g);
553 }
554 
555 template <typename PointSource, typename PointTarget, typename Scalar>
556 inline void
559 {
560  Matrix4 transformation_matrix = gicp_->base_transformation_;
561  gicp_->applyState(transformation_matrix, x);
562  f = 0;
563  g.setZero();
564  // the transpose of the derivative of the cost function w.r.t rotation matrix
565  Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
566  const int m = static_cast<int>(gicp_->tmp_idx_src_->size());
567  for (int i = 0; i < m; ++i) {
568  // The last coordinate, p_src[3] is guaranteed to be set to 1.0 in registration.hpp
569  Vector4fMapConst p_src =
570  (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
571  // The last coordinate, p_tgt[3] is guaranteed to be set to 1.0 in registration.hpp
572  Vector4fMapConst p_tgt =
573  (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
574  Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
575  // The last coordinate is still guaranteed to be set to 1.0
576  // The d here is the negative of the d in the paper
577  Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
578  p_trans_src[1] - p_tgt[1],
579  p_trans_src[2] - p_tgt[2]);
580  // Md = M*d
581  Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
582  // Increment total error
583  f += static_cast<double>(d.transpose() * Md);
584  // Increment translation gradient
585  // g.head<3> ()+= 2*M*d/num_matches (we postpone 2/num_matches after the loop
586  // closes)
587  g.head<3>() += Md;
588  p_trans_src = gicp_->base_transformation_.template cast<float>() * p_src;
589  Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
590  // Increment rotation gradient
591  dCost_dR_T += p_base_src * Md.transpose();
592  }
593  f /= static_cast<double>(m);
594  g.head<3>() *= (2.0 / m);
595  dCost_dR_T *= 2.0 / m;
596  gicp_->computeRDerivative(x, dCost_dR_T, g);
597 }
598 
599 template <typename PointSource, typename PointTarget, typename Scalar>
600 inline void
603  Vector6d& gradient,
604  Matrix6d& hessian)
605 {
606  Matrix4 transformation_matrix = gicp_->base_transformation_;
607  gicp_->applyState(transformation_matrix, x);
608  const Eigen::Matrix4f transformation_matrix_float =
609  transformation_matrix.template cast<float>();
610  const Eigen::Matrix4f base_transformation_float =
611  gicp_->base_transformation_.template cast<float>();
612  // Zero out gradient and hessian
613  gradient.setZero();
614  hessian.setZero();
615  // Helper matrices
616  Eigen::Matrix3d dR_dPhi;
617  Eigen::Matrix3d dR_dTheta;
618  Eigen::Matrix3d dR_dPsi;
619  gicp_->getRDerivatives(x[3], x[4], x[5], dR_dPhi, dR_dTheta, dR_dPsi);
620  Eigen::Matrix3d ddR_dPhi_dPhi;
621  Eigen::Matrix3d ddR_dPhi_dTheta;
622  Eigen::Matrix3d ddR_dPhi_dPsi;
623  Eigen::Matrix3d ddR_dTheta_dTheta;
624  Eigen::Matrix3d ddR_dTheta_dPsi;
625  Eigen::Matrix3d ddR_dPsi_dPsi;
626  gicp_->getR2ndDerivatives(x[3],
627  x[4],
628  x[5],
629  ddR_dPhi_dPhi,
630  ddR_dPhi_dTheta,
631  ddR_dPhi_dPsi,
632  ddR_dTheta_dTheta,
633  ddR_dTheta_dPsi,
634  ddR_dPsi_dPsi);
635  Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
636  Eigen::Matrix3d dCost_dR_T1 = Eigen::Matrix3d::Zero();
637  Eigen::Matrix3d dCost_dR_T2 = Eigen::Matrix3d::Zero();
638  Eigen::Matrix3d dCost_dR_T3 = Eigen::Matrix3d::Zero();
639  Eigen::Matrix3d dCost_dR_T1b = Eigen::Matrix3d::Zero();
640  Eigen::Matrix3d dCost_dR_T2b = Eigen::Matrix3d::Zero();
641  Eigen::Matrix3d dCost_dR_T3b = Eigen::Matrix3d::Zero();
642  Eigen::Matrix3d hessian_rot_phi = Eigen::Matrix3d::Zero();
643  Eigen::Matrix3d hessian_rot_theta = Eigen::Matrix3d::Zero();
644  Eigen::Matrix3d hessian_rot_psi = Eigen::Matrix3d::Zero();
645  Eigen::Matrix<double, 9, 6> hessian_rot_tmp = Eigen::Matrix<double, 9, 6>::Zero();
646 
647  int m = static_cast<int>(gicp_->tmp_idx_src_->size());
648  for (int i = 0; i < m; ++i) {
649  // The last coordinate, p_src[3] is guaranteed to be set to 1.0 in registration.hpp
650  const auto& src_idx = (*gicp_->tmp_idx_src_)[i];
651  Vector4fMapConst p_src = (*gicp_->tmp_src_)[src_idx].getVector4fMap();
652  // The last coordinate, p_tgt[3] is guaranteed to be set to 1.0 in registration.hpp
653  Vector4fMapConst p_tgt =
654  (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
655  Eigen::Vector4f p_trans_src(transformation_matrix_float * p_src);
656  // The last coordinate is still guaranteed to be set to 1.0
657  // The d here is the negative of the d in the paper
658  const Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
659  p_trans_src[1] - p_tgt[1],
660  p_trans_src[2] - p_tgt[2]);
661  const Eigen::Matrix3d& M = gicp_->mahalanobis(src_idx);
662  const Eigen::Vector3d Md(M * d); // Md = M*d
663  gradient.head<3>() += Md; // translation gradient
664  hessian.topLeftCorner<3, 3>() += M; // translation-translation hessian
665  p_trans_src = base_transformation_float * p_src;
666  const Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
667  dCost_dR_T.noalias() += p_base_src * Md.transpose();
668  dCost_dR_T1b += p_base_src[0] * M;
669  dCost_dR_T2b += p_base_src[1] * M;
670  dCost_dR_T3b += p_base_src[2] * M;
671  hessian_rot_tmp.noalias() +=
672  Eigen::Map<const Eigen::Matrix<double, 9, 1>>{M.data()} *
673  (Eigen::Matrix<double, 1, 6>() << p_base_src[0] * p_base_src[0],
674  p_base_src[0] * p_base_src[1],
675  p_base_src[0] * p_base_src[2],
676  p_base_src[1] * p_base_src[1],
677  p_base_src[1] * p_base_src[2],
678  p_base_src[2] * p_base_src[2])
679  .finished();
680  }
681  gradient.head<3>() *= 2.0 / m; // translation gradient
682  dCost_dR_T *= 2.0 / m;
683  gicp_->computeRDerivative(x, dCost_dR_T, gradient); // rotation gradient
684  hessian.topLeftCorner<3, 3>() *= 2.0 / m; // translation-translation hessian
685  // translation-rotation hessian
686  dCost_dR_T1.row(0) = dCost_dR_T1b.col(0);
687  dCost_dR_T1.row(1) = dCost_dR_T2b.col(0);
688  dCost_dR_T1.row(2) = dCost_dR_T3b.col(0);
689  dCost_dR_T2.row(0) = dCost_dR_T1b.col(1);
690  dCost_dR_T2.row(1) = dCost_dR_T2b.col(1);
691  dCost_dR_T2.row(2) = dCost_dR_T3b.col(1);
692  dCost_dR_T3.row(0) = dCost_dR_T1b.col(2);
693  dCost_dR_T3.row(1) = dCost_dR_T2b.col(2);
694  dCost_dR_T3.row(2) = dCost_dR_T3b.col(2);
695  dCost_dR_T1 *= 2.0 / m;
696  dCost_dR_T2 *= 2.0 / m;
697  dCost_dR_T3 *= 2.0 / m;
698  hessian(3, 0) = (dR_dPhi * dCost_dR_T1).trace();
699  hessian(4, 0) = (dR_dTheta * dCost_dR_T1).trace();
700  hessian(5, 0) = (dR_dPsi * dCost_dR_T1).trace();
701  hessian(3, 1) = (dR_dPhi * dCost_dR_T2).trace();
702  hessian(4, 1) = (dR_dTheta * dCost_dR_T2).trace();
703  hessian(5, 1) = (dR_dPsi * dCost_dR_T2).trace();
704  hessian(3, 2) = (dR_dPhi * dCost_dR_T3).trace();
705  hessian(4, 2) = (dR_dTheta * dCost_dR_T3).trace();
706  hessian(5, 2) = (dR_dPsi * dCost_dR_T3).trace();
707  hessian.block<3, 3>(0, 3) = hessian.block<3, 3>(3, 0).transpose();
708  // rotation-rotation hessian
709  int lookup[3][3] = {{0, 1, 2}, {1, 3, 4}, {2, 4, 5}};
710  for (int l = 0; l < 3; ++l) {
711  for (int i = 0; i < 3; ++i) {
712  double phi_tmp = 0.0, theta_tmp = 0.0, psi_tmp = 0.0;
713  for (int j = 0; j < 3; ++j) {
714  for (int k = 0; k < 3; ++k) {
715  phi_tmp += hessian_rot_tmp(3 * j + i, lookup[l][k]) * dR_dPhi(j, k);
716  theta_tmp += hessian_rot_tmp(3 * j + i, lookup[l][k]) * dR_dTheta(j, k);
717  psi_tmp += hessian_rot_tmp(3 * j + i, lookup[l][k]) * dR_dPsi(j, k);
718  }
719  }
720  hessian_rot_phi(i, l) = phi_tmp;
721  hessian_rot_theta(i, l) = theta_tmp;
722  hessian_rot_psi(i, l) = psi_tmp;
723  }
724  }
725  hessian_rot_phi *= 2.0 / m;
726  hessian_rot_theta *= 2.0 / m;
727  hessian_rot_psi *= 2.0 / m;
728  hessian(3, 3) = (dR_dPhi.transpose() * hessian_rot_phi).trace() +
729  (ddR_dPhi_dPhi * dCost_dR_T).trace();
730  hessian(3, 4) = (dR_dPhi.transpose() * hessian_rot_theta).trace() +
731  (ddR_dPhi_dTheta * dCost_dR_T).trace();
732  hessian(3, 5) = (dR_dPhi.transpose() * hessian_rot_psi).trace() +
733  (ddR_dPhi_dPsi * dCost_dR_T).trace();
734  hessian(4, 4) = (dR_dTheta.transpose() * hessian_rot_theta).trace() +
735  (ddR_dTheta_dTheta * dCost_dR_T).trace();
736  hessian(4, 5) = (dR_dTheta.transpose() * hessian_rot_psi).trace() +
737  (ddR_dTheta_dPsi * dCost_dR_T).trace();
738  hessian(5, 5) = (dR_dPsi.transpose() * hessian_rot_psi).trace() +
739  (ddR_dPsi_dPsi * dCost_dR_T).trace();
740  hessian(4, 3) = hessian(3, 4);
741  hessian(5, 3) = hessian(3, 5);
742  hessian(5, 4) = hessian(4, 5);
743 }
744 
745 template <typename PointSource, typename PointTarget, typename Scalar>
746 inline BFGSSpace::Status
749 {
750  auto translation_epsilon = gicp_->translation_gradient_tolerance_;
751  auto rotation_epsilon = gicp_->rotation_gradient_tolerance_;
752 
753  if ((translation_epsilon < 0.) || (rotation_epsilon < 0.))
755 
756  // express translation gradient as norm of translation parameters
757  auto translation_grad = g.head<3>().norm();
758 
759  // express rotation gradient as a norm of rotation parameters
760  auto rotation_grad = g.tail<3>().norm();
761 
762  if ((translation_grad < translation_epsilon) && (rotation_grad < rotation_epsilon))
763  return BFGSSpace::Success;
764 
765  return BFGSSpace::Running;
766 }
767 
768 template <typename PointSource, typename PointTarget, typename Scalar>
769 inline void
771  computeTransformation(PointCloudSource& output, const Matrix4& guess)
772 {
774  // Difference between consecutive transforms
775  double delta = 0;
776  // Get the size of the source point cloud
777  const std::size_t N = indices_->size();
778  // Set the mahalanobis matrices to identity
779  mahalanobis_.resize(N, Eigen::Matrix3d::Identity());
780  // Compute target cloud covariance matrices
781  if ((!target_covariances_) || (target_covariances_->empty())) {
782  target_covariances_.reset(new MatricesVector);
783  computeCovariances<PointTarget>(target_, tree_, *target_covariances_);
784  }
785  // Compute input cloud covariance matrices
786  if ((!input_covariances_) || (input_covariances_->empty())) {
787  input_covariances_.reset(new MatricesVector);
788  computeCovariances<PointSource>(input_, tree_reciprocal_, *input_covariances_);
789  }
790 
791  base_transformation_ = Matrix4::Identity();
792  nr_iterations_ = 0;
793  converged_ = false;
794  double dist_threshold = corr_dist_threshold_ * corr_dist_threshold_;
795  pcl::Indices nn_indices(1);
796  std::vector<float> nn_dists(1);
797 
798  pcl::transformPointCloud(output, output, guess);
799 
800  while (!converged_) {
801  std::size_t cnt = 0;
802  pcl::Indices source_indices(indices_->size());
803  pcl::Indices target_indices(indices_->size());
804 
805  // guess corresponds to base_t and transformation_ to t
806  Eigen::Matrix4d transform_R = Eigen::Matrix4d::Zero();
807  for (std::size_t i = 0; i < 4; i++)
808  for (std::size_t j = 0; j < 4; j++)
809  for (std::size_t k = 0; k < 4; k++)
810  transform_R(i, j) += static_cast<double>(transformation_(i, k)) *
811  static_cast<double>(guess(k, j));
812 
813  Eigen::Matrix3d R = transform_R.topLeftCorner<3, 3>();
814 
815  for (std::size_t i = 0; i < N; i++) {
816  PointSource query = output[i];
817  query.getVector4fMap() =
818  transformation_.template cast<float>() * query.getVector4fMap();
819 
820  if (!searchForNeighbors(query, nn_indices, nn_dists)) {
821  PCL_ERROR("[pcl::%s::computeTransformation] Unable to find a nearest neighbor "
822  "in the target dataset for point %d in the source!\n",
823  getClassName().c_str(),
824  (*indices_)[i]);
825  return;
826  }
827 
828  // Check if the distance to the nearest neighbor is smaller than the user imposed
829  // threshold
830  if (nn_dists[0] < dist_threshold) {
831  Eigen::Matrix3d& C1 = (*input_covariances_)[i];
832  Eigen::Matrix3d& C2 = (*target_covariances_)[nn_indices[0]];
833  Eigen::Matrix3d& M = mahalanobis_[i];
834  // M = R*C1
835  M = R * C1;
836  // temp = M*R' + C2 = R*C1*R' + C2
837  Eigen::Matrix3d temp = M * R.transpose();
838  temp += C2;
839  // M = temp^-1
840  M = temp.inverse();
841  source_indices[cnt] = static_cast<int>(i);
842  target_indices[cnt] = nn_indices[0];
843  cnt++;
844  }
845  }
846  // Resize to the actual number of valid correspondences
847  source_indices.resize(cnt);
848  target_indices.resize(cnt);
849  /* optimize transformation using the current assignment and Mahalanobis metrics*/
850  previous_transformation_ = transformation_;
851  // optimization right here
852  try {
853  rigid_transformation_estimation_(
854  output, source_indices, *target_, target_indices, transformation_);
855  /* compute the delta from this iteration */
856  delta = 0.;
857  for (int k = 0; k < 4; k++) {
858  for (int l = 0; l < 4; l++) {
859  double ratio = 1;
860  if (k < 3 && l < 3) // rotation part of the transform
861  ratio = 1. / rotation_epsilon_;
862  else
863  ratio = 1. / transformation_epsilon_;
864  double c_delta =
865  ratio * std::abs(previous_transformation_(k, l) - transformation_(k, l));
866  if (c_delta > delta)
867  delta = c_delta;
868  }
869  }
870  } catch (PCLException& e) {
871  PCL_DEBUG("[pcl::%s::computeTransformation] Optimization issue %s\n",
872  getClassName().c_str(),
873  e.what());
874  break;
875  }
876  nr_iterations_++;
877 
878  if (update_visualizer_ != nullptr) {
879  PointCloudSourcePtr input_transformed(new PointCloudSource);
880  pcl::transformPointCloud(output, *input_transformed, transformation_);
881  update_visualizer_(*input_transformed, source_indices, *target_, target_indices);
882  }
883 
884  // Check for convergence
885  if (nr_iterations_ >= max_iterations_ || delta < 1) {
886  converged_ = true;
887  PCL_DEBUG("[pcl::%s::computeTransformation] Convergence reached. Number of "
888  "iterations: %d out of %d. Transformation difference: %f\n",
889  getClassName().c_str(),
890  nr_iterations_,
891  max_iterations_,
892  (transformation_ - previous_transformation_).array().abs().sum());
893  previous_transformation_ = transformation_;
894  }
895  else
896  PCL_DEBUG("[pcl::%s::computeTransformation] Convergence failed\n",
897  getClassName().c_str());
898  }
899  final_transformation_ = previous_transformation_ * guess;
900 
901  PCL_DEBUG("Transformation "
902  "is:\n\t%5f\t%5f\t%5f\t%5f\n\t%5f\t%5f\t%5f\t%5f\n\t%5f\t%5f\t%5f\t%5f\n\t%"
903  "5f\t%5f\t%5f\t%5f\n",
904  final_transformation_(0, 0),
905  final_transformation_(0, 1),
906  final_transformation_(0, 2),
907  final_transformation_(0, 3),
908  final_transformation_(1, 0),
909  final_transformation_(1, 1),
910  final_transformation_(1, 2),
911  final_transformation_(1, 3),
912  final_transformation_(2, 0),
913  final_transformation_(2, 1),
914  final_transformation_(2, 2),
915  final_transformation_(2, 3),
916  final_transformation_(3, 0),
917  final_transformation_(3, 1),
918  final_transformation_(3, 2),
919  final_transformation_(3, 3));
920 
921  // Transform the point cloud
922  pcl::transformPointCloud(*input_, output, final_transformation_);
923 }
924 
925 template <typename PointSource, typename PointTarget, typename Scalar>
926 void
928  Matrix4& t, const Vector6d& x) const
929 {
930  // Z Y X euler angles convention
931  Matrix3 R = (AngleAxis(static_cast<Scalar>(x[5]), Vector3::UnitZ()) *
932  AngleAxis(static_cast<Scalar>(x[4]), Vector3::UnitY()) *
933  AngleAxis(static_cast<Scalar>(x[3]), Vector3::UnitX()))
934  .toRotationMatrix();
935  Matrix4 T = Matrix4::Identity();
936  T.template block<3, 3>(0, 0) = R;
937  T.template block<3, 1>(0, 3) = Vector3(
938  static_cast<Scalar>(x[0]), static_cast<Scalar>(x[1]), static_cast<Scalar>(x[2]));
939  t = T * t;
940 }
941 
942 } // namespace pcl
943 
944 #endif // PCL_REGISTRATION_IMPL_GICP_HPP_
BFGS stands for Broyden–Fletcher–Goldfarb–Shanno (BFGS) method for solving unconstrained nonlinear op...
Definition: bfgs.h:121
BFGSSpace::Status testGradient()
Definition: bfgs.h:476
BFGSSpace::Status minimizeInit(FVectorType &x)
Definition: bfgs.h:361
BFGSSpace::Status minimizeOneStep(FVectorType &x)
Definition: bfgs.h:393
Parameters parameters
Definition: bfgs.h:169
GeneralizedIterativeClosestPoint is an ICP variant that implements the generalized iterative closest ...
Definition: gicp.h:76
void estimateRigidTransformationBFGS(const PointCloudSource &cloud_src, const pcl::Indices &indices_src, const PointCloudTarget &cloud_tgt, const pcl::Indices &indices_tgt, Matrix4 &transformation_matrix)
Estimate a rigid rotation transformation between a source and a target point cloud using an iterative...
Definition: gicp.hpp:298
typename IterativeClosestPoint< PointSource, PointTarget, Scalar >::Matrix4 Matrix4
Definition: gicp.h:130
void applyState(Matrix4 &t, const Vector6d &x) const
compute transformation matrix from transformation matrix
Definition: gicp.hpp:927
Eigen::Matrix< double, 6, 6 > Matrix6d
Definition: gicp.h:131
typename Eigen::Matrix< Scalar, 3, 1 > Vector3
Definition: gicp.h:125
std::vector< Eigen::Matrix3d, Eigen::aligned_allocator< Eigen::Matrix3d > > MatricesVector
Definition: gicp.h:112
void computeCovariances(typename pcl::PointCloud< PointT >::ConstPtr cloud, const typename pcl::search::KdTree< PointT >::Ptr tree, MatricesVector &cloud_covariances)
compute points covariances matrices according to the K nearest neighbors.
Definition: gicp.hpp:73
void estimateRigidTransformationNewton(const PointCloudSource &cloud_src, const pcl::Indices &indices_src, const PointCloudTarget &cloud_tgt, const pcl::Indices &indices_tgt, Matrix4 &transformation_matrix)
Estimate a rigid rotation transformation between a source and a target point cloud using an iterative...
Definition: gicp.hpp:374
typename PointCloudSource::Ptr PointCloudSourcePtr
Definition: gicp.h:101
void setNumberOfThreads(unsigned int nr_threads=0)
Initialize the scheduler and set the number of threads to use.
Definition: gicp.hpp:50
void computeRDerivative(const Vector6d &x, const Eigen::Matrix3d &dCost_dR_T, Vector6d &g) const
Computes the derivative of the cost function w.r.t rotation angles.
Definition: gicp.hpp:206
typename Eigen::AngleAxis< Scalar > AngleAxis
Definition: gicp.h:132
void computeTransformation(PointCloudSource &output, const Matrix4 &guess) override
Rigid transformation computation method with initial guess.
Definition: gicp.hpp:771
typename Eigen::Matrix< Scalar, 3, 3 > Matrix3
Definition: gicp.h:128
Eigen::Matrix< double, 6, 1 > Vector6d
Definition: gicp.h:127
An exception that is thrown when the number of correspondents is not equal to the minimum required.
Definition: exceptions.h:63
A base class for all pcl exceptions which inherits from std::runtime_error.
Definition: exceptions.h:66
std::size_t size() const
Definition: point_cloud.h:443
shared_ptr< const PointCloud< PointT > > ConstPtr
Definition: point_cloud.h:414
bool initComputeReciprocal()
Internal computation when reciprocal lookup is needed.
An exception that is thrown when the non linear solver didn't converge.
Definition: exceptions.h:49
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
shared_ptr< KdTree< PointT, Tree > > Ptr
Definition: kdtree.h:75
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
Status
Definition: bfgs.h:70
@ NoProgress
Definition: bfgs.h:75
@ Running
Definition: bfgs.h:73
@ Success
Definition: bfgs.h:74
@ NegativeGradientEpsilon
Definition: bfgs.h:71
const Eigen::Map< const Eigen::Vector4f, Eigen::Aligned > Vector4fMapConst
IndicesAllocator<> Indices
Type used for indices in PCL.
Definition: types.h:133
Scalar sigma
Definition: bfgs.h:147
Scalar tau2
Definition: bfgs.h:149
Scalar rho
Definition: bfgs.h:146
Scalar tau1
Definition: bfgs.h:148
Scalar tau3
Definition: bfgs.h:150
Index order
Definition: bfgs.h:152
void dfddf(const Vector6d &x, Vector6d &df, Matrix6d &ddf)
Definition: gicp.hpp:602
void df(const Vector6d &x, Vector6d &df) override
Definition: gicp.hpp:515
BFGSSpace::Status checkGradient(const Vector6d &g) override
Definition: gicp.hpp:748
void fdf(const Vector6d &x, double &f, Vector6d &df) override
Definition: gicp.hpp:558
A point structure representing Euclidean xyz coordinates, and the RGB color.