41 #ifndef PCL_REGISTRATION_NDT_IMPL_H_
42 #define PCL_REGISTRATION_NDT_IMPL_H_
46 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
52 , outlier_ratio_(0.55)
57 reg_name_ =
"NormalDistributionsTransform";
62 const double gauss_d3 = -std::log(gauss_c2);
63 gauss_d1_ = -std::log(gauss_c1 + gauss_c2) - gauss_d3;
65 -2 * std::log((-std::log(gauss_c1 * std::exp(-0.5) + gauss_c2) - gauss_d3) /
72 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
81 const double gauss_c1 = 10 * (1 - outlier_ratio_);
82 const double gauss_c2 = outlier_ratio_ / pow(resolution_, 3);
83 const double gauss_d3 = -std::log(gauss_c2);
84 gauss_d1_ = -std::log(gauss_c1 + gauss_c2) - gauss_d3;
86 -2 * std::log((-std::log(gauss_c1 * std::exp(-0.5) + gauss_c2) - gauss_d3) /
89 if (guess != Matrix4::Identity()) {
91 final_transformation_ = guess;
97 point_jacobian_.setZero();
98 point_jacobian_.block<3, 3>(0, 0).setIdentity();
99 point_hessian_.setZero();
101 Eigen::Transform<Scalar, 3, Eigen::Affine, Eigen::ColMajor> eig_transformation;
102 eig_transformation.matrix() = final_transformation_;
105 Eigen::Matrix<double, 6, 1> transform, score_gradient;
106 Vector3 init_translation = eig_transformation.translation();
107 Vector3 init_rotation = eig_transformation.rotation().eulerAngles(0, 1, 2);
108 transform << init_translation.template cast<double>(),
109 init_rotation.template cast<double>();
111 Eigen::Matrix<double, 6, 6> hessian;
115 double score = computeDerivatives(score_gradient, hessian, output, transform);
117 while (!converged_) {
119 previous_transformation_ = transformation_;
123 Eigen::JacobiSVD<Eigen::Matrix<double, 6, 6>> sv(
124 hessian, Eigen::ComputeFullU | Eigen::ComputeFullV);
126 Eigen::Matrix<double, 6, 1> delta = sv.solve(-score_gradient);
129 double delta_norm = delta.norm();
131 if (delta_norm == 0 || std::isnan(delta_norm)) {
132 trans_likelihood_ = score /
static_cast<double>(input_->size());
133 converged_ = delta_norm == 0;
138 delta_norm = computeStepLengthMT(transform,
142 transformation_epsilon_ / 2,
150 convertTransform(delta, transformation_);
155 if (update_visualizer_)
158 const double cos_angle =
159 0.5 * (transformation_.template block<3, 3>(0, 0).trace() - 1);
160 const double translation_sqr =
161 transformation_.template block<3, 1>(0, 3).squaredNorm();
165 if (nr_iterations_ >= max_iterations_ ||
166 ((transformation_epsilon_ > 0 && translation_sqr <= transformation_epsilon_) &&
167 (transformation_rotation_epsilon_ > 0 &&
168 cos_angle >= transformation_rotation_epsilon_)) ||
169 ((transformation_epsilon_ <= 0) &&
170 (transformation_rotation_epsilon_ > 0 &&
171 cos_angle >= transformation_rotation_epsilon_)) ||
172 ((transformation_epsilon_ > 0 && translation_sqr <= transformation_epsilon_) &&
173 (transformation_rotation_epsilon_ <= 0))) {
181 trans_likelihood_ = score /
static_cast<double>(input_->size());
184 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
187 Eigen::Matrix<double, 6, 1>& score_gradient,
188 Eigen::Matrix<double, 6, 6>& hessian,
190 const Eigen::Matrix<double, 6, 1>& transform,
191 bool compute_hessian)
193 score_gradient.setZero();
198 computeAngleDerivatives(transform);
201 for (std::size_t idx = 0; idx < input_->size(); idx++) {
203 const auto& x_trans_pt = trans_cloud[idx];
207 std::vector<TargetGridLeafConstPtr> neighborhood;
208 std::vector<float> distances;
209 target_cells_.radiusSearch(x_trans_pt, resolution_, neighborhood, distances);
211 for (
const auto& cell : neighborhood) {
213 const auto& x_pt = (*input_)[idx];
214 const Eigen::Vector3d x = x_pt.getVector3fMap().template cast<double>();
217 const Eigen::Vector3d x_trans =
218 x_trans_pt.getVector3fMap().template cast<double>() - cell->getMean();
221 const Eigen::Matrix3d c_inv = cell->getInverseCov();
225 computePointDerivatives(x);
229 updateDerivatives(score_gradient, hessian, x_trans, c_inv, compute_hessian);
235 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
238 const Eigen::Matrix<double, 6, 1>& transform,
bool compute_hessian)
241 const auto calculate_cos_sin = [](
double angle,
double& c,
double& s) {
242 if (std::abs(angle) < 10e-5) {
252 double cx, cy, cz, sx, sy, sz;
253 calculate_cos_sin(transform(3), cx, sx);
254 calculate_cos_sin(transform(4), cy, sy);
255 calculate_cos_sin(transform(5), cz, sz);
259 angular_jacobian_.setZero();
260 angular_jacobian_.row(0).noalias() = Eigen::Vector4d(
261 (-sx * sz + cx * sy * cz), (-sx * cz - cx * sy * sz), (-cx * cy), 1.0);
262 angular_jacobian_.row(1).noalias() = Eigen::Vector4d(
263 (cx * sz + sx * sy * cz), (cx * cz - sx * sy * sz), (-sx * cy), 1.0);
264 angular_jacobian_.row(2).noalias() =
265 Eigen::Vector4d((-sy * cz), sy * sz, cy, 1.0);
266 angular_jacobian_.row(3).noalias() =
267 Eigen::Vector4d(sx * cy * cz, (-sx * cy * sz), sx * sy, 1.0);
268 angular_jacobian_.row(4).noalias() =
269 Eigen::Vector4d((-cx * cy * cz), cx * cy * sz, (-cx * sy), 1.0);
270 angular_jacobian_.row(5).noalias() =
271 Eigen::Vector4d((-cy * sz), (-cy * cz), 0, 1.0);
272 angular_jacobian_.row(6).noalias() =
273 Eigen::Vector4d((cx * cz - sx * sy * sz), (-cx * sz - sx * sy * cz), 0, 1.0);
274 angular_jacobian_.row(7).noalias() =
275 Eigen::Vector4d((sx * cz + cx * sy * sz), (cx * sy * cz - sx * sz), 0, 1.0);
277 if (compute_hessian) {
280 angular_hessian_.setZero();
281 angular_hessian_.row(0).noalias() = Eigen::Vector4d(
282 (-cx * sz - sx * sy * cz), (-cx * cz + sx * sy * sz), sx * cy, 0.0f);
283 angular_hessian_.row(1).noalias() = Eigen::Vector4d(
284 (-sx * sz + cx * sy * cz), (-cx * sy * sz - sx * cz), (-cx * cy), 0.0f);
286 angular_hessian_.row(2).noalias() =
287 Eigen::Vector4d((cx * cy * cz), (-cx * cy * sz), (cx * sy), 0.0f);
288 angular_hessian_.row(3).noalias() =
289 Eigen::Vector4d((sx * cy * cz), (-sx * cy * sz), (sx * sy), 0.0f);
292 angular_hessian_.row(4).noalias() = Eigen::Vector4d(
293 (-sx * cz - cx * sy * sz), (sx * sz - cx * sy * cz), 0, 0.0f);
294 angular_hessian_.row(5).noalias() = Eigen::Vector4d(
295 (cx * cz - sx * sy * sz), (-sx * sy * cz - cx * sz), 0, 0.0f);
297 angular_hessian_.row(6).noalias() =
298 Eigen::Vector4d((-cy * cz), (cy * sz), (-sy), 0.0f);
299 angular_hessian_.row(7).noalias() =
300 Eigen::Vector4d((-sx * sy * cz), (sx * sy * sz), (sx * cy), 0.0f);
301 angular_hessian_.row(8).noalias() =
302 Eigen::Vector4d((cx * sy * cz), (-cx * sy * sz), (-cx * cy), 0.0f);
304 angular_hessian_.row(9).noalias() =
305 Eigen::Vector4d((sy * sz), (sy * cz), 0, 0.0f);
306 angular_hessian_.row(10).noalias() =
307 Eigen::Vector4d((-sx * cy * sz), (-sx * cy * cz), 0, 0.0f);
308 angular_hessian_.row(11).noalias() =
309 Eigen::Vector4d((cx * cy * sz), (cx * cy * cz), 0, 0.0f);
311 angular_hessian_.row(12).noalias() =
312 Eigen::Vector4d((-cy * cz), (cy * sz), 0, 0.0f);
313 angular_hessian_.row(13).noalias() = Eigen::Vector4d(
314 (-cx * sz - sx * sy * cz), (-cx * cz + sx * sy * sz), 0, 0.0f);
315 angular_hessian_.row(14).noalias() = Eigen::Vector4d(
316 (-sx * sz + cx * sy * cz), (-cx * sy * sz - sx * cz), 0, 0.0f);
320 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
323 const Eigen::Vector3d& x,
bool compute_hessian)
328 Eigen::Matrix<double, 8, 1> point_angular_jacobian =
329 angular_jacobian_ * Eigen::Vector4d(x[0], x[1], x[2], 0.0);
330 point_jacobian_(1, 3) = point_angular_jacobian[0];
331 point_jacobian_(2, 3) = point_angular_jacobian[1];
332 point_jacobian_(0, 4) = point_angular_jacobian[2];
333 point_jacobian_(1, 4) = point_angular_jacobian[3];
334 point_jacobian_(2, 4) = point_angular_jacobian[4];
335 point_jacobian_(0, 5) = point_angular_jacobian[5];
336 point_jacobian_(1, 5) = point_angular_jacobian[6];
337 point_jacobian_(2, 5) = point_angular_jacobian[7];
339 if (compute_hessian) {
340 Eigen::Matrix<double, 15, 1> point_angular_hessian =
341 angular_hessian_ * Eigen::Vector4d(x[0], x[1], x[2], 0.0);
344 const Eigen::Vector3d a(0, point_angular_hessian[0], point_angular_hessian[1]);
345 const Eigen::Vector3d b(0, point_angular_hessian[2], point_angular_hessian[3]);
346 const Eigen::Vector3d c(0, point_angular_hessian[4], point_angular_hessian[5]);
347 const Eigen::Vector3d d = point_angular_hessian.block<3, 1>(6, 0);
348 const Eigen::Vector3d e = point_angular_hessian.block<3, 1>(9, 0);
349 const Eigen::Vector3d f = point_angular_hessian.block<3, 1>(12, 0);
354 point_hessian_.block<3, 1>(9, 3) = a;
355 point_hessian_.block<3, 1>(12, 3) = b;
356 point_hessian_.block<3, 1>(15, 3) = c;
357 point_hessian_.block<3, 1>(9, 4) = b;
358 point_hessian_.block<3, 1>(12, 4) = d;
359 point_hessian_.block<3, 1>(15, 4) = e;
360 point_hessian_.block<3, 1>(9, 5) = c;
361 point_hessian_.block<3, 1>(12, 5) = e;
362 point_hessian_.block<3, 1>(15, 5) = f;
366 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
369 Eigen::Matrix<double, 6, 1>& score_gradient,
370 Eigen::Matrix<double, 6, 6>& hessian,
371 const Eigen::Vector3d& x_trans,
372 const Eigen::Matrix3d& c_inv,
373 bool compute_hessian)
const
376 double e_x_cov_x = std::exp(-gauss_d2_ * x_trans.dot(c_inv * x_trans) / 2);
379 const double score_inc = -gauss_d1_ * e_x_cov_x;
381 e_x_cov_x = gauss_d2_ * e_x_cov_x;
384 if (e_x_cov_x > 1 || e_x_cov_x < 0 || std::isnan(e_x_cov_x)) {
389 e_x_cov_x *= gauss_d1_;
391 for (
int i = 0; i < 6; i++) {
394 const Eigen::Vector3d cov_dxd_pi = c_inv * point_jacobian_.col(i);
397 score_gradient(i) += x_trans.dot(cov_dxd_pi) * e_x_cov_x;
399 if (compute_hessian) {
400 for (Eigen::Index j = 0; j < hessian.cols(); j++) {
403 e_x_cov_x * (-gauss_d2_ * x_trans.dot(cov_dxd_pi) *
404 x_trans.dot(c_inv * point_jacobian_.col(j)) +
405 x_trans.dot(c_inv * point_hessian_.block<3, 1>(3 * i, j)) +
406 point_jacobian_.col(j).dot(cov_dxd_pi));
414 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
424 for (std::size_t idx = 0; idx < input_->size(); idx++) {
426 const auto& x_trans_pt = trans_cloud[idx];
430 std::vector<TargetGridLeafConstPtr> neighborhood;
431 std::vector<float> distances;
432 target_cells_.radiusSearch(x_trans_pt, resolution_, neighborhood, distances);
434 for (
const auto& cell : neighborhood) {
436 const auto& x_pt = (*input_)[idx];
437 const Eigen::Vector3d x = x_pt.getVector3fMap().template cast<double>();
440 const Eigen::Vector3d x_trans =
441 x_trans_pt.getVector3fMap().template cast<double>() - cell->getMean();
444 const Eigen::Matrix3d c_inv = cell->getInverseCov();
448 computePointDerivatives(x);
451 updateHessian(hessian, x_trans, c_inv);
456 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
459 Eigen::Matrix<double, 6, 6>& hessian,
460 const Eigen::Vector3d& x_trans,
461 const Eigen::Matrix3d& c_inv)
const
465 gauss_d2_ * std::exp(-gauss_d2_ * x_trans.dot(c_inv * x_trans) / 2);
468 if (e_x_cov_x > 1 || e_x_cov_x < 0 || std::isnan(e_x_cov_x)) {
473 e_x_cov_x *= gauss_d1_;
475 for (
int i = 0; i < 6; i++) {
478 const Eigen::Vector3d cov_dxd_pi = c_inv * point_jacobian_.col(i);
480 for (Eigen::Index j = 0; j < hessian.cols(); j++) {
483 e_x_cov_x * (-gauss_d2_ * x_trans.dot(cov_dxd_pi) *
484 x_trans.dot(c_inv * point_jacobian_.col(j)) +
485 x_trans.dot(c_inv * point_hessian_.block<3, 1>(3 * i, j)) +
486 point_jacobian_.col(j).dot(cov_dxd_pi));
491 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
514 if (g_t * (a_l - a_t) > 0) {
522 if (g_t * (a_l - a_t) < 0) {
536 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
549 if (a_t == a_l && a_t == a_u) {
554 enum class EndpointsCondition { Case1, Case2, Case3, Case4 };
555 EndpointsCondition condition;
558 condition = EndpointsCondition::Case4;
560 else if (f_t > f_l) {
561 condition = EndpointsCondition::Case1;
563 else if (g_t * g_l < 0) {
564 condition = EndpointsCondition::Case2;
566 else if (std::fabs(g_t) <= std::fabs(g_l)) {
567 condition = EndpointsCondition::Case3;
570 condition = EndpointsCondition::Case4;
574 case EndpointsCondition::Case1: {
577 const double z = 3 * (f_t - f_l) / (a_t - a_l) - g_t - g_l;
578 const double w = std::sqrt(z * z - g_t * g_l);
580 const double a_c = a_l + (a_t - a_l) * (w - g_l - z) / (g_t - g_l + 2 * w);
585 a_l - 0.5 * (a_l - a_t) * g_l / (g_l - (f_l - f_t) / (a_l - a_t));
587 if (std::fabs(a_c - a_l) < std::fabs(a_q - a_l)) {
590 return 0.5 * (a_q + a_c);
593 case EndpointsCondition::Case2: {
596 const double z = 3 * (f_t - f_l) / (a_t - a_l) - g_t - g_l;
597 const double w = std::sqrt(z * z - g_t * g_l);
599 const double a_c = a_l + (a_t - a_l) * (w - g_l - z) / (g_t - g_l + 2 * w);
603 const double a_s = a_l - (a_l - a_t) / (g_l - g_t) * g_l;
605 if (std::fabs(a_c - a_t) >= std::fabs(a_s - a_t)) {
611 case EndpointsCondition::Case3: {
614 const double z = 3 * (f_t - f_l) / (a_t - a_l) - g_t - g_l;
615 const double w = std::sqrt(z * z - g_t * g_l);
616 const double a_c = a_l + (a_t - a_l) * (w - g_l - z) / (g_t - g_l + 2 * w);
620 const double a_s = a_l - (a_l - a_t) / (g_l - g_t) * g_l;
624 if (std::fabs(a_c - a_t) < std::fabs(a_s - a_t)) {
632 return std::min(a_t + 0.66 * (a_u - a_t), a_t_next);
634 return std::max(a_t + 0.66 * (a_u - a_t), a_t_next);
638 case EndpointsCondition::Case4: {
641 const double z = 3 * (f_t - f_u) / (a_t - a_u) - g_t - g_u;
642 const double w = std::sqrt(z * z - g_t * g_u);
644 return a_u + (a_t - a_u) * (w - g_u - z) / (g_t - g_u + 2 * w);
649 template <
typename Po
intSource,
typename Po
intTarget,
typename Scalar>
652 const Eigen::Matrix<double, 6, 1>& x,
653 Eigen::Matrix<double, 6, 1>& step_dir,
658 Eigen::Matrix<double, 6, 1>& score_gradient,
659 Eigen::Matrix<double, 6, 6>& hessian,
663 const double phi_0 = -score;
665 double d_phi_0 = -(score_gradient.dot(step_dir));
679 const int max_step_iterations = 10;
680 int step_iterations = 0;
683 const double mu = 1.e-4;
685 const double nu = 0.9;
688 double a_l = 0, a_u = 0;
692 double f_l = auxilaryFunction_PsiMT(a_l, phi_0, phi_0, d_phi_0, mu);
693 double g_l = auxilaryFunction_dPsiMT(d_phi_0, d_phi_0, mu);
695 double f_u = auxilaryFunction_PsiMT(a_u, phi_0, phi_0, d_phi_0, mu);
696 double g_u = auxilaryFunction_dPsiMT(d_phi_0, d_phi_0, mu);
700 bool interval_converged = (step_max - step_min) < 0, open_interval =
true;
702 double a_t = step_init;
703 a_t = std::min(a_t, step_max);
704 a_t = std::max(a_t, step_min);
706 Eigen::Matrix<double, 6, 1> x_t = x + step_dir * a_t;
709 convertTransform(x_t, final_transformation_);
718 score = computeDerivatives(score_gradient, hessian, trans_cloud, x_t,
true);
721 double phi_t = -score;
723 double d_phi_t = -(score_gradient.dot(step_dir));
726 double psi_t = auxilaryFunction_PsiMT(a_t, phi_t, phi_0, d_phi_0, mu);
728 double d_psi_t = auxilaryFunction_dPsiMT(d_phi_t, d_phi_0, mu);
733 while (!interval_converged && step_iterations < max_step_iterations &&
735 d_phi_t <= -nu * d_phi_0 )) {
738 a_t = trialValueSelectionMT(a_l, f_l, g_l, a_u, f_u, g_u, a_t, psi_t, d_psi_t);
741 a_t = trialValueSelectionMT(a_l, f_l, g_l, a_u, f_u, g_u, a_t, phi_t, d_phi_t);
744 a_t = std::min(a_t, step_max);
745 a_t = std::max(a_t, step_min);
747 x_t = x + step_dir * a_t;
750 convertTransform(x_t, final_transformation_);
757 score = computeDerivatives(score_gradient, hessian, trans_cloud, x_t,
false);
762 d_phi_t = -(score_gradient.dot(step_dir));
765 psi_t = auxilaryFunction_PsiMT(a_t, phi_t, phi_0, d_phi_0, mu);
767 d_psi_t = auxilaryFunction_dPsiMT(d_phi_t, d_phi_0, mu);
770 if (open_interval && (psi_t <= 0 && d_psi_t >= 0)) {
771 open_interval =
false;
774 f_l += phi_0 - mu * d_phi_0 * a_l;
778 f_u += phi_0 - mu * d_phi_0 * a_u;
785 updateIntervalMT(a_l, f_l, g_l, a_u, f_u, g_u, a_t, psi_t, d_psi_t);
791 updateIntervalMT(a_l, f_l, g_l, a_u, f_u, g_u, a_t, phi_t, d_phi_t);
800 if (step_iterations) {
801 computeHessian(hessian, trans_cloud);
std::string reg_name_
The registration method name.
int max_iterations_
The maximum number of iterations the internal optimization should run for.
double transformation_epsilon_
The maximum difference between two consecutive transformations in order to consider convergence (user...
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.
IndicesAllocator<> Indices
Type used for indices in PCL.