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 std::vector<float> nn_dist_sq(k_correspondences_);
69 if (cloud_covariances.size() < cloud->
size())
70 cloud_covariances.resize(cloud->
size());
72 auto matrices_iterator = cloud_covariances.begin();
73 for (
auto points_iterator = cloud->
begin(); points_iterator != cloud->
end();
74 ++points_iterator, ++matrices_iterator) {
75 const PointT& query_point = *points_iterator;
76 Eigen::Matrix3d& cov = *matrices_iterator;
82 kdtree->
nearestKSearch(query_point, k_correspondences_, nn_indices, nn_dist_sq);
85 for (
int j = 0; j < k_correspondences_; j++) {
87 const double ptx = (*cloud)[nn_indices[j]].x - query_point.x,
88 pty = (*cloud)[nn_indices[j]].y - query_point.y,
89 ptz = (*cloud)[nn_indices[j]].z - query_point.z;
95 cov(0, 0) += ptx * ptx;
97 cov(1, 0) += pty * ptx;
98 cov(1, 1) += pty * pty;
100 cov(2, 0) += ptz * ptx;
101 cov(2, 1) += ptz * pty;
102 cov(2, 2) += ptz * ptz;
105 mean /=
static_cast<double>(k_correspondences_);
107 for (
int k = 0; k < 3; k++)
108 for (
int l = 0; l <= k; l++) {
109 cov(k, l) /=
static_cast<double>(k_correspondences_);
110 cov(k, l) -= mean[k] * mean[l];
111 cov(l, k) = cov(k, l);
115 Eigen::JacobiSVD<Eigen::Matrix3d> svd(cov, Eigen::ComputeFullU);
117 Eigen::Matrix3d U = svd.matrixU();
120 for (
int k = 0; k < 3; k++) {
121 Eigen::Vector3d col = U.col(k);
125 cov += v * col * col.transpose();
130 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
135 Eigen::Matrix3d dR_dPhi;
136 Eigen::Matrix3d dR_dTheta;
137 Eigen::Matrix3d dR_dPsi;
139 double phi = x[3], theta = x[4], psi = x[5];
141 double cphi = std::cos(phi), sphi = sin(phi);
142 double ctheta = std::cos(theta), stheta = sin(theta);
143 double cpsi = std::cos(psi), spsi = sin(psi);
149 dR_dPhi(0, 1) = sphi * spsi + cphi * cpsi * stheta;
150 dR_dPhi(1, 1) = -cpsi * sphi + cphi * spsi * stheta;
151 dR_dPhi(2, 1) = cphi * ctheta;
153 dR_dPhi(0, 2) = cphi * spsi - cpsi * sphi * stheta;
154 dR_dPhi(1, 2) = -cphi * cpsi - sphi * spsi * stheta;
155 dR_dPhi(2, 2) = -ctheta * sphi;
157 dR_dTheta(0, 0) = -cpsi * stheta;
158 dR_dTheta(1, 0) = -spsi * stheta;
159 dR_dTheta(2, 0) = -ctheta;
161 dR_dTheta(0, 1) = cpsi * ctheta * sphi;
162 dR_dTheta(1, 1) = ctheta * sphi * spsi;
163 dR_dTheta(2, 1) = -sphi * stheta;
165 dR_dTheta(0, 2) = cphi * cpsi * ctheta;
166 dR_dTheta(1, 2) = cphi * ctheta * spsi;
167 dR_dTheta(2, 2) = -cphi * stheta;
169 dR_dPsi(0, 0) = -ctheta * spsi;
170 dR_dPsi(1, 0) = cpsi * ctheta;
173 dR_dPsi(0, 1) = -cphi * cpsi - sphi * spsi * stheta;
174 dR_dPsi(1, 1) = -cphi * spsi + cpsi * sphi * stheta;
177 dR_dPsi(0, 2) = cpsi * sphi - cphi * spsi * stheta;
178 dR_dPsi(1, 2) = sphi * spsi + cphi * cpsi * stheta;
181 g[3] = matricesInnerProd(dR_dPhi, dCost_dR_T);
182 g[4] = matricesInnerProd(dR_dTheta, dCost_dR_T);
183 g[5] = matricesInnerProd(dR_dPsi, dCost_dR_T);
186 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
193 Matrix4& transformation_matrix)
196 if (indices_src.size() < min_number_correspondences_) {
199 "[pcl::GeneralizedIterativeClosestPoint::estimateRigidTransformationBFGS] Need "
201 << min_number_correspondences_
202 <<
" points to estimate a transform! "
203 "Source and target have "
204 << indices_src.size() <<
" points!");
210 x[0] = transformation_matrix(0, 3);
211 x[1] = transformation_matrix(1, 3);
212 x[2] = transformation_matrix(2, 3);
215 x[3] = std::atan2(transformation_matrix(2, 1), transformation_matrix(2, 2));
216 x[4] = asin(-transformation_matrix(2, 0));
217 x[5] = std::atan2(transformation_matrix(1, 0), transformation_matrix(0, 0));
220 tmp_src_ = &cloud_src;
221 tmp_tgt_ = &cloud_tgt;
222 tmp_idx_src_ = &indices_src;
223 tmp_idx_tgt_ = &indices_tgt;
226 OptimizationFunctorWithIndices functor(
this);
228 bfgs.parameters.sigma = 0.01;
229 bfgs.parameters.rho = 0.01;
230 bfgs.parameters.tau1 = 9;
231 bfgs.parameters.tau2 = 0.05;
232 bfgs.parameters.tau3 = 0.5;
233 bfgs.parameters.order = 3;
235 int inner_iterations_ = 0;
236 int result = bfgs.minimizeInit(x);
240 result = bfgs.minimizeOneStep(x);
244 result = bfgs.testGradient();
247 inner_iterations_ == max_inner_iterations_) {
248 PCL_DEBUG(
"[pcl::registration::TransformationEstimationBFGS::"
249 "estimateRigidTransformation]");
250 PCL_DEBUG(
"BFGS solver finished with exit code %i \n", result);
251 transformation_matrix.setIdentity();
252 applyState(transformation_matrix, x);
256 SolverDidntConvergeException,
257 "[pcl::" << getClassName()
258 <<
"::TransformationEstimationBFGS::estimateRigidTransformation] BFGS "
259 "solver didn't converge!");
262 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
267 Matrix4 transformation_matrix = gicp_->base_transformation_;
268 gicp_->applyState(transformation_matrix, x);
270 int m =
static_cast<int>(gicp_->tmp_idx_src_->size());
271 for (
int i = 0; i < m; ++i) {
274 (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
277 (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
278 Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
282 Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
283 p_trans_src[1] - p_tgt[1],
284 p_trans_src[2] - p_tgt[2]);
285 Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
288 f +=
static_cast<double>(d.transpose() * Md);
293 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
298 Matrix4 transformation_matrix = gicp_->base_transformation_;
299 gicp_->applyState(transformation_matrix, x);
304 Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
305 int m =
static_cast<int>(gicp_->tmp_idx_src_->size());
306 for (
int i = 0; i < m; ++i) {
309 (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
312 (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
314 Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
317 Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
318 p_trans_src[1] - p_tgt[1],
319 p_trans_src[2] - p_tgt[2]);
321 Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
327 p_trans_src = gicp_->base_transformation_.template cast<float>() * p_src;
328 Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
329 dCost_dR_T += p_base_src * Md.transpose();
331 g.head<3>() *= 2.0 / m;
332 dCost_dR_T *= 2.0 / m;
333 gicp_->computeRDerivative(x, dCost_dR_T, g);
336 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
341 Matrix4 transformation_matrix = gicp_->base_transformation_;
342 gicp_->applyState(transformation_matrix, x);
346 Eigen::Matrix3d dCost_dR_T = Eigen::Matrix3d::Zero();
347 const int m =
static_cast<int>(gicp_->tmp_idx_src_->size());
348 for (
int i = 0; i < m; ++i) {
351 (*gicp_->tmp_src_)[(*gicp_->tmp_idx_src_)[i]].getVector4fMap();
354 (*gicp_->tmp_tgt_)[(*gicp_->tmp_idx_tgt_)[i]].getVector4fMap();
355 Eigen::Vector4f p_trans_src(transformation_matrix.template cast<float>() * p_src);
358 Eigen::Vector3d d(p_trans_src[0] - p_tgt[0],
359 p_trans_src[1] - p_tgt[1],
360 p_trans_src[2] - p_tgt[2]);
362 Eigen::Vector3d Md(gicp_->mahalanobis((*gicp_->tmp_idx_src_)[i]) * d);
364 f +=
static_cast<double>(d.transpose() * Md);
369 p_trans_src = gicp_->base_transformation_.template cast<float>() * p_src;
370 Eigen::Vector3d p_base_src(p_trans_src[0], p_trans_src[1], p_trans_src[2]);
372 dCost_dR_T += p_base_src * Md.transpose();
374 f /=
static_cast<double>(m);
375 g.head<3>() *= (2.0 / m);
376 dCost_dR_T *= 2.0 / m;
377 gicp_->computeRDerivative(x, dCost_dR_T, g);
380 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
385 auto translation_epsilon = gicp_->translation_gradient_tolerance_;
386 auto rotation_epsilon = gicp_->rotation_gradient_tolerance_;
388 if ((translation_epsilon < 0.) || (rotation_epsilon < 0.))
392 auto translation_grad = g.head<3>().norm();
395 auto rotation_grad = g.tail<3>().norm();
397 if ((translation_grad < translation_epsilon) && (rotation_grad < rotation_epsilon))
403 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
412 const std::size_t N = indices_->size();
414 mahalanobis_.resize(N, Eigen::Matrix3d::Identity());
416 if ((!target_covariances_) || (target_covariances_->empty())) {
418 computeCovariances<PointTarget>(target_, tree_, *target_covariances_);
421 if ((!input_covariances_) || (input_covariances_->empty())) {
422 input_covariances_.reset(
new MatricesVector);
423 computeCovariances<PointSource>(input_, tree_reciprocal_, *input_covariances_);
426 base_transformation_ = Matrix4::Identity();
429 double dist_threshold = corr_dist_threshold_ * corr_dist_threshold_;
431 std::vector<float> nn_dists(1);
435 while (!converged_) {
441 Eigen::Matrix4d transform_R = Eigen::Matrix4d::Zero();
442 for (std::size_t i = 0; i < 4; i++)
443 for (std::size_t j = 0; j < 4; j++)
444 for (std::size_t k = 0; k < 4; k++)
445 transform_R(i, j) +=
static_cast<double>(transformation_(i, k)) *
446 static_cast<double>(guess(k, j));
448 Eigen::Matrix3d R = transform_R.topLeftCorner<3, 3>();
450 for (std::size_t i = 0; i < N; i++) {
451 PointSource query = output[i];
452 query.getVector4fMap() =
453 transformation_.template cast<float>() * query.getVector4fMap();
455 if (!searchForNeighbors(query, nn_indices, nn_dists)) {
456 PCL_ERROR(
"[pcl::%s::computeTransformation] Unable to find a nearest neighbor "
457 "in the target dataset for point %d in the source!\n",
458 getClassName().c_str(),
465 if (nn_dists[0] < dist_threshold) {
466 Eigen::Matrix3d& C1 = (*input_covariances_)[i];
467 Eigen::Matrix3d& C2 = (*target_covariances_)[nn_indices[0]];
468 Eigen::Matrix3d& M = mahalanobis_[i];
472 Eigen::Matrix3d temp = M * R.transpose();
476 source_indices[cnt] =
static_cast<int>(i);
477 target_indices[cnt] = nn_indices[0];
482 source_indices.resize(cnt);
483 target_indices.resize(cnt);
485 previous_transformation_ = transformation_;
488 rigid_transformation_estimation_(
489 output, source_indices, *target_, target_indices, transformation_);
492 for (
int k = 0; k < 4; k++) {
493 for (
int l = 0; l < 4; l++) {
496 ratio = 1. / rotation_epsilon_;
498 ratio = 1. / transformation_epsilon_;
500 ratio * std::abs(previous_transformation_(k, l) - transformation_(k, l));
505 }
catch (PCLException& e) {
506 PCL_DEBUG(
"[pcl::%s::computeTransformation] Optimization issue %s\n",
507 getClassName().c_str(),
513 if (update_visualizer_ !=
nullptr) {
514 PointCloudSourcePtr input_transformed(
new PointCloudSource);
516 update_visualizer_(*input_transformed, source_indices, *target_, target_indices);
520 if (nr_iterations_ >= max_iterations_ || delta < 1) {
522 PCL_DEBUG(
"[pcl::%s::computeTransformation] Convergence reached. Number of "
523 "iterations: %d out of %d. Transformation difference: %f\n",
524 getClassName().c_str(),
527 (transformation_ - previous_transformation_).array().abs().sum());
528 previous_transformation_ = transformation_;
531 PCL_DEBUG(
"[pcl::%s::computeTransformation] Convergence failed\n",
532 getClassName().c_str());
534 final_transformation_ = previous_transformation_ * guess;
536 PCL_DEBUG(
"Transformation "
537 "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%"
538 "5f\t%5f\t%5f\t%5f\n",
539 final_transformation_(0, 0),
540 final_transformation_(0, 1),
541 final_transformation_(0, 2),
542 final_transformation_(0, 3),
543 final_transformation_(1, 0),
544 final_transformation_(1, 1),
545 final_transformation_(1, 2),
546 final_transformation_(1, 3),
547 final_transformation_(2, 0),
548 final_transformation_(2, 1),
549 final_transformation_(2, 2),
550 final_transformation_(2, 3),
551 final_transformation_(3, 0),
552 final_transformation_(3, 1),
553 final_transformation_(3, 2),
554 final_transformation_(3, 3));
560 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
567 AngleAxis(
static_cast<Scalar
>(x[4]), Vector3::UnitY()) *
568 AngleAxis(
static_cast<Scalar
>(x[3]), Vector3::UnitX()))
570 Matrix4 T = Matrix4::Identity();
571 T.template block<3, 3>(0, 0) = R;
572 T.template block<3, 1>(0, 3) =
Vector3(
573 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.