Point Cloud Library (PCL)  1.11.1-dev
sac_model_cone.hpp
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2009-2012, Willow Garage, Inc.
6  * Copyright (c) 2012-, Open Perception, Inc.
7  *
8  * All rights reserved.
9  *
10  * Redistribution and use in source and binary forms, with or without
11  * modification, are permitted provided that the following conditions
12  * are met:
13  *
14  * * Redistributions of source code must retain the above copyright
15  * notice, this list of conditions and the following disclaimer.
16  * * Redistributions in binary form must reproduce the above
17  * copyright notice, this list of conditions and the following
18  * disclaimer in the documentation and/or other materials provided
19  * with the distribution.
20  * * Neither the name of the copyright holder(s) nor the names of its
21  * contributors may be used to endorse or promote products derived
22  * from this software without specific prior written permission.
23  *
24  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
28  * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
29  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
30  * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
31  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
32  * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
33  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
34  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
35  * POSSIBILITY OF SUCH DAMAGE.
36  *
37  */
38 
39 #ifndef PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_CONE_H_
40 #define PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_CONE_H_
41 
42 #include <unsupported/Eigen/NonLinearOptimization> // for LevenbergMarquardt
43 #include <pcl/sample_consensus/sac_model_cone.h>
44 #include <pcl/common/common.h> // for getAngle3D
45 #include <pcl/common/concatenate.h>
46 
47 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
48 template <typename PointT, typename PointNT> bool
50 {
51  if (samples.size () != sample_size_)
52  {
53  PCL_ERROR ("[pcl::SampleConsensusModelCone::isSampleGood] Wrong number of samples (is %lu, should be %lu)!\n", samples.size (), sample_size_);
54  return (false);
55  }
56  return (true);
57 }
58 
59 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
60 template <typename PointT, typename PointNT> bool
62  const Indices &samples, Eigen::VectorXf &model_coefficients) const
63 {
64  // Need 3 samples
65  if (samples.size () != sample_size_)
66  {
67  PCL_ERROR ("[pcl::SampleConsensusModelCone::computeModelCoefficients] Invalid set of samples given (%lu)!\n", samples.size ());
68  return (false);
69  }
70 
71  if (!normals_)
72  {
73  PCL_ERROR ("[pcl::SampleConsensusModelCone::computeModelCoefficients] No input dataset containing normals was given!\n");
74  return (false);
75  }
76 
77  Eigen::Vector4f p1 ((*input_)[samples[0]].x, (*input_)[samples[0]].y, (*input_)[samples[0]].z, 0.0f);
78  Eigen::Vector4f p2 ((*input_)[samples[1]].x, (*input_)[samples[1]].y, (*input_)[samples[1]].z, 0.0f);
79  Eigen::Vector4f p3 ((*input_)[samples[2]].x, (*input_)[samples[2]].y, (*input_)[samples[2]].z, 0.0f);
80 
81  Eigen::Vector4f n1 ((*normals_)[samples[0]].normal[0], (*normals_)[samples[0]].normal[1], (*normals_)[samples[0]].normal[2], 0.0f);
82  Eigen::Vector4f n2 ((*normals_)[samples[1]].normal[0], (*normals_)[samples[1]].normal[1], (*normals_)[samples[1]].normal[2], 0.0f);
83  Eigen::Vector4f n3 ((*normals_)[samples[2]].normal[0], (*normals_)[samples[2]].normal[1], (*normals_)[samples[2]].normal[2], 0.0f);
84 
85  //calculate apex (intersection of the three planes defined by points and belonging normals
86  Eigen::Vector4f ortho12 = n1.cross3(n2);
87  Eigen::Vector4f ortho23 = n2.cross3(n3);
88  Eigen::Vector4f ortho31 = n3.cross3(n1);
89 
90  float denominator = n1.dot(ortho23);
91 
92  float d1 = p1.dot (n1);
93  float d2 = p2.dot (n2);
94  float d3 = p3.dot (n3);
95 
96  Eigen::Vector4f apex = (d1 * ortho23 + d2 * ortho31 + d3 * ortho12) / denominator;
97 
98  //compute axis (normal of plane defined by: { apex+(p1-apex)/(||p1-apex||), apex+(p2-apex)/(||p2-apex||), apex+(p3-apex)/(||p3-apex||)}
99  Eigen::Vector4f ap1 = p1 - apex;
100  Eigen::Vector4f ap2 = p2 - apex;
101  Eigen::Vector4f ap3 = p3 - apex;
102 
103  Eigen::Vector4f np1 = apex + (ap1/ap1.norm ());
104  Eigen::Vector4f np2 = apex + (ap2/ap2.norm ());
105  Eigen::Vector4f np3 = apex + (ap3/ap3.norm ());
106 
107  Eigen::Vector4f np1np2 = np2 - np1;
108  Eigen::Vector4f np1np3 = np3 - np1;
109 
110  Eigen::Vector4f axis_dir = np1np2.cross3 (np1np3);
111  axis_dir.normalize ();
112 
113  // normalize the vector (apex->p) for opening angle calculation
114  ap1.normalize ();
115  ap2.normalize ();
116  ap3.normalize ();
117 
118  //compute opening angle
119  float opening_angle = ( std::acos (ap1.dot (axis_dir)) + std::acos (ap2.dot (axis_dir)) + std::acos (ap3.dot (axis_dir)) ) / 3.0f;
120 
121  model_coefficients.resize (model_size_);
122  // model_coefficients.template head<3> () = line_pt.template head<3> ();
123  model_coefficients[0] = apex[0];
124  model_coefficients[1] = apex[1];
125  model_coefficients[2] = apex[2];
126  // model_coefficients.template segment<3> (3) = line_dir.template head<3> ();
127  model_coefficients[3] = axis_dir[0];
128  model_coefficients[4] = axis_dir[1];
129  model_coefficients[5] = axis_dir[2];
130  // cone radius
131  model_coefficients[6] = opening_angle;
132 
133  if (model_coefficients[6] != -std::numeric_limits<double>::max() && model_coefficients[6] < min_angle_)
134  return (false);
135  if (model_coefficients[6] != std::numeric_limits<double>::max() && model_coefficients[6] > max_angle_)
136  return (false);
137 
138  PCL_DEBUG ("[pcl::SampleConsensusModelCone::computeModelCoefficients] Model is (%g,%g,%g,%g,%g,%g,%g).\n",
139  model_coefficients[0], model_coefficients[1], model_coefficients[2], model_coefficients[3],
140  model_coefficients[4], model_coefficients[5], model_coefficients[6]);
141  return (true);
142 }
143 
144 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
145 template <typename PointT, typename PointNT> void
147  const Eigen::VectorXf &model_coefficients, std::vector<double> &distances) const
148 {
149  // Check if the model is valid given the user constraints
150  if (!isModelValid (model_coefficients))
151  {
152  distances.clear ();
153  return;
154  }
155 
156  distances.resize (indices_->size ());
157 
158  Eigen::Vector4f apex (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
159  Eigen::Vector4f axis_dir (model_coefficients[3], model_coefficients[4], model_coefficients[5], 0.0f);
160  float opening_angle = model_coefficients[6];
161 
162  float apexdotdir = apex.dot (axis_dir);
163  float dirdotdir = 1.0f / axis_dir.dot (axis_dir);
164  // Iterate through the 3d points and calculate the distances from them to the cone
165  for (std::size_t i = 0; i < indices_->size (); ++i)
166  {
167  Eigen::Vector4f pt ((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z, 0.0f);
168 
169  // Calculate the point's projection on the cone axis
170  float k = (pt.dot (axis_dir) - apexdotdir) * dirdotdir;
171  Eigen::Vector4f pt_proj = apex + k * axis_dir;
172 
173  // Calculate the actual radius of the cone at the level of the projected point
174  Eigen::Vector4f height = apex - pt_proj;
175  float actual_cone_radius = tanf (opening_angle) * height.norm ();
176 
177  // Approximate the distance from the point to the cone as the difference between
178  // dist(point,cone_axis) and actual cone radius
179  const double weighted_euclid_dist = (1.0 - normal_distance_weight_) * std::abs (pointToAxisDistance (pt, model_coefficients) - actual_cone_radius);
180 
181  // Calculate the direction of the point from center
182  Eigen::Vector4f dir = pt - pt_proj;
183  dir.normalize ();
184 
185  // Calculate the cones perfect normals
186  height.normalize ();
187  Eigen::Vector4f cone_normal = sinf (opening_angle) * height + std::cos (opening_angle) * dir;
188 
189  // Calculate the angular distance between the point normal and the (dir=pt_proj->pt) vector
190  Eigen::Vector4f n ((*normals_)[(*indices_)[i]].normal[0], (*normals_)[(*indices_)[i]].normal[1], (*normals_)[(*indices_)[i]].normal[2], 0.0f);
191  double d_normal = std::abs (getAngle3D (n, cone_normal));
192  d_normal = (std::min) (d_normal, M_PI - d_normal);
193 
194  distances[i] = std::abs (normal_distance_weight_ * d_normal + weighted_euclid_dist);
195  }
196 }
197 
198 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
199 template <typename PointT, typename PointNT> void
201  const Eigen::VectorXf &model_coefficients, const double threshold, Indices &inliers)
202 {
203  // Check if the model is valid given the user constraints
204  if (!isModelValid (model_coefficients))
205  {
206  inliers.clear ();
207  return;
208  }
209 
210  inliers.clear ();
211  error_sqr_dists_.clear ();
212  inliers.reserve (indices_->size ());
213  error_sqr_dists_.reserve (indices_->size ());
214 
215  Eigen::Vector4f apex (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
216  Eigen::Vector4f axis_dir (model_coefficients[3], model_coefficients[4], model_coefficients[5], 0.0f);
217  float opening_angle = model_coefficients[6];
218 
219  float apexdotdir = apex.dot (axis_dir);
220  float dirdotdir = 1.0f / axis_dir.dot (axis_dir);
221  // Iterate through the 3d points and calculate the distances from them to the cone
222  for (std::size_t i = 0; i < indices_->size (); ++i)
223  {
224  Eigen::Vector4f pt ((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z, 0.0f);
225 
226  // Calculate the point's projection on the cone axis
227  float k = (pt.dot (axis_dir) - apexdotdir) * dirdotdir;
228  Eigen::Vector4f pt_proj = apex + k * axis_dir;
229 
230  // Calculate the actual radius of the cone at the level of the projected point
231  Eigen::Vector4f height = apex - pt_proj;
232  double actual_cone_radius = tan(opening_angle) * height.norm ();
233 
234  // Approximate the distance from the point to the cone as the difference between
235  // dist(point,cone_axis) and actual cone radius
236  const double weighted_euclid_dist = (1.0 - normal_distance_weight_) * std::abs (pointToAxisDistance (pt, model_coefficients) - actual_cone_radius);
237  if (weighted_euclid_dist > threshold) // Early termination: cannot be an inlier
238  continue;
239 
240  // Calculate the direction of the point from center
241  Eigen::Vector4f pp_pt_dir = pt - pt_proj;
242  pp_pt_dir.normalize ();
243 
244  // Calculate the cones perfect normals
245  height.normalize ();
246  Eigen::Vector4f cone_normal = sinf (opening_angle) * height + std::cos (opening_angle) * pp_pt_dir;
247 
248  // Calculate the angular distance between the point normal and the (dir=pt_proj->pt) vector
249  Eigen::Vector4f n ((*normals_)[(*indices_)[i]].normal[0], (*normals_)[(*indices_)[i]].normal[1], (*normals_)[(*indices_)[i]].normal[2], 0.0f);
250  double d_normal = std::abs (getAngle3D (n, cone_normal));
251  d_normal = (std::min) (d_normal, M_PI - d_normal);
252 
253  double distance = std::abs (normal_distance_weight_ * d_normal + weighted_euclid_dist);
254 
255  if (distance < threshold)
256  {
257  // Returns the indices of the points whose distances are smaller than the threshold
258  inliers.push_back ((*indices_)[i]);
259  error_sqr_dists_.push_back (distance);
260  }
261  }
262 }
263 
264 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
265 template <typename PointT, typename PointNT> std::size_t
267  const Eigen::VectorXf &model_coefficients, const double threshold) const
268 {
269 
270  // Check if the model is valid given the user constraints
271  if (!isModelValid (model_coefficients))
272  return (0);
273 
274  std::size_t nr_p = 0;
275 
276  Eigen::Vector4f apex (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
277  Eigen::Vector4f axis_dir (model_coefficients[3], model_coefficients[4], model_coefficients[5], 0.0f);
278  float opening_angle = model_coefficients[6];
279 
280  float apexdotdir = apex.dot (axis_dir);
281  float dirdotdir = 1.0f / axis_dir.dot (axis_dir);
282  // Iterate through the 3d points and calculate the distances from them to the cone
283  for (std::size_t i = 0; i < indices_->size (); ++i)
284  {
285  Eigen::Vector4f pt ((*input_)[(*indices_)[i]].x, (*input_)[(*indices_)[i]].y, (*input_)[(*indices_)[i]].z, 0.0f);
286 
287  // Calculate the point's projection on the cone axis
288  float k = (pt.dot (axis_dir) - apexdotdir) * dirdotdir;
289  Eigen::Vector4f pt_proj = apex + k * axis_dir;
290 
291  // Calculate the actual radius of the cone at the level of the projected point
292  Eigen::Vector4f height = apex - pt_proj;
293  double actual_cone_radius = tan(opening_angle) * height.norm ();
294 
295  // Approximate the distance from the point to the cone as the difference between
296  // dist(point,cone_axis) and actual cone radius
297  const double weighted_euclid_dist = (1.0 - normal_distance_weight_) * std::abs (pointToAxisDistance (pt, model_coefficients) - actual_cone_radius);
298  if (weighted_euclid_dist > threshold) // Early termination: cannot be an inlier
299  continue;
300 
301  // Calculate the direction of the point from center
302  Eigen::Vector4f pp_pt_dir = pt - pt_proj;
303  pp_pt_dir.normalize ();
304 
305  // Calculate the cones perfect normals
306  height.normalize ();
307  Eigen::Vector4f cone_normal = sinf (opening_angle) * height + std::cos (opening_angle) * pp_pt_dir;
308 
309  // Calculate the angular distance between the point normal and the (dir=pt_proj->pt) vector
310  Eigen::Vector4f n ((*normals_)[(*indices_)[i]].normal[0], (*normals_)[(*indices_)[i]].normal[1], (*normals_)[(*indices_)[i]].normal[2], 0.0f);
311  double d_normal = std::abs (getAngle3D (n, cone_normal));
312  d_normal = (std::min) (d_normal, M_PI - d_normal);
313 
314  if (std::abs (normal_distance_weight_ * d_normal + weighted_euclid_dist) < threshold)
315  nr_p++;
316  }
317  return (nr_p);
318 }
319 
320 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
321 template <typename PointT, typename PointNT> void
323  const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const
324 {
325  optimized_coefficients = model_coefficients;
326 
327  // Needs a set of valid model coefficients
328  if (!isModelValid (model_coefficients))
329  {
330  PCL_ERROR ("[pcl::SampleConsensusModelCone::optimizeModelCoefficients] Given model is invalid!\n");
331  return;
332  }
333 
334  // Need more than the minimum sample size to make a difference
335  if (inliers.size () <= sample_size_)
336  {
337  PCL_ERROR ("[pcl::SampleConsensusModelCone:optimizeModelCoefficients] Not enough inliers found to optimize model coefficients (%lu)! Returning the same coefficients.\n", inliers.size ());
338  return;
339  }
340 
341  OptimizationFunctor functor (this, inliers);
342  Eigen::NumericalDiff<OptimizationFunctor > num_diff (functor);
343  Eigen::LevenbergMarquardt<Eigen::NumericalDiff<OptimizationFunctor>, float> lm (num_diff);
344  int info = lm.minimize (optimized_coefficients);
345 
346  // Compute the L2 norm of the residuals
347  PCL_DEBUG ("[pcl::SampleConsensusModelCone::optimizeModelCoefficients] LM solver finished with exit code %i, having a residual norm of %g. \nInitial solution: %g %g %g %g %g %g %g \nFinal solution: %g %g %g %g %g %g %g\n",
348  info, lm.fvec.norm (), model_coefficients[0], model_coefficients[1], model_coefficients[2], model_coefficients[3],
349  model_coefficients[4], model_coefficients[5], model_coefficients[6], optimized_coefficients[0], optimized_coefficients[1], optimized_coefficients[2], optimized_coefficients[3], optimized_coefficients[4], optimized_coefficients[5], optimized_coefficients[6]);
350 
351  Eigen::Vector3f line_dir (optimized_coefficients[3], optimized_coefficients[4], optimized_coefficients[5]);
352  line_dir.normalize ();
353  optimized_coefficients[3] = line_dir[0];
354  optimized_coefficients[4] = line_dir[1];
355  optimized_coefficients[5] = line_dir[2];
356 }
357 
358 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
359 template <typename PointT, typename PointNT> void
361  const Indices &inliers, const Eigen::VectorXf &model_coefficients, PointCloud &projected_points, bool copy_data_fields) const
362 {
363  // Needs a valid set of model coefficients
364  if (!isModelValid (model_coefficients))
365  {
366  PCL_ERROR ("[pcl::SampleConsensusModelCone::projectPoints] Given model is invalid!\n");
367  return;
368  }
369 
370  projected_points.header = input_->header;
371  projected_points.is_dense = input_->is_dense;
372 
373  Eigen::Vector4f apex (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
374  Eigen::Vector4f axis_dir (model_coefficients[3], model_coefficients[4], model_coefficients[5], 0.0f);
375  float opening_angle = model_coefficients[6];
376 
377  float apexdotdir = apex.dot (axis_dir);
378  float dirdotdir = 1.0f / axis_dir.dot (axis_dir);
379 
380  // Copy all the data fields from the input cloud to the projected one?
381  if (copy_data_fields)
382  {
383  // Allocate enough space and copy the basics
384  projected_points.resize (input_->size ());
385  projected_points.width = input_->width;
386  projected_points.height = input_->height;
387 
388  using FieldList = typename pcl::traits::fieldList<PointT>::type;
389  // Iterate over each point
390  for (std::size_t i = 0; i < projected_points.size (); ++i)
391  // Iterate over each dimension
392  pcl::for_each_type <FieldList> (NdConcatenateFunctor <PointT, PointT> ((*input_)[i], projected_points[i]));
393 
394  // Iterate through the 3d points and calculate the distances from them to the cone
395  for (const auto &inlier : inliers)
396  {
397  Eigen::Vector4f pt ((*input_)[inlier].x,
398  (*input_)[inlier].y,
399  (*input_)[inlier].z,
400  1);
401 
402  float k = (pt.dot (axis_dir) - apexdotdir) * dirdotdir;
403 
404  pcl::Vector4fMap pp = projected_points[inlier].getVector4fMap ();
405  pp.matrix () = apex + k * axis_dir;
406 
407  Eigen::Vector4f dir = pt - pp;
408  dir.normalize ();
409 
410  // Calculate the actual radius of the cone at the level of the projected point
411  Eigen::Vector4f height = apex - pp;
412  float actual_cone_radius = tanf (opening_angle) * height.norm ();
413 
414  // Calculate the projection of the point onto the cone
415  pp += dir * actual_cone_radius;
416  }
417  }
418  else
419  {
420  // Allocate enough space and copy the basics
421  projected_points.resize (inliers.size ());
422  projected_points.width = inliers.size ();
423  projected_points.height = 1;
424 
425  using FieldList = typename pcl::traits::fieldList<PointT>::type;
426  // Iterate over each point
427  for (std::size_t i = 0; i < inliers.size (); ++i)
428  // Iterate over each dimension
429  pcl::for_each_type <FieldList> (NdConcatenateFunctor <PointT, PointT> ((*input_)[inliers[i]], projected_points[i]));
430 
431  // Iterate through the 3d points and calculate the distances from them to the cone
432  for (std::size_t i = 0; i < inliers.size (); ++i)
433  {
434  pcl::Vector4fMap pp = projected_points[i].getVector4fMap ();
435  pcl::Vector4fMapConst pt = (*input_)[inliers[i]].getVector4fMap ();
436 
437  float k = (pt.dot (axis_dir) - apexdotdir) * dirdotdir;
438  // Calculate the projection of the point on the line
439  pp.matrix () = apex + k * axis_dir;
440 
441  Eigen::Vector4f dir = pt - pp;
442  dir.normalize ();
443 
444  // Calculate the actual radius of the cone at the level of the projected point
445  Eigen::Vector4f height = apex - pp;
446  float actual_cone_radius = tanf (opening_angle) * height.norm ();
447 
448  // Calculate the projection of the point onto the cone
449  pp += dir * actual_cone_radius;
450  }
451  }
452 }
453 
454 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
455 template <typename PointT, typename PointNT> bool
457  const std::set<index_t> &indices, const Eigen::VectorXf &model_coefficients, const double threshold) const
458 {
459  // Needs a valid model coefficients
460  if (!isModelValid (model_coefficients))
461  {
462  PCL_ERROR ("[pcl::SampleConsensusModelCone::doSamplesVerifyModel] Given model is invalid!\n");
463  return (false);
464  }
465 
466  Eigen::Vector4f apex (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
467  Eigen::Vector4f axis_dir (model_coefficients[3], model_coefficients[4], model_coefficients[5], 0.0f);
468  float openning_angle = model_coefficients[6];
469 
470  float apexdotdir = apex.dot (axis_dir);
471  float dirdotdir = 1.0f / axis_dir.dot (axis_dir);
472 
473  // Iterate through the 3d points and calculate the distances from them to the cone
474  for (const auto &index : indices)
475  {
476  Eigen::Vector4f pt ((*input_)[index].x, (*input_)[index].y, (*input_)[index].z, 0.0f);
477 
478  // Calculate the point's projection on the cone axis
479  float k = (pt.dot (axis_dir) - apexdotdir) * dirdotdir;
480  Eigen::Vector4f pt_proj = apex + k * axis_dir;
481  Eigen::Vector4f dir = pt - pt_proj;
482  dir.normalize ();
483 
484  // Calculate the actual radius of the cone at the level of the projected point
485  Eigen::Vector4f height = apex - pt_proj;
486  double actual_cone_radius = tan (openning_angle) * height.norm ();
487 
488  // Approximate the distance from the point to the cone as the difference between
489  // dist(point,cone_axis) and actual cone radius
490  if (std::abs (static_cast<double>(pointToAxisDistance (pt, model_coefficients) - actual_cone_radius)) > threshold)
491  return (false);
492  }
493 
494  return (true);
495 }
496 
497 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
498 template <typename PointT, typename PointNT> double
500  const Eigen::Vector4f &pt, const Eigen::VectorXf &model_coefficients) const
501 {
502  Eigen::Vector4f apex (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
503  Eigen::Vector4f axis_dir (model_coefficients[3], model_coefficients[4], model_coefficients[5], 0.0f);
504  return sqrt(pcl::sqrPointToLineDistance (pt, apex, axis_dir));
505 }
506 
507 //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
508 template <typename PointT, typename PointNT> bool
509 pcl::SampleConsensusModelCone<PointT, PointNT>::isModelValid (const Eigen::VectorXf &model_coefficients) const
510 {
511  if (!SampleConsensusModel<PointT>::isModelValid (model_coefficients))
512  return (false);
513 
514  // Check against template, if given
515  if (eps_angle_ > 0.0)
516  {
517  // Obtain the cone direction
518  const Eigen::Vector3f coeff(model_coefficients[3], model_coefficients[4], model_coefficients[5]);
519 
520  double angle_diff = std::abs (getAngle3D (axis_, coeff));
521  angle_diff = (std::min) (angle_diff, M_PI - angle_diff);
522  // Check whether the current cone model satisfies our angle threshold criterion with respect to the given axis
523  if (angle_diff > eps_angle_)
524  {
525  PCL_DEBUG ("[pcl::SampleConsensusModelCone::isModelValid] Angle between cone direction and given axis is too large.\n");
526  return (false);
527  }
528  }
529 
530  if (model_coefficients[6] != -std::numeric_limits<double>::max() && model_coefficients[6] < min_angle_)
531  {
532  PCL_DEBUG ("[pcl::SampleConsensusModelCone::isModelValid] The opening angle is too small: should be larger than %g, but is %g.\n",
533  min_angle_, model_coefficients[6]);
534  return (false);
535  }
536  if (model_coefficients[6] != std::numeric_limits<double>::max() && model_coefficients[6] > max_angle_)
537  {
538  PCL_DEBUG ("[pcl::SampleConsensusModelCone::isModelValid] The opening angle is too big: should be smaller than %g, but is %g.\n",
539  max_angle_, model_coefficients[6]);
540  return (false);
541  }
542 
543  return (true);
544 }
545 
546 #define PCL_INSTANTIATE_SampleConsensusModelCone(PointT, PointNT) template class PCL_EXPORTS pcl::SampleConsensusModelCone<PointT, PointNT>;
547 
548 #endif // PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_CONE_H_
549 
pcl::SampleConsensusModelCone::countWithinDistance
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.
Definition: sac_model_cone.hpp:266
pcl::PointCloud::height
std::uint32_t height
The point cloud height (if organized as an image-structure).
Definition: point_cloud.h:394
pcl::SampleConsensusModelCone::isModelValid
bool isModelValid(const Eigen::VectorXf &model_coefficients) const override
Check whether a model is valid given the user constraints.
Definition: sac_model_cone.hpp:509
common.h
pcl::geometry::distance
float distance(const PointT &p1, const PointT &p2)
Definition: geometry.h:60
pcl::SampleConsensusModelCone::doSamplesVerifyModel
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 cone model coefficients.
Definition: sac_model_cone.hpp:456
pcl::NdConcatenateFunctor
Helper functor structure for concatenate.
Definition: concatenate.h:49
pcl::PointCloud< pcl::PointXYZRGB >
pcl::SampleConsensusModelCone::projectPoints
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 cone model.
Definition: sac_model_cone.hpp:360
pcl::PointCloud::width
std::uint32_t width
The point cloud width (if organized as an image-structure).
Definition: point_cloud.h:392
pcl::SampleConsensusModelCone::pointToAxisDistance
double pointToAxisDistance(const Eigen::Vector4f &pt, const Eigen::VectorXf &model_coefficients) const
Get the distance from a point to a line (represented by a point and a direction)
Definition: sac_model_cone.hpp:499
pcl::SampleConsensusModelCone::getDistancesToModel
void getDistancesToModel(const Eigen::VectorXf &model_coefficients, std::vector< double > &distances) const override
Compute all distances from the cloud data to a given cone model.
Definition: sac_model_cone.hpp:146
pcl::getAngle3D
double getAngle3D(const Eigen::Vector4f &v1, const Eigen::Vector4f &v2, const bool in_degree=false)
Compute the smallest angle between two 3D vectors in radians (default) or degree.
Definition: common.hpp:47
M_PI
#define M_PI
Definition: pcl_macros.h:201
pcl::SampleConsensusModelCone::isSampleGood
bool isSampleGood(const Indices &samples) const override
Check if a sample of indices results in a good sample of points indices.
Definition: sac_model_cone.hpp:49
pcl::PointCloud::is_dense
bool is_dense
True if no points are invalid (e.g., have NaN or Inf values in any of their floating point fields).
Definition: point_cloud.h:397
pcl::PointCloud::resize
void resize(std::size_t count)
Resizes the container to contain count elements.
Definition: point_cloud.h:456
pcl::PointCloud::header
pcl::PCLHeader header
The point cloud header.
Definition: point_cloud.h:386
pcl::Indices
IndicesAllocator<> Indices
Type used for indices in PCL.
Definition: types.h:133
pcl::PointCloud::size
std::size_t size() const
Definition: point_cloud.h:437
pcl::SampleConsensusModel
SampleConsensusModel represents the base model class.
Definition: sac_model.h:69
pcl::Vector4fMap
Eigen::Map< Eigen::Vector4f, Eigen::Aligned > Vector4fMap
Definition: point_types.hpp:184
pcl::SampleConsensusModelCone::selectWithinDistance
void selectWithinDistance(const Eigen::VectorXf &model_coefficients, const double threshold, Indices &inliers) override
Select all the points which respect the given model coefficients as inliers.
Definition: sac_model_cone.hpp:200
pcl::SampleConsensusModelCone::computeModelCoefficients
bool computeModelCoefficients(const Indices &samples, Eigen::VectorXf &model_coefficients) const override
Check whether the given index samples can form a valid cone model, compute the model coefficients fro...
Definition: sac_model_cone.hpp:61
pcl::sqrPointToLineDistance
double sqrPointToLineDistance(const Eigen::Vector4f &pt, const Eigen::Vector4f &line_pt, const Eigen::Vector4f &line_dir)
Get the square distance from a point to a line (represented by a point and a direction)
Definition: distances.h:75
pcl::SampleConsensusModelCone::optimizeModelCoefficients
void optimizeModelCoefficients(const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const override
Recompute the cone coefficients using the given inlier set and return them to the user.
Definition: sac_model_cone.hpp:322
pcl::Vector4fMapConst
const Eigen::Map< const Eigen::Vector4f, Eigen::Aligned > Vector4fMapConst
Definition: point_types.hpp:185