41 #ifndef PCL_REGISTRATION_IMPL_GICP_HPP_
42 #define PCL_REGISTRATION_IMPL_GICP_HPP_
44 #include <pcl/registration/exceptions.h>
48 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
49 template <
typename Po
intT>
56 if (k_correspondences_ >
static_cast<int>(cloud->
size())) {
57 PCL_ERROR(
"[pcl::GeneralizedIterativeClosestPoint::computeCovariances] Number or "
58 "points in cloud (%lu) is less than k_correspondences_ (%lu)!\n",
66 nn_indices.reserve(k_correspondences_);
67 std::vector<float> nn_dist_sq;
68 nn_dist_sq.reserve(k_correspondences_);
71 if (cloud_covariances.size() < cloud->
size())
72 cloud_covariances.resize(cloud->
size());
74 auto matrices_iterator = cloud_covariances.begin();
75 for (
auto points_iterator = cloud->
begin(); points_iterator != cloud->
end();
76 ++points_iterator, ++matrices_iterator) {
77 const PointT& query_point = *points_iterator;
78 Eigen::Matrix3d& cov = *matrices_iterator;
84 kdtree->
nearestKSearch(query_point, k_correspondences_, nn_indices, nn_dist_sq);
87 for (
int j = 0; j < k_correspondences_; j++) {
88 const PointT& pt = (*cloud)[nn_indices[j]];
94 cov(0, 0) += pt.x * pt.x;
96 cov(1, 0) += pt.y * pt.x;
97 cov(1, 1) += pt.y * pt.y;
99 cov(2, 0) += pt.z * pt.x;
100 cov(2, 1) += pt.z * pt.y;
101 cov(2, 2) += pt.z * pt.z;
104 mean /=
static_cast<double>(k_correspondences_);
106 for (
int k = 0; k < 3; k++)
107 for (
int l = 0; l <= k; l++) {
108 cov(k, l) /=
static_cast<double>(k_correspondences_);
109 cov(k, l) -= mean[k] * mean[l];
110 cov(l, k) = cov(k, l);
114 Eigen::JacobiSVD<Eigen::Matrix3d> svd(cov, Eigen::ComputeFullU);
116 Eigen::Matrix3d U = svd.matrixU();
119 for (
int k = 0; k < 3; k++) {
120 Eigen::Vector3d col = U.col(k);
124 cov += v * col * col.transpose();
129 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
134 Eigen::Matrix3d dR_dPhi;
135 Eigen::Matrix3d dR_dTheta;
136 Eigen::Matrix3d dR_dPsi;
138 double phi = x[3], theta = x[4], psi = x[5];
140 double cphi = std::cos(phi), sphi = sin(phi);
141 double ctheta = std::cos(theta), stheta = sin(theta);
142 double cpsi = std::cos(psi), spsi = sin(psi);
148 dR_dPhi(0, 1) = sphi * spsi + cphi * cpsi * stheta;
149 dR_dPhi(1, 1) = -cpsi * sphi + cphi * spsi * stheta;
150 dR_dPhi(2, 1) = cphi * ctheta;
152 dR_dPhi(0, 2) = cphi * spsi - cpsi * sphi * stheta;
153 dR_dPhi(1, 2) = -cphi * cpsi - sphi * spsi * stheta;
154 dR_dPhi(2, 2) = -ctheta * sphi;
156 dR_dTheta(0, 0) = -cpsi * stheta;
157 dR_dTheta(1, 0) = -spsi * stheta;
158 dR_dTheta(2, 0) = -ctheta;
160 dR_dTheta(0, 1) = cpsi * ctheta * sphi;
161 dR_dTheta(1, 1) = ctheta * sphi * spsi;
162 dR_dTheta(2, 1) = -sphi * stheta;
164 dR_dTheta(0, 2) = cphi * cpsi * ctheta;
165 dR_dTheta(1, 2) = cphi * ctheta * spsi;
166 dR_dTheta(2, 2) = -cphi * stheta;
168 dR_dPsi(0, 0) = -ctheta * spsi;
169 dR_dPsi(1, 0) = cpsi * ctheta;
172 dR_dPsi(0, 1) = -cphi * cpsi - sphi * spsi * stheta;
173 dR_dPsi(1, 1) = -cphi * spsi + cpsi * sphi * stheta;
176 dR_dPsi(0, 2) = cpsi * sphi - cphi * spsi * stheta;
177 dR_dPsi(1, 2) = sphi * spsi + cphi * cpsi * stheta;
180 g[3] = matricesInnerProd(dR_dPhi, dCost_dR_T);
181 g[4] = matricesInnerProd(dR_dTheta, dCost_dR_T);
182 g[5] = matricesInnerProd(dR_dPsi, dCost_dR_T);
185 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
192 Matrix4& transformation_matrix)
195 if (indices_src.size() < min_number_correspondences_) {
198 "[pcl::GeneralizedIterativeClosestPoint::estimateRigidTransformationBFGS] Need "
200 << min_number_correspondences_
201 <<
" points to estimate a transform! "
202 "Source and target have "
203 << indices_src.size() <<
" points!");
209 x[0] = transformation_matrix(0, 3);
210 x[1] = transformation_matrix(1, 3);
211 x[2] = transformation_matrix(2, 3);
214 x[3] = std::atan2(transformation_matrix(2, 1), transformation_matrix(2, 2));
215 x[4] = asin(-transformation_matrix(2, 0));
216 x[5] = std::atan2(transformation_matrix(1, 0), transformation_matrix(0, 0));
219 tmp_src_ = &cloud_src;
220 tmp_tgt_ = &cloud_tgt;
221 tmp_idx_src_ = &indices_src;
222 tmp_idx_tgt_ = &indices_tgt;
225 OptimizationFunctorWithIndices functor(
this);
227 bfgs.parameters.sigma = 0.01;
228 bfgs.parameters.rho = 0.01;
229 bfgs.parameters.tau1 = 9;
230 bfgs.parameters.tau2 = 0.05;
231 bfgs.parameters.tau3 = 0.5;
232 bfgs.parameters.order = 3;
234 int inner_iterations_ = 0;
235 int result = bfgs.minimizeInit(x);
239 result = bfgs.minimizeOneStep(x);
243 result = bfgs.testGradient();
246 inner_iterations_ == max_inner_iterations_) {
247 PCL_DEBUG(
"[pcl::registration::TransformationEstimationBFGS::"
248 "estimateRigidTransformation]");
249 PCL_DEBUG(
"BFGS solver finished with exit code %i \n", result);
250 transformation_matrix.setIdentity();
251 applyState(transformation_matrix, x);
255 SolverDidntConvergeException,
256 "[pcl::" << getClassName()
257 <<
"::TransformationEstimationBFGS::estimateRigidTransformation] BFGS "
258 "solver didn't converge!");
261 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
266 Matrix4 transformation_matrix = gicp_->base_transformation_;
267 gicp_->applyState(transformation_matrix, x);
269 int m =
static_cast<int>(gicp_->tmp_idx_src_->size());
270 for (
int i = 0; i < m; ++i) {
273 (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
276 (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
277 Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
281 Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
282 p_trans_src[1] - p_tgt[1],
283 p_trans_src[2] - p_tgt[2]);
284 Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
287 f +=
static_cast<double>(d.transpose() * Md);
292 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
297 Matrix4 transformation_matrix = gicp_->base_transformation_;
298 gicp_->applyState(transformation_matrix, x);
303 Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
304 int m =
static_cast<int>(gicp_->tmp_idx_src_->size());
305 for (
int i = 0; i < m; ++i) {
308 (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
311 (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
313 Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
316 Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
317 p_trans_src[1] - p_tgt[1],
318 p_trans_src[2] - p_tgt[2]);
320 Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
326 p_trans_src = gicp_->base_transformation_.template cast<float>() * p_src;
327 Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
328 dCost_dR_T += p_base_src * Md.transpose();
330 g.head<3>() *= 2.0 / m;
331 dCost_dR_T *= 2.0 / m;
332 gicp_->computeRDerivative(x, dCost_dR_T, g);
335 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
340 Matrix4 transformation_matrix = gicp_->base_transformation_;
341 gicp_->applyState(transformation_matrix, x);
345 Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
346 const int m =
static_cast<int>(gicp_->tmp_idx_src_->size());
347 for (
int i = 0; i < m; ++i) {
350 (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
353 (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
354 Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
357 Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
358 p_trans_src[1] - p_tgt[1],
359 p_trans_src[2] - p_tgt[2]);
361 Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
363 f +=
static_cast<double>(d.transpose() * Md);
368 p_trans_src = gicp_->base_transformation_.template cast<float>() * p_src;
369 Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
371 dCost_dR_T += p_base_src * Md.transpose();
373 f /=
static_cast<double>(m);
374 g.head<3>() *= (2.0 / m);
375 dCost_dR_T *= 2.0 / m;
376 gicp_->computeRDerivative(x, dCost_dR_T, g);
379 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
384 auto translation_epsilon = gicp_->translation_gradient_tolerance_;
385 auto rotation_epsilon = gicp_->rotation_gradient_tolerance_;
387 if ((translation_epsilon < 0.) || (rotation_epsilon < 0.))
391 auto translation_grad = g.head<3>().norm();
394 auto rotation_grad = g.tail<3>().norm();
396 if ((translation_grad < translation_epsilon) && (rotation_grad < rotation_epsilon))
402 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
411 const std::size_t N = indices_->size();
413 mahalanobis_.resize(N, Eigen::Matrix3d::Identity());
415 if ((!target_covariances_) || (target_covariances_->empty())) {
417 computeCovariances<PointTarget>(target_, tree_, *target_covariances_);
420 if ((!input_covariances_) || (input_covariances_->empty())) {
421 input_covariances_.reset(
new MatricesVector);
422 computeCovariances<PointSource>(input_, tree_reciprocal_, *input_covariances_);
425 base_transformation_ = Matrix4::Identity();
428 double dist_threshold = corr_dist_threshold_ * corr_dist_threshold_;
430 std::vector<float> nn_dists(1);
434 while (!converged_) {
440 Eigen::Matrix4d transform_R = Eigen::Matrix4d::Zero();
441 for (std::size_t i = 0; i < 4; i++)
442 for (std::size_t j = 0; j < 4; j++)
443 for (std::size_t k = 0; k < 4; k++)
444 transform_R(i, j) +=
static_cast<double>(transformation_(i, k)) *
445 static_cast<double>(guess(k, j));
447 Eigen::Matrix3d R = transform_R.topLeftCorner<3, 3>();
449 for (std::size_t i = 0; i < N; i++) {
450 PointSource query = output[i];
451 query.getVector4fMap() =
452 transformation_.template cast<float>() * query.getVector4fMap();
454 if (!searchForNeighbors(query, nn_indices, nn_dists)) {
455 PCL_ERROR(
"[pcl::%s::computeTransformation] Unable to find a nearest neighbor "
456 "in the target dataset for point %d in the source!\n",
457 getClassName().c_str(),
464 if (nn_dists[0] < dist_threshold) {
465 Eigen::Matrix3d& C1 = (*input_covariances_)[i];
466 Eigen::Matrix3d& C2 = (*target_covariances_)[nn_indices[0]];
467 Eigen::Matrix3d& M = mahalanobis_[i];
471 Eigen::Matrix3d temp = M * R.transpose();
475 source_indices[cnt] =
static_cast<int>(i);
476 target_indices[cnt] = nn_indices[0];
481 source_indices.resize(cnt);
482 target_indices.resize(cnt);
484 previous_transformation_ = transformation_;
487 rigid_transformation_estimation_(
488 output, source_indices, *target_, target_indices, transformation_);
491 for (
int k = 0; k < 4; k++) {
492 for (
int l = 0; l < 4; l++) {
495 ratio = 1. / rotation_epsilon_;
497 ratio = 1. / transformation_epsilon_;
499 ratio * std::abs(previous_transformation_(k, l) - transformation_(k, l));
504 }
catch (PCLException& e) {
505 PCL_DEBUG(
"[pcl::%s::computeTransformation] Optimization issue %s\n",
506 getClassName().c_str(),
512 if (update_visualizer_ !=
nullptr) {
513 PointCloudSourcePtr input_transformed(
new PointCloudSource);
515 update_visualizer_(*input_transformed, source_indices, *target_, target_indices);
519 if (nr_iterations_ >= max_iterations_ || delta < 1) {
521 PCL_DEBUG(
"[pcl::%s::computeTransformation] Convergence reached. Number of "
522 "iterations: %d out of %d. Transformation difference: %f\n",
523 getClassName().c_str(),
526 (transformation_ - previous_transformation_).array().abs().sum());
527 previous_transformation_ = transformation_;
530 PCL_DEBUG(
"[pcl::%s::computeTransformation] Convergence failed\n",
531 getClassName().c_str());
533 final_transformation_ = previous_transformation_ * guess;
535 PCL_DEBUG(
"Transformation "
536 "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%"
537 "5f\t%5f\t%5f\t%5f\n",
538 final_transformation_(0, 0),
539 final_transformation_(0, 1),
540 final_transformation_(0, 2),
541 final_transformation_(0, 3),
542 final_transformation_(1, 0),
543 final_transformation_(1, 1),
544 final_transformation_(1, 2),
545 final_transformation_(1, 3),
546 final_transformation_(2, 0),
547 final_transformation_(2, 1),
548 final_transformation_(2, 2),
549 final_transformation_(2, 3),
550 final_transformation_(3, 0),
551 final_transformation_(3, 1),
552 final_transformation_(3, 2),
553 final_transformation_(3, 3));
559 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
566 AngleAxis(
static_cast<Scalar
>(x[4]), Vector3::UnitY()) *
567 AngleAxis(
static_cast<Scalar
>(x[3]), Vector3::UnitX()))
569 Matrix4 T = Matrix4::Identity();
570 T.template block<3, 3>(0, 0) = R;
571 T.template block<3, 1>(0, 3) =
Vector3(
572 static_cast<Scalar
>(x[0]),
static_cast<Scalar
>(x[1]),
static_cast<Scalar
>(x[2]));
BFGS stands for Broyden–Fletcher–Goldfarb–Shanno (BFGS) method for solving unconstrained nonlinear op...
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...
typename IterativeClosestPoint< PointSource, PointTarget, Scalar >::Matrix4 Matrix4
void applyState(Matrix4 &t, const Vector6d &x) const
compute transformation matrix from transformation matrix
typename Eigen::Matrix< Scalar, 3, 1 > Vector3
std::vector< Eigen::Matrix3d, Eigen::aligned_allocator< Eigen::Matrix3d > > MatricesVector
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.
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.
typename Eigen::AngleAxis< Scalar > AngleAxis
void computeTransformation(PointCloudSource &output, const Matrix4 &guess) override
Rigid transformation computation method with initial guess.
typename Eigen::Matrix< Scalar, 3, 3 > Matrix3
Eigen::Matrix< double, 6, 1 > Vector6d
An exception that is thrown when the number of correspondents is not equal to the minimum required.
iterator begin() noexcept
shared_ptr< const PointCloud< PointT > > ConstPtr
bool initComputeReciprocal()
Internal computation when reciprocal lookup is needed.
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.
shared_ptr< KdTree< PointT, Tree > > Ptr
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.
@ NegativeGradientEpsilon
const Eigen::Map< const Eigen::Vector4f, Eigen::Aligned > Vector4fMapConst
IndicesAllocator<> Indices
Type used for indices in PCL.
double operator()(const Vector6d &x) override
void df(const Vector6d &x, Vector6d &df) override
BFGSSpace::Status checkGradient(const Vector6d &g) override
void fdf(const Vector6d &x, double &f, Vector6d &df) override
A point structure representing Euclidean xyz coordinates, and the RGB color.