14 #include <unsupported/Eigen/NonLinearOptimization>
15 #include <pcl/sample_consensus/sac_model_ellipse3d.h>
16 #include <pcl/common/concatenate.h>
18 #include <Eigen/Eigenvalues>
22 template <
typename Po
intT>
bool
26 if (samples.size () != sample_size_)
28 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::isSampleGood] Wrong number of samples (is %lu, should be %lu)!\n", samples.size (), sample_size_);
33 const Eigen::Vector3d p0 ((*input_)[samples[0]].x, (*input_)[samples[0]].y, (*input_)[samples[0]].z);
34 const Eigen::Vector3d p1 ((*input_)[samples[1]].x, (*input_)[samples[1]].y, (*input_)[samples[1]].z);
35 const Eigen::Vector3d p2 ((*input_)[samples[2]].x, (*input_)[samples[2]].y, (*input_)[samples[2]].z);
40 if ((p1 - p0).cross(p1 - p2).squaredNorm() < Eigen::NumTraits<float>::dummy_precision ())
42 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::isSampleGood] Sample points too similar or collinear!\n");
50 template <
typename Po
intT>
bool
54 if (samples.size () != sample_size_)
56 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Invalid set of samples given (%lu)!\n", samples.size ());
60 model_coefficients.resize (model_size_);
62 const Eigen::Vector3f p0((*input_)[samples[0]].x, (*input_)[samples[0]].y, (*input_)[samples[0]].z);
63 const Eigen::Vector3f p1((*input_)[samples[1]].x, (*input_)[samples[1]].y, (*input_)[samples[1]].z);
64 const Eigen::Vector3f p2((*input_)[samples[2]].x, (*input_)[samples[2]].y, (*input_)[samples[2]].z);
65 const Eigen::Vector3f p3((*input_)[samples[3]].x, (*input_)[samples[3]].y, (*input_)[samples[3]].z);
66 const Eigen::Vector3f p4((*input_)[samples[4]].x, (*input_)[samples[4]].y, (*input_)[samples[4]].z);
67 const Eigen::Vector3f p5((*input_)[samples[5]].x, (*input_)[samples[5]].y, (*input_)[samples[5]].z);
69 const Eigen::Vector3f common_helper_vec = (p1 - p0).cross(p1 - p2);
70 const Eigen::Vector3f ellipse_normal = common_helper_vec.normalized();
74 if (common_helper_vec.squaredNorm() < Eigen::NumTraits<float>::dummy_precision ())
76 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Sample points too similar or collinear!\n");
81 Eigen::Vector3f x_axis = (p1 - p0).normalized();
82 const Eigen::Vector3f z_axis = ellipse_normal.normalized();
83 const Eigen::Vector3f y_axis = z_axis.cross(x_axis).normalized();
86 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
87 << x_axis(0), y_axis(0), z_axis(0),
88 x_axis(1), y_axis(1), z_axis(1),
89 x_axis(2), y_axis(2), z_axis(2))
91 const Eigen::Matrix3f Rot_T = Rot.transpose();
94 const Eigen::Vector3f p0_ = Rot_T * (p0 - p0);
95 const Eigen::Vector3f p1_ = Rot_T * (p1 - p0);
96 const Eigen::Vector3f p2_ = Rot_T * (p2 - p0);
97 const Eigen::Vector3f p3_ = Rot_T * (p3 - p0);
98 const Eigen::Vector3f p4_ = Rot_T * (p4 - p0);
99 const Eigen::Vector3f p5_ = Rot_T * (p5 - p0);
107 const Eigen::VectorXf X = (Eigen::VectorXf(6) << p0_(0), p1_(0), p2_(0), p3_(0), p4_(0), p5_(0)).finished();
108 const Eigen::VectorXf Y = (Eigen::VectorXf(6) << p0_(1), p1_(1), p2_(1), p3_(1), p4_(1), p5_(1)).finished();
111 const Eigen::MatrixXf D = (Eigen::MatrixXf(6,6)
112 << X(0) * X(0), X(0) * Y(0), Y(0) * Y(0), X(0), Y(0), 1.0,
113 X(1) * X(1), X(1) * Y(1), Y(1) * Y(1), X(1), Y(1), 1.0,
114 X(2) * X(2), X(2) * Y(2), Y(2) * Y(2), X(2), Y(2), 1.0,
115 X(3) * X(3), X(3) * Y(3), Y(3) * Y(3), X(3), Y(3), 1.0,
116 X(4) * X(4), X(4) * Y(4), Y(4) * Y(4), X(4), Y(4), 1.0,
117 X(5) * X(5), X(5) * Y(5), Y(5) * Y(5), X(5), Y(5), 1.0)
121 const Eigen::MatrixXf S = D.transpose() * D;
124 const Eigen::MatrixXf C = (Eigen::MatrixXf(6,6)
125 << 0.0, 0.0, -2.0, 0.0, 0.0, 0.0,
126 0.0, 1.0, 0.0, 0.0, 0.0, 0.0,
127 -2.0, 0.0, 0.0, 0.0, 0.0, 0.0,
128 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
129 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
130 0.0, 0.0, 0.0, 0.0, 0.0, 0.0)
134 Eigen::GeneralizedEigenSolver<Eigen::MatrixXf> solver;
135 solver.compute(S, C);
136 const Eigen::VectorXf eigvals = solver.eigenvalues().real();
141 for (
size_t i(0); i < static_cast<size_t>(eigvals.size()); ++i) {
142 if (eigvals(i) < absmin && !std::isinf(eigvals(i))) {
148 PCL_DEBUG(
"[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Failed to find the negative eigenvalue in the GES.\n");
151 const Eigen::VectorXf neigvec = solver.eigenvectors().real().col(idx).normalized();
157 const float con_A(neigvec(0));
158 const float con_B(neigvec(1));
159 const float con_C(neigvec(2));
160 const float con_D(neigvec(3));
161 const float con_E(neigvec(4));
162 const float con_F(neigvec(5));
165 const Eigen::Matrix3f M0 = (Eigen::Matrix3f()
166 << con_F, con_D/2.0, con_E/2.0,
167 con_D/2.0, con_A, con_B/2.0,
168 con_E/2.0, con_B/2.0, con_C)
172 const Eigen::Matrix2f M = (Eigen::Matrix2f()
178 const Eigen::SelfAdjointEigenSolver<Eigen::Matrix2f> solver_M(M, Eigen::EigenvaluesOnly);
180 Eigen::Vector2f eigvals_M = solver_M.eigenvalues();
183 float aux_eigval(0.0);
184 if (std::abs(eigvals_M(0) - con_A) > std::abs(eigvals_M(0) - con_C)) {
185 aux_eigval = eigvals_M(0);
186 eigvals_M(0) = eigvals_M(1);
187 eigvals_M(1) = aux_eigval;
191 float par_a = std::sqrt(-M0.determinant() / (M.determinant() * eigvals_M(0)));
192 float par_b = std::sqrt(-M0.determinant() / (M.determinant() * eigvals_M(1)));
193 const float par_h = (con_B * con_E - 2.0 * con_C * con_D) / (4.0 * con_A * con_C - std::pow(con_B, 2));
194 const float par_k = (con_B * con_D - 2.0 * con_A * con_E) / (4.0 * con_A * con_C - std::pow(con_B, 2));
195 const float par_t = (
M_PI / 2.0 - std::atan((con_A - con_C) / con_B)) / 2.0;
199 Eigen::Vector3f p_ctr;
202 p_ctr = p0 + Rot * Eigen::Vector3f(par_h, par_k, 0.0);
207 p_ctr = p0 + Rot * Eigen::Vector3f(par_k, par_h, 0.0);
211 model_coefficients[0] =
static_cast<float>(p_ctr(0));
212 model_coefficients[1] =
static_cast<float>(p_ctr(1));
213 model_coefficients[2] =
static_cast<float>(p_ctr(2));
216 model_coefficients[3] =
static_cast<float>(par_a);
218 model_coefficients[4] =
static_cast<float>(par_b);
221 model_coefficients[5] =
static_cast<float>(ellipse_normal[0]);
222 model_coefficients[6] =
static_cast<float>(ellipse_normal[1]);
223 model_coefficients[7] =
static_cast<float>(ellipse_normal[2]);
226 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, par_h, par_k, par_t).finished();
227 Eigen::Vector3f p_th_(0.0, 0.0, 0.0);
228 get_ellipse_point(params, par_t, p_th_(0), p_th_(1));
231 x_axis = (Rot * p_th_).normalized();
232 model_coefficients[8] =
static_cast<float>(x_axis[0]);
233 model_coefficients[9] =
static_cast<float>(x_axis[1]);
234 model_coefficients[10] =
static_cast<float>(x_axis[2]);
237 PCL_DEBUG (
"[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Model is (%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g).\n",
238 model_coefficients[0], model_coefficients[1], model_coefficients[2], model_coefficients[3],
239 model_coefficients[4], model_coefficients[5], model_coefficients[6], model_coefficients[7],
240 model_coefficients[8], model_coefficients[9], model_coefficients[10]);
246 template <
typename Po
intT>
void
250 if (!isModelValid (model_coefficients))
255 distances.resize (indices_->size ());
258 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
260 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
262 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
264 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
266 const float par_a(model_coefficients[3]);
268 const float par_b(model_coefficients[4]);
271 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
272 << x_axis(0), y_axis(0), n_axis(0),
273 x_axis(1), y_axis(1), n_axis(1),
274 x_axis(2), y_axis(2), n_axis(2))
276 const Eigen::Matrix3f Rot_T = Rot.transpose();
279 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
283 for (std::size_t i = 0; i < indices_->size (); ++i)
292 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
295 const Eigen::Vector3f p_ = Rot_T * (p - c);
301 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
303 distances[i] = distanceVector.norm();
308 template <
typename Po
intT>
void
310 const Eigen::VectorXf &model_coefficients,
const double threshold,
314 error_sqr_dists_.clear();
316 if (!isModelValid (model_coefficients))
320 inliers.reserve (indices_->size ());
321 error_sqr_dists_.reserve (indices_->size ());
324 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
326 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
328 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
330 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
332 const float par_a(model_coefficients[3]);
334 const float par_b(model_coefficients[4]);
337 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
338 << x_axis(0), y_axis(0), n_axis(0),
339 x_axis(1), y_axis(1), n_axis(1),
340 x_axis(2), y_axis(2), n_axis(2))
342 const Eigen::Matrix3f Rot_T = Rot.transpose();
344 const auto squared_threshold = threshold * threshold;
346 for (std::size_t i = 0; i < indices_->size (); ++i)
349 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
352 const Eigen::Vector3f p_ = Rot_T * (p - c);
358 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
360 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
362 const double sqr_dist = distanceVector.squaredNorm();
363 if (sqr_dist < squared_threshold)
366 inliers.push_back ((*indices_)[i]);
367 error_sqr_dists_.push_back (sqr_dist);
373 template <
typename Po
intT> std::size_t
375 const Eigen::VectorXf &model_coefficients,
const double threshold)
const
378 if (!isModelValid (model_coefficients))
380 std::size_t nr_p = 0;
383 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
385 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
387 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
389 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
391 const float par_a(model_coefficients[3]);
393 const float par_b(model_coefficients[4]);
396 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
397 << x_axis(0), y_axis(0), n_axis(0),
398 x_axis(1), y_axis(1), n_axis(1),
399 x_axis(2), y_axis(2), n_axis(2))
401 const Eigen::Matrix3f Rot_T = Rot.transpose();
403 const auto squared_threshold = threshold * threshold;
405 for (std::size_t i = 0; i < indices_->size (); ++i)
408 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
411 const Eigen::Vector3f p_ = Rot_T * (p - c);
417 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
419 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
421 if (distanceVector.squaredNorm() < squared_threshold)
428 template <
typename Po
intT>
void
431 const Eigen::VectorXf &model_coefficients,
432 Eigen::VectorXf &optimized_coefficients)
const
434 optimized_coefficients = model_coefficients;
437 if (!isModelValid (model_coefficients))
439 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::optimizeModelCoefficients] Given model is invalid!\n");
444 if (inliers.size () <= sample_size_)
446 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::optimizeModelCoefficients] Not enough inliers to refine/optimize the model's coefficients (%lu)! Returning the same coefficients.\n", inliers.size ());
450 OptimizationFunctor functor(
this, inliers);
451 Eigen::NumericalDiff<OptimizationFunctor> num_diff(functor);
452 Eigen::LevenbergMarquardt<Eigen::NumericalDiff<OptimizationFunctor>,
double> lm(num_diff);
453 Eigen::VectorXd coeff = model_coefficients.cast<
double>();
454 int info = lm.minimize(coeff);
455 optimized_coefficients = coeff.cast<
float>();
458 PCL_DEBUG (
"[pcl::SampleConsensusModelEllipse3D::optimizeModelCoefficients] LM solver finished with exit code %i, having a residual norm of %g. \nInitial solution: %g %g %g %g %g %g %g %g %g %g %g\nFinal solution: %g %g %g %g %g %g %g %g %g %g %g\n",
459 info, lm.fvec.norm (),
461 model_coefficients[0],
462 model_coefficients[1],
463 model_coefficients[2],
464 model_coefficients[3],
465 model_coefficients[4],
466 model_coefficients[5],
467 model_coefficients[6],
468 model_coefficients[7],
469 model_coefficients[8],
470 model_coefficients[9],
471 model_coefficients[10],
473 optimized_coefficients[0],
474 optimized_coefficients[1],
475 optimized_coefficients[2],
476 optimized_coefficients[3],
477 optimized_coefficients[4],
478 optimized_coefficients[5],
479 optimized_coefficients[6],
480 optimized_coefficients[7],
481 optimized_coefficients[8],
482 optimized_coefficients[9],
483 optimized_coefficients[10]);
487 template <
typename Po
intT>
void
489 const Indices &inliers,
const Eigen::VectorXf &model_coefficients,
490 PointCloud &projected_points,
bool copy_data_fields)
const
493 if (!isModelValid (model_coefficients))
495 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::projectPoints] Given model is invalid!\n");
499 projected_points.
header = input_->header;
500 projected_points.
is_dense = input_->is_dense;
503 if (copy_data_fields)
506 projected_points.
resize (input_->size ());
507 projected_points.
width = input_->width;
508 projected_points.
height = input_->height;
510 using FieldList =
typename pcl::traits::fieldList<PointT>::type;
512 for (std::size_t i = 0; i < projected_points.
size(); ++i)
519 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
521 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
523 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
525 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
527 const float par_a(model_coefficients[3]);
529 const float par_b(model_coefficients[4]);
532 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
533 << x_axis(0), y_axis(0), n_axis(0),
534 x_axis(1), y_axis(1), n_axis(1),
535 x_axis(2), y_axis(2), n_axis(2))
537 const Eigen::Matrix3f Rot_T = Rot.transpose();
540 for (std::size_t i = 0; i < inliers.size (); ++i)
543 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
546 const Eigen::Vector3f p_ = Rot_T * (p - c);
552 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
554 dvec2ellipse(params, p_(0), p_(1), th_opt);
557 Eigen::Vector3f k_(0.0, 0.0, 0.0);
558 get_ellipse_point(params, th_opt, k_[0], k_[1]);
560 const Eigen::Vector3f k = c + Rot * k_;
562 projected_points[i].x =
static_cast<float> (k[0]);
563 projected_points[i].y =
static_cast<float> (k[1]);
564 projected_points[i].z =
static_cast<float> (k[2]);
570 projected_points.
resize (inliers.size ());
571 projected_points.
width = inliers.size ();
572 projected_points.
height = 1;
574 using FieldList =
typename pcl::traits::fieldList<PointT>::type;
576 for (std::size_t i = 0; i < inliers.size (); ++i)
581 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
583 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
585 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
587 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
589 const float par_a(model_coefficients[3]);
591 const float par_b(model_coefficients[4]);
594 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
595 << x_axis(0), y_axis(0), n_axis(0),
596 x_axis(1), y_axis(1), n_axis(1),
597 x_axis(2), y_axis(2), n_axis(2))
599 const Eigen::Matrix3f Rot_T = Rot.transpose();
602 for (std::size_t i = 0; i < inliers.size (); ++i)
605 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
608 const Eigen::Vector3f p_ = Rot_T * (p - c);
614 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
616 dvec2ellipse(params, p_(0), p_(1), th_opt);
620 Eigen::Vector3f k_(0.0, 0.0, 0.0);
621 get_ellipse_point(params, th_opt, k_[0], k_[1]);
623 const Eigen::Vector3f k = c + Rot * k_;
625 projected_points[i].x =
static_cast<float> (k[0]);
626 projected_points[i].y =
static_cast<float> (k[1]);
627 projected_points[i].z =
static_cast<float> (k[2]);
633 template <
typename Po
intT>
bool
635 const std::set<index_t> &indices,
636 const Eigen::VectorXf &model_coefficients,
637 const double threshold)
const
640 if (!isModelValid (model_coefficients))
642 PCL_ERROR (
"[pcl::SampleConsensusModelEllipse3D::doSamplesVerifyModel] Given model is invalid!\n");
647 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
649 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
651 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
653 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
655 const float par_a(model_coefficients[3]);
657 const float par_b(model_coefficients[4]);
660 const Eigen::Matrix3f Rot = (Eigen::Matrix3f(3,3)
661 << x_axis(0), y_axis(0), n_axis(0),
662 x_axis(1), y_axis(1), n_axis(1),
663 x_axis(2), y_axis(2), n_axis(2))
665 const Eigen::Matrix3f Rot_T = Rot.transpose();
667 const auto squared_threshold = threshold * threshold;
668 for (
const auto &index : indices)
671 const Eigen::Vector3f p((*input_)[index].x, (*input_)[index].y, (*input_)[index].z);
674 const Eigen::Vector3f p_ = Rot_T * (p - c);
680 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
682 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
684 if (distanceVector.squaredNorm() > squared_threshold)
691 template <
typename Po
intT>
bool
697 if (radius_min_ != std::numeric_limits<double>::lowest() && (model_coefficients[3] < radius_min_ || model_coefficients[4] < radius_min_))
699 PCL_DEBUG (
"[pcl::SampleConsensusModelEllipse3D::isModelValid] Semi-minor axis OR semi-major axis (radii) of ellipse is/are too small: should be larger than %g, but are {%g, %g}.\n",
700 radius_min_, model_coefficients[3], model_coefficients[4]);
703 if (radius_max_ != std::numeric_limits<double>::max() && (model_coefficients[3] > radius_max_ || model_coefficients[4] > radius_max_))
705 PCL_DEBUG (
"[pcl::SampleConsensusModelEllipse3D::isModelValid] Semi-minor axis OR semi-major axis (radii) of ellipse is/are too big: should be smaller than %g, but are {%g, %g}.\n",
706 radius_max_, model_coefficients[3], model_coefficients[4]);
716 template <
typename Po
intT>
718 const Eigen::VectorXf& par,
float th,
float& x,
float& y)
725 const float par_a(par[0]);
726 const float par_b(par[1]);
727 const float par_h(par[2]);
728 const float par_k(par[3]);
729 const float par_t(par[4]);
731 x = par_h + std::cos(par_t) * par_a * std::cos(th) -
732 std::sin(par_t) * par_b * std::sin(th);
733 y = par_k + std::sin(par_t) * par_a * std::cos(th) +
734 std::cos(par_t) * par_b * std::sin(th);
740 template <
typename Po
intT>
742 const Eigen::VectorXf& par,
float u,
float v,
float& th_opt)
750 const float par_h = par[2];
751 const float par_k = par[3];
753 const Eigen::Vector2f center(par_h, par_k);
754 Eigen::Vector2f p(u, v);
758 Eigen::Vector2f x_axis(0.0, 0.0);
759 get_ellipse_point(par, 0.0, x_axis(0), x_axis(1));
763 Eigen::Vector2f y_axis(0.0, 0.0);
764 get_ellipse_point(par,
M_PI / 2.0, y_axis(0), y_axis(1));
768 const float x_proj = p.dot(x_axis) / x_axis.norm();
769 const float y_proj = p.dot(y_axis) / y_axis.norm();
773 float th_min(0.0), th_max(0.0);
774 const float th = std::atan2(y_proj, x_proj);
776 if (-
M_PI <= th && th < -
M_PI / 2.0) {
779 th_max = -
M_PI / 2.0;
781 if (-
M_PI / 2.0 <= th && th < 0.0) {
783 th_min = -
M_PI / 2.0;
786 if (0.0 <= th && th <
M_PI / 2.0) {
791 if (
M_PI / 2.0 <= th && th <=
M_PI) {
798 th_opt = golden_section_search(par, u, v, th_min, th_max, 1.e-3);
801 float x(0.0), y(0.0);
802 get_ellipse_point(par, th_opt, x, y);
803 Eigen::Vector2f distanceVector(u - x, v - y);
804 return distanceVector;
808 template <
typename Po
intT>
810 const Eigen::VectorXf& par,
821 constexpr
float phi(1.61803398874989484820f);
824 float tl(th_min), tu(th_max);
825 float ta = tl + (tu - tl) * (1 - 1 / phi);
826 float tb = tl + (tu - tl) * 1 / phi;
828 while ((tu - tl) > epsilon) {
831 float x_a(0.0), y_a(0.0);
832 get_ellipse_point(par, ta, x_a, y_a);
833 float squared_dist_ta = (u - x_a) * (u - x_a) + (v - y_a) * (v - y_a);
836 float x_b(0.0), y_b(0.0);
837 get_ellipse_point(par, tb, x_b, y_b);
838 float squared_dist_tb = (u - x_b) * (u - x_b) + (v - y_b) * (v - y_b);
840 if (squared_dist_ta < squared_dist_tb) {
843 ta = tl + (tu - tl) * (1 - 1 / phi);
845 else if (squared_dist_ta > squared_dist_tb) {
848 tb = tl + (tu - tl) * 1 / phi;
853 ta = tl + (tu - tl) * (1 - 1 / phi);
854 tb = tl + (tu - tl) * 1 / phi;
857 return (tl + tu) / 2.0;
861 #define PCL_INSTANTIATE_SampleConsensusModelEllipse3D(T) template class PCL_EXPORTS pcl::SampleConsensusModelEllipse3D<T>;
PointCloud represents the base class in PCL for storing collections of 3D points.
bool is_dense
True if no points are invalid (e.g., have NaN or Inf values in any of their floating point fields).
void resize(std::size_t count)
Resizes the container to contain count elements.
std::uint32_t width
The point cloud width (if organized as an image-structure).
pcl::PCLHeader header
The point cloud header.
std::uint32_t height
The point cloud height (if organized as an image-structure).
SampleConsensusModelEllipse3D defines a model for 3D ellipse segmentation.
void optimizeModelCoefficients(const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const override
Recompute the 3d ellipse coefficients using the given inlier set and return them to the user.
std::size_t countWithinDistance(const Eigen::VectorXf &model_coefficients, const double threshold) const override
Count all the points which respect the given model coefficients as inliers.
void projectPoints(const Indices &inliers, const Eigen::VectorXf &model_coefficients, PointCloud &projected_points, bool copy_data_fields=true) const override
Create a new point cloud with inliers projected onto the 3d ellipse model.
bool computeModelCoefficients(const Indices &samples, Eigen::VectorXf &model_coefficients) const override
Check whether the given index samples can form a valid 3D ellipse model, compute the model coefficien...
bool doSamplesVerifyModel(const std::set< index_t > &indices, const Eigen::VectorXf &model_coefficients, const double threshold) const override
Verify whether a subset of indices verifies the given 3d ellipse model coefficients.
void selectWithinDistance(const Eigen::VectorXf &model_coefficients, const double threshold, Indices &inliers) override
Compute all distances from the cloud data to a given 3D ellipse model.
bool isSampleGood(const Indices &samples) const override
Check if a sample of indices results in a good sample of points indices.
bool isModelValid(const Eigen::VectorXf &model_coefficients) const override
Check whether a model is valid given the user constraints.
void getDistancesToModel(const Eigen::VectorXf &model_coefficients, std::vector< double > &distances) const override
Compute all distances from the cloud data to a given 3D ellipse model.
SampleConsensusModel represents the base model class.
IndicesAllocator<> Indices
Type used for indices in PCL.
Helper functor structure for concatenate.