Point Cloud Library (PCL) 1.15.1
Loading...
Searching...
No Matches
sac_model_ellipse3d.hpp
1/*
2 * SPDX-License-Identifier: BSD-3-Clause
3 *
4 * Point Cloud Library (PCL) - www.pointclouds.org
5 * Copyright (c) 2014-, Open Perception Inc.
6 *
7 * All rights reserved
8 */
9
10#pragma once
11
12#include <limits>
13
14#include <unsupported/Eigen/NonLinearOptimization> // for LevenbergMarquardt
15#include <pcl/sample_consensus/sac_model_ellipse3d.h>
16#include <pcl/common/concatenate.h>
17
18#include <Eigen/Eigenvalues>
19
20
21//////////////////////////////////////////////////////////////////////////
22template <typename PointT> bool
24 const Indices &samples) const
25{
26 if (samples.size () != sample_size_)
27 {
28 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::isSampleGood] Wrong number of samples (is %lu, should be %lu)!\n", samples.size (), sample_size_);
29 return (false);
30 }
31
32 // Use three points out of the 6 samples
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);
36
37 // Check if the squared norm of the cross-product is non-zero, otherwise
38 // common_helper_vec, which plays an important role in computeModelCoefficients,
39 // would likely be ill-formed.
40 if ((p1 - p0).cross(p1 - p2).squaredNorm() < Eigen::NumTraits<float>::dummy_precision ())
41 {
42 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::isSampleGood] Sample points too similar or collinear!\n");
43 return (false);
44 }
45
46 return (true);
47}
48
49//////////////////////////////////////////////////////////////////////////
50template <typename PointT> bool
51pcl::SampleConsensusModelEllipse3D<PointT>::computeModelCoefficients (const Indices &samples, Eigen::VectorXf &model_coefficients) const
52{
53 // Uses 6 samples
54 if (samples.size () != sample_size_)
55 {
56 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Invalid set of samples given (%lu)!\n", samples.size ());
57 return (false);
58 }
59
60 model_coefficients.resize (model_size_); // 11 coefs
61
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);
68
69 const Eigen::Vector3f common_helper_vec = (p1 - p0).cross(p1 - p2);
70 const Eigen::Vector3f ellipse_normal = common_helper_vec.normalized();
71
72 // The same check is implemented in isSampleGood, so be sure to look there too
73 // if you find the need to change something here.
74 if (common_helper_vec.squaredNorm() < Eigen::NumTraits<float>::dummy_precision ())
75 {
76 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Sample points too similar or collinear!\n");
77 return (false);
78 }
79
80 // Definition of the local reference frame of the ellipse
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();
84
85 // Compute the rotation matrix and its transpose
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))
90 .finished();
91 const Eigen::Matrix3f Rot_T = Rot.transpose();
92
93 // Convert the points to local coordinates
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);
100
101
102 // Fit an ellipse to the samples to obtain its conic equation parameters
103 // (this implementation follows the manuscript "Direct Least Square Fitting of Ellipses"
104 // A. Fitzgibbon, M. Pilu and R. Fisher, IEEE TPAMI, 21(5) : 476–480, May 1999).
105
106 // xOy projections only
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();
109
110 // Design matrix D
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)
118 .finished();
119
120 // Scatter matrix S
121 const Eigen::MatrixXf S = D.transpose() * D;
123 // Constraint matrix C
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)
131 .finished();
132
133 // Solve the Generalized Eigensystem: S*a = lambda*C*a
134 Eigen::GeneralizedEigenSolver<Eigen::MatrixXf> solver;
135 solver.compute(S, C);
136 const Eigen::VectorXf eigvals = solver.eigenvalues().real();
137
138 // Find the negative eigenvalue 'neigvec' (the largest, if many exist)
139 int idx(-1);
140 float absmin(0.0);
141 for (size_t i(0); i < static_cast<size_t>(eigvals.size()); ++i) {
142 if (eigvals(i) < absmin && !std::isinf(eigvals(i))) {
143 idx = i;
144 }
145 }
146 // Return "false" in case the negative eigenvalue was not found
147 if (idx == -1) {
148 PCL_DEBUG("[pcl::SampleConsensusModelEllipse3D::computeModelCoefficients] Failed to find the negative eigenvalue in the GES.\n");
149 return (false);
150 }
151 const Eigen::VectorXf neigvec = solver.eigenvectors().real().col(idx).normalized();
153
154 // Convert the conic model parameters to parametric ones
155
156 // Conic parameters
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));
164 // Build matrix M0
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)
169 .finished();
170
171 // Build matrix M
172 const Eigen::Matrix2f M = (Eigen::Matrix2f()
173 << con_A, con_B/2.0,
174 con_B/2.0, con_C)
175 .finished();
176
177 // Calculate the eigenvalues and eigenvectors of matrix M
178 const Eigen::SelfAdjointEigenSolver<Eigen::Matrix2f> solver_M(M, Eigen::EigenvaluesOnly);
179
180 Eigen::Vector2f eigvals_M = solver_M.eigenvalues();
181
182 // Order the eigenvalues so that |lambda_0 - con_A| <= |lambda_0 - con_C|
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;
188 }
189
190 // Parametric parameters of the ellipse
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; // equivalent to: acot((con_A - con_C) / con_B) / 2.0;
197 // Convert the center point of the ellipse to global coordinates
198 // (the if statement ensures that 'par_a' always refers to the semi-major axis length)
199 Eigen::Vector3f p_ctr;
200 float aux_par(0.0);
201 if (par_a > par_b) {
202 p_ctr.noalias() = p0 + Rot * Eigen::Vector3f(par_h, par_k, 0.0);
203 } else {
204 aux_par = par_a;
205 par_a = par_b;
206 par_b = aux_par;
207 p_ctr.noalias() = p0 + Rot * Eigen::Vector3f(par_k, par_h, 0.0);
208 }
209
210 // Center (x, y, z)
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));
214
215 // Semi-major axis length 'a' (along the local x-axis)
216 model_coefficients[3] = static_cast<float>(par_a);
217 // Semi-minor axis length 'b' (along the local y-axis)
218 model_coefficients[4] = static_cast<float>(par_b);
219
220 // Ellipse normal
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]);
224
225 // Retrieve the ellipse point at the tilt angle t (par_t), along the local x-axis
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));
229
230 // Redefinition of the x-axis of the ellipse's local reference frame
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]);
235
236
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]);
241 return (true);
242}
243
244
245//////////////////////////////////////////////////////////////////////////
246template <typename PointT> void
247pcl::SampleConsensusModelEllipse3D<PointT>::getDistancesToModel (const Eigen::VectorXf &model_coefficients, std::vector<double> &distances) const
248{
249 // Check if the model is valid given the user constraints
250 if (!isModelValid (model_coefficients))
251 {
252 distances.clear ();
253 return;
254 }
255 distances.resize (indices_->size ());
256
257 // c : Ellipse Center
258 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
259 // n : Ellipse (Plane) Normal
260 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
261 // x : Ellipse (Plane) X-Axis
262 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
263 // y : Ellipse (Plane) Y-Axis
264 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
265 // a : Ellipse semi-major axis (X) length
266 const float par_a(model_coefficients[3]);
267 // b : Ellipse semi-minor axis (Y) length
268 const float par_b(model_coefficients[4]);
269
270 // Compute the rotation matrix and its transpose
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))
275 .finished();
276 const Eigen::Matrix3f Rot_T = Rot.transpose();
277
278 // Ellipse parameters
279 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
280 float th_opt;
281
282 // Iterate through the 3D points and calculate the distances from them to the ellipse
283 for (std::size_t i = 0; i < indices_->size (); ++i)
284 // Calculate the distance from the point to the ellipse:
285 // 1. calculate intersection point of the plane in which the ellipse lies and the
286 // line from the sample point with the direction of the plane normal (projected point)
287 // 2. calculate the intersection point of the line from the ellipse center to the projected point
288 // with the ellipse
289 // 3. calculate distance from corresponding point on the ellipse to the sample point
290 {
291 // p : Sample Point
292 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
293
294 // Local coordinates of sample point p
295 const Eigen::Vector3f p_ = Rot_T * (p - c);
296
297 // k : Point on Ellipse
298 // Calculate the shortest distance from the point to the ellipse which is given by
299 // the norm of a vector that is normal to the ellipse tangent calculated at the
300 // point it intersects the tangent.
301 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
302
303 distances[i] = distanceVector.norm();
304 }
305}
306
307//////////////////////////////////////////////////////////////////////////
308template <typename PointT> void
310 const Eigen::VectorXf &model_coefficients, const double threshold,
311 Indices &inliers)
312{
313 inliers.clear();
314 error_sqr_dists_.clear();
315 // Check if the model is valid given the user constraints
316 if (!isModelValid (model_coefficients))
317 {
318 return;
319 }
320 inliers.reserve (indices_->size ());
321 error_sqr_dists_.reserve (indices_->size ());
322
323 // c : Ellipse Center
324 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
325 // n : Ellipse (Plane) Normal
326 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
327 // x : Ellipse (Plane) X-Axis
328 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
329 // y : Ellipse (Plane) Y-Axis
330 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
331 // a : Ellipse semi-major axis (X) length
332 const float par_a(model_coefficients[3]);
333 // b : Ellipse semi-minor axis (Y) length
334 const float par_b(model_coefficients[4]);
335
336 // Compute the rotation matrix and its transpose
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))
341 .finished();
342 const Eigen::Matrix3f Rot_T = Rot.transpose();
343
344 const auto squared_threshold = threshold * threshold;
345 // Iterate through the 3d points and calculate the distances from them to the ellipse
346 for (std::size_t i = 0; i < indices_->size (); ++i)
347 {
348 // p : Sample Point
349 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
350
351 // Local coordinates of sample point p
352 const Eigen::Vector3f p_ = Rot_T * (p - c);
353
354 // k : Point on Ellipse
355 // Calculate the shortest distance from the point to the ellipse which is given by
356 // the norm of a vector that is normal to the ellipse tangent calculated at the
357 // point it intersects the tangent.
358 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
359 float th_opt;
360 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
361
362 const double sqr_dist = distanceVector.squaredNorm();
363 if (sqr_dist < squared_threshold)
364 {
365 // Returns the indices of the points whose distances are smaller than the threshold
366 inliers.push_back ((*indices_)[i]);
367 error_sqr_dists_.push_back (sqr_dist);
368 }
369 }
370}
371
372//////////////////////////////////////////////////////////////////////////
373template <typename PointT> std::size_t
375 const Eigen::VectorXf &model_coefficients, const double threshold) const
376{
377 // Check if the model is valid given the user constraints
378 if (!isModelValid (model_coefficients))
379 return (0);
380 std::size_t nr_p = 0;
381
382 // c : Ellipse Center
383 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
384 // n : Ellipse (Plane) Normal
385 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
386 // x : Ellipse (Plane) X-Axis
387 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
388 // y : Ellipse (Plane) Y-Axis
389 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
390 // a : Ellipse semi-major axis (X) length
391 const float par_a(model_coefficients[3]);
392 // b : Ellipse semi-minor axis (Y) length
393 const float par_b(model_coefficients[4]);
394
395 // Compute the rotation matrix and its transpose
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))
400 .finished();
401 const Eigen::Matrix3f Rot_T = Rot.transpose();
402
403 const auto squared_threshold = threshold * threshold;
404 // Iterate through the 3d points and calculate the distances from them to the ellipse
405 for (std::size_t i = 0; i < indices_->size (); ++i)
406 {
407 // p : Sample Point
408 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
409
410 // Local coordinates of sample point p
411 const Eigen::Vector3f p_ = Rot_T * (p - c);
412
413 // k : Point on Ellipse
414 // Calculate the shortest distance from the point to the ellipse which is given by
415 // the norm of a vector that is normal to the ellipse tangent calculated at the
416 // point it intersects the tangent.
417 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
418 float th_opt;
419 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
420
421 if (distanceVector.squaredNorm() < squared_threshold)
422 nr_p++;
423 }
424 return (nr_p);
425}
426
427//////////////////////////////////////////////////////////////////////////
428template <typename PointT> void
430 const Indices &inliers,
431 const Eigen::VectorXf &model_coefficients,
432 Eigen::VectorXf &optimized_coefficients) const
433{
434 optimized_coefficients = model_coefficients;
435
436 // Needs a set of valid model coefficients
437 if (!isModelValid (model_coefficients))
438 {
439 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::optimizeModelCoefficients] Given model is invalid!\n");
440 return;
441 }
442
443 // Need more than the minimum sample size to make a difference
444 if (inliers.size () <= sample_size_)
445 {
446 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::optimizeModelCoefficients] Not enough inliers to refine/optimize the model's coefficients (%lu)! Returning the same coefficients.\n", inliers.size ());
447 return;
448 }
449
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>();
456
457 // Compute the L2 norm of the residuals
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 (),
460
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],
472
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]);
484}
485
486//////////////////////////////////////////////////////////////////////////
487template <typename PointT> void
489 const Indices &inliers, const Eigen::VectorXf &model_coefficients,
490 PointCloud &projected_points, bool copy_data_fields) const
491{
492 // Needs a valid set of model coefficients
493 if (!isModelValid (model_coefficients))
494 {
495 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::projectPoints] Given model is invalid!\n");
496 return;
497 }
498
499 projected_points.header = input_->header;
500 projected_points.is_dense = input_->is_dense;
501
502 // Copy all the data fields from the input cloud to the projected one?
503 if (copy_data_fields)
504 {
505 // Allocate enough space and copy the basics
506 projected_points.resize (input_->size ());
507 projected_points.width = input_->width;
508 projected_points.height = input_->height;
509
510 using FieldList = typename pcl::traits::fieldList<PointT>::type;
511 // Iterate over each point
512 for (std::size_t i = 0; i < projected_points.size(); ++i)
513 {
514 // Iterate over each dimension
516 }
517
518 // c : Ellipse Center
519 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
520 // n : Ellipse (Plane) Normal
521 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
522 // x : Ellipse (Plane) X-Axis
523 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
524 // y : Ellipse (Plane) Y-Axis
525 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
526 // a : Ellipse semi-major axis (X) length
527 const float par_a(model_coefficients[3]);
528 // b : Ellipse semi-minor axis (Y) length
529 const float par_b(model_coefficients[4]);
530
531 // Compute the rotation matrix and its transpose
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))
536 .finished();
537 const Eigen::Matrix3f Rot_T = Rot.transpose();
538
539 // Iterate through the 3d points and calculate the distances from them to the plane
540 for (std::size_t i = 0; i < inliers.size (); ++i)
541 {
542 // p : Sample Point
543 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
544
545 // Local coordinates of sample point p
546 const Eigen::Vector3f p_ = Rot_T * (p - c);
547
548 // k : Point on Ellipse
549 // Calculate the shortest distance from the point to the ellipse which is given by
550 // the norm of a vector that is normal to the ellipse tangent calculated at the
551 // point it intersects the tangent.
552 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
553 float th_opt;
554 dvec2ellipse(params, p_(0), p_(1), th_opt);
555
556 // Retrieve the ellipse point at the tilt angle t, along the local x-axis
557 Eigen::Vector3f k_(0.0, 0.0, 0.0);
558 get_ellipse_point(params, th_opt, k_[0], k_[1]);
559
560 const Eigen::Vector3f k = c + Rot * k_;
561
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]);
565 }
566 }
567 else
568 {
569 // Allocate enough space and copy the basics
570 projected_points.resize (inliers.size ());
571 projected_points.width = inliers.size ();
572 projected_points.height = 1;
573
574 using FieldList = typename pcl::traits::fieldList<PointT>::type;
575 // Iterate over each point
576 for (std::size_t i = 0; i < inliers.size (); ++i)
577 // Iterate over each dimension
578 pcl::for_each_type <FieldList> (NdConcatenateFunctor <PointT, PointT> ((*input_)[inliers[i]], projected_points[i]));
579
580 // c : Ellipse Center
581 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
582 // n : Ellipse (Plane) Normal
583 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
584 // x : Ellipse (Plane) X-Axis
585 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
586 // y : Ellipse (Plane) Y-Axis
587 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
588 // a : Ellipse semi-major axis (X) length
589 const float par_a(model_coefficients[3]);
590 // b : Ellipse semi-minor axis (Y) length
591 const float par_b(model_coefficients[4]);
592
593 // Compute the rotation matrix and its transpose
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))
598 .finished();
599 const Eigen::Matrix3f Rot_T = Rot.transpose();
600
601 // Iterate through the 3d points and calculate the distances from them to the plane
602 for (std::size_t i = 0; i < inliers.size (); ++i)
603 {
604 // p : Sample Point
605 const Eigen::Vector3f p((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z);
606
607 // Local coordinates of sample point p
608 const Eigen::Vector3f p_ = Rot_T * (p - c);
609
610 // k : Point on Ellipse
611 // Calculate the shortest distance from the point to the ellipse which is given by
612 // the norm of a vector that is normal to the ellipse tangent calculated at the
613 // point it intersects the tangent.
614 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
615 float th_opt;
616 dvec2ellipse(params, p_(0), p_(1), th_opt);
617
618 // Retrieve the ellipse point at the tilt angle t, along the local x-axis
619 //// model_coefficients[5] = static_cast<float>(par_t);
620 Eigen::Vector3f k_(0.0, 0.0, 0.0);
621 get_ellipse_point(params, th_opt, k_[0], k_[1]);
622
623 const Eigen::Vector3f k = c + Rot * k_;
624
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]);
628 }
629 }
630}
631
632//////////////////////////////////////////////////////////////////////////
633template <typename PointT> bool
635 const std::set<index_t> &indices,
636 const Eigen::VectorXf &model_coefficients,
637 const double threshold) const
638{
639 // Needs a valid model coefficients
640 if (!isModelValid (model_coefficients))
641 {
642 PCL_ERROR ("[pcl::SampleConsensusModelEllipse3D::doSamplesVerifyModel] Given model is invalid!\n");
643 return (false);
644 }
645
646 // c : Ellipse Center
647 const Eigen::Vector3f c(model_coefficients[0], model_coefficients[1], model_coefficients[2]);
648 // n : Ellipse (Plane) Normal
649 const Eigen::Vector3f n_axis(model_coefficients[5], model_coefficients[6], model_coefficients[7]);
650 // x : Ellipse (Plane) X-Axis
651 const Eigen::Vector3f x_axis(model_coefficients[8], model_coefficients[9], model_coefficients[10]);
652 // y : Ellipse (Plane) Y-Axis
653 const Eigen::Vector3f y_axis = n_axis.cross(x_axis).normalized();
654 // a : Ellipse semi-major axis (X) length
655 const float par_a(model_coefficients[3]);
656 // b : Ellipse semi-minor axis (Y) length
657 const float par_b(model_coefficients[4]);
658
659 // Compute the rotation matrix and its transpose
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))
664 .finished();
665 const Eigen::Matrix3f Rot_T = Rot.transpose();
666
667 const auto squared_threshold = threshold * threshold;
668 for (const auto &index : indices)
669 {
670 // p : Sample Point
671 const Eigen::Vector3f p((*input_)[index].x, (*input_)[index].y, (*input_)[index].z);
672
673 // Local coordinates of sample point p
674 const Eigen::Vector3f p_ = Rot_T * (p - c);
675
676 // k : Point on Ellipse
677 // Calculate the shortest distance from the point to the ellipse which is given by
678 // the norm of a vector that is normal to the ellipse tangent calculated at the
679 // point it intersects the tangent.
680 const Eigen::VectorXf params = (Eigen::VectorXf(5) << par_a, par_b, 0.0, 0.0, 0.0).finished();
681 float th_opt;
682 const Eigen::Vector2f distanceVector = dvec2ellipse(params, p_(0), p_(1), th_opt);
683
684 if (distanceVector.squaredNorm() > squared_threshold)
685 return (false);
686 }
687 return (true);
688}
689
690//////////////////////////////////////////////////////////////////////////
691template <typename PointT> bool
692pcl::SampleConsensusModelEllipse3D<PointT>::isModelValid (const Eigen::VectorXf &model_coefficients) const
693{
694 if (!SampleConsensusModel<PointT>::isModelValid (model_coefficients))
695 return (false);
696
697 if (radius_min_ != std::numeric_limits<double>::lowest() && (model_coefficients[3] < radius_min_ || model_coefficients[4] < radius_min_))
698 {
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]);
701 return (false);
702 }
703 if (radius_max_ != std::numeric_limits<double>::max() && (model_coefficients[3] > radius_max_ || model_coefficients[4] > radius_max_))
704 {
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]);
707 return (false);
708 }
709
710 return (true);
711}
712
713
714
715//////////////////////////////////////////////////////////////////////////
716template <typename PointT>
717void inline pcl::SampleConsensusModelEllipse3D<PointT>::get_ellipse_point(
718 const Eigen::VectorXf& par, float th, float& x, float& y)
719{
720 /*
721 * Calculates a point on the ellipse model 'par' using the angle 'th'.
722 */
723
724 // Parametric eq.params
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]);
730
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);
735
736 return;
737}
738
739//////////////////////////////////////////////////////////////////////////
740template <typename PointT>
741Eigen::Vector2f inline pcl::SampleConsensusModelEllipse3D<PointT>::dvec2ellipse(
742 const Eigen::VectorXf& par, float u, float v, float& th_opt)
743{
744 /*
745 * Minimum distance vector from point p=(u,v) to the ellipse model 'par'.
746 */
747
748 // Parametric eq.params
749 // (par_a, par_b, and par_t do not need to be declared)
750 const float par_h = par[2];
751 const float par_k = par[3];
752
753 const Eigen::Vector2f center(par_h, par_k);
754 Eigen::Vector2f p(u, v);
755 p -= center;
756
757 // Local x-axis of the ellipse
758 Eigen::Vector2f x_axis(0.0, 0.0);
759 get_ellipse_point(par, 0.0, x_axis(0), x_axis(1));
760 x_axis -= center;
761
762 // Local y-axis of the ellipse
763 Eigen::Vector2f y_axis(0.0, 0.0);
764 get_ellipse_point(par, M_PI / 2.0, y_axis(0), y_axis(1));
765 y_axis -= center;
766
767 // Convert the point p=(u,v) to local ellipse coordinates
768 const float x_proj = p.dot(x_axis) / x_axis.norm();
769 const float y_proj = p.dot(y_axis) / y_axis.norm();
770
771 // Find the ellipse quandrant to where the point p=(u,v) belongs,
772 // and limit the search interval to 'th_min' and 'th_max'.
773 float th_min(0.0), th_max(0.0);
774 const float th = std::atan2(y_proj, x_proj);
775
776 if (-M_PI <= th && th < -M_PI / 2.0) {
777 // Q3
778 th_min = -M_PI;
779 th_max = -M_PI / 2.0;
780 }
781 if (-M_PI / 2.0 <= th && th < 0.0) {
782 // Q4
783 th_min = -M_PI / 2.0;
784 th_max = 0.0;
785 }
786 if (0.0 <= th && th < M_PI / 2.0) {
787 // Q1
788 th_min = 0.0;
789 th_max = M_PI / 2.0;
790 }
791 if (M_PI / 2.0 <= th && th <= M_PI) {
792 // Q2
793 th_min = M_PI / 2.0;
794 th_max = M_PI;
795 }
796
797 // Use an unconstrained line search optimizer to find the optimal th_opt
798 th_opt = golden_section_search(par, u, v, th_min, th_max, 1.e-3);
799
800 // Distance vector from a point (u,v) to a given point in the ellipse model 'par' at an angle 'th_opt'.
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;
805}
806
807//////////////////////////////////////////////////////////////////////////
808template <typename PointT>
809float inline pcl::SampleConsensusModelEllipse3D<PointT>::golden_section_search(
810 const Eigen::VectorXf& par,
811 float u,
812 float v,
813 float th_min,
814 float th_max,
815 float epsilon)
816{
817 /*
818 * Golden section search
819 */
820
821 constexpr float phi(1.61803398874989484820f); // Golden ratio
822
823 // tl (theta lower bound), tu (theta upper bound)
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;
827
828 while ((tu - tl) > epsilon) {
829
830 // theta a
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);
834
835 // theta b
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);
839
840 if (squared_dist_ta < squared_dist_tb) {
841 tu = tb;
842 tb = ta;
843 ta = tl + (tu - tl) * (1 - 1 / phi);
844 }
845 else if (squared_dist_ta > squared_dist_tb) {
846 tl = ta;
847 ta = tb;
848 tb = tl + (tu - tl) * 1 / phi;
849 }
850 else {
851 tl = ta;
852 tu = tb;
853 ta = tl + (tu - tl) * (1 - 1 / phi);
854 tb = tl + (tu - tl) * 1 / phi;
855 }
856 }
857 return (tl + tu) / 2.0;
858}
859
860
861#define PCL_INSTANTIATE_SampleConsensusModelEllipse3D(T) template class PCL_EXPORTS pcl::SampleConsensusModelEllipse3D<T>;
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...
typename SampleConsensusModel< PointT >::PointCloud PointCloud
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.
double radius_min_
The minimum and maximum radius limits for the model.
Definition sac_model.h:565
unsigned int sample_size_
The size of a sample from which the model is computed.
Definition sac_model.h:589
IndicesPtr indices_
A pointer to the vector of point indices to use.
Definition sac_model.h:557
PointCloudConstPtr input_
A boost shared pointer to the point cloud data array.
Definition sac_model.h:554
virtual bool isModelValid(const Eigen::VectorXf &model_coefficients) const
Check whether a model is valid given the user constraints.
Definition sac_model.h:528
unsigned int model_size_
The number of coefficients in the model.
Definition sac_model.h:592
std::vector< double > error_sqr_dists_
A vector holding the distances to the computed model.
Definition sac_model.h:586
void for_each_type(F f)
IndicesAllocator<> Indices
Type used for indices in PCL.
Definition types.h:133
#define M_PI
Definition pcl_macros.h:203
Helper functor structure for concatenate.
Definition concatenate.h:50