Point Cloud Library (PCL)  1.11.1-dev
sac_model_sphere.hpp
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2009, 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  * $Id$
38  *
39  */
40 
41 #ifndef PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_SPHERE_H_
42 #define PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_SPHERE_H_
43 
44 #include <pcl/sample_consensus/eigen.h>
45 #include <pcl/sample_consensus/sac_model_sphere.h>
46 
47 //////////////////////////////////////////////////////////////////////////
48 template <typename PointT> bool
50 {
51  if (samples.size () != sample_size_)
52  {
53  PCL_ERROR ("[pcl::SampleConsensusModelSphere::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> bool
62  const Indices &samples, Eigen::VectorXf &model_coefficients) const
63 {
64  // Need 4 samples
65  if (samples.size () != sample_size_)
66  {
67  PCL_ERROR ("[pcl::SampleConsensusModelSphere::computeModelCoefficients] Invalid set of samples given (%lu)!\n", samples.size ());
68  return (false);
69  }
70 
71  Eigen::Matrix4f temp;
72  for (int i = 0; i < 4; i++)
73  {
74  temp (i, 0) = (*input_)[samples[i]].x;
75  temp (i, 1) = (*input_)[samples[i]].y;
76  temp (i, 2) = (*input_)[samples[i]].z;
77  temp (i, 3) = 1;
78  }
79  float m11 = temp.determinant ();
80  if (m11 == 0)
81  {
82  return (false); // the points don't define a sphere!
83  }
84 
85  for (int i = 0; i < 4; ++i)
86  {
87  temp (i, 0) = ((*input_)[samples[i]].x) * ((*input_)[samples[i]].x) +
88  ((*input_)[samples[i]].y) * ((*input_)[samples[i]].y) +
89  ((*input_)[samples[i]].z) * ((*input_)[samples[i]].z);
90  }
91  float m12 = temp.determinant ();
92 
93  for (int i = 0; i < 4; ++i)
94  {
95  temp (i, 1) = temp (i, 0);
96  temp (i, 0) = (*input_)[samples[i]].x;
97  }
98  float m13 = temp.determinant ();
99 
100  for (int i = 0; i < 4; ++i)
101  {
102  temp (i, 2) = temp (i, 1);
103  temp (i, 1) = (*input_)[samples[i]].y;
104  }
105  float m14 = temp.determinant ();
106 
107  for (int i = 0; i < 4; ++i)
108  {
109  temp (i, 0) = temp (i, 2);
110  temp (i, 1) = (*input_)[samples[i]].x;
111  temp (i, 2) = (*input_)[samples[i]].y;
112  temp (i, 3) = (*input_)[samples[i]].z;
113  }
114  float m15 = temp.determinant ();
115 
116  // Center (x , y, z)
117  model_coefficients.resize (model_size_);
118  model_coefficients[0] = 0.5f * m12 / m11;
119  model_coefficients[1] = 0.5f * m13 / m11;
120  model_coefficients[2] = 0.5f * m14 / m11;
121  // Radius
122  model_coefficients[3] = std::sqrt (model_coefficients[0] * model_coefficients[0] +
123  model_coefficients[1] * model_coefficients[1] +
124  model_coefficients[2] * model_coefficients[2] - m15 / m11);
125 
126  return (true);
127 }
128 
129 #define AT(POS) ((*input_)[(*indices_)[(POS)]])
130 
131 #ifdef __AVX__
132 // This function computes the squared distances (i.e. the distances without the square root) of 8 points to the center of the sphere
133 template <typename PointT> inline __m256 pcl::SampleConsensusModelSphere<PointT>::sqr_dist8 (const std::size_t i, const __m256 a_vec, const __m256 b_vec, const __m256 c_vec) const
134 {
135  const __m256 tmp1 = _mm256_sub_ps (_mm256_set_ps (AT(i ).x, AT(i+1).x, AT(i+2).x, AT(i+3).x, AT(i+4).x, AT(i+5).x, AT(i+6).x, AT(i+7).x), a_vec);
136  const __m256 tmp2 = _mm256_sub_ps (_mm256_set_ps (AT(i ).y, AT(i+1).y, AT(i+2).y, AT(i+3).y, AT(i+4).y, AT(i+5).y, AT(i+6).y, AT(i+7).y), b_vec);
137  const __m256 tmp3 = _mm256_sub_ps (_mm256_set_ps (AT(i ).z, AT(i+1).z, AT(i+2).z, AT(i+3).z, AT(i+4).z, AT(i+5).z, AT(i+6).z, AT(i+7).z), c_vec);
138  return _mm256_add_ps (_mm256_add_ps (_mm256_mul_ps (tmp1, tmp1), _mm256_mul_ps (tmp2, tmp2)), _mm256_mul_ps(tmp3, tmp3));
139 }
140 #endif // ifdef __AVX__
141 
142 #ifdef __SSE__
143 // This function computes the squared distances (i.e. the distances without the square root) of 4 points to the center of the sphere
144 template <typename PointT> inline __m128 pcl::SampleConsensusModelSphere<PointT>::sqr_dist4 (const std::size_t i, const __m128 a_vec, const __m128 b_vec, const __m128 c_vec) const
145 {
146  const __m128 tmp1 = _mm_sub_ps (_mm_set_ps (AT(i ).x, AT(i+1).x, AT(i+2).x, AT(i+3).x), a_vec);
147  const __m128 tmp2 = _mm_sub_ps (_mm_set_ps (AT(i ).y, AT(i+1).y, AT(i+2).y, AT(i+3).y), b_vec);
148  const __m128 tmp3 = _mm_sub_ps (_mm_set_ps (AT(i ).z, AT(i+1).z, AT(i+2).z, AT(i+3).z), c_vec);
149  return _mm_add_ps (_mm_add_ps (_mm_mul_ps (tmp1, tmp1), _mm_mul_ps (tmp2, tmp2)), _mm_mul_ps(tmp3, tmp3));
150 }
151 #endif // ifdef __SSE__
152 
153 #undef AT
154 
155 //////////////////////////////////////////////////////////////////////////
156 template <typename PointT> void
158  const Eigen::VectorXf &model_coefficients, std::vector<double> &distances) const
159 {
160  // Check if the model is valid given the user constraints
161  if (!isModelValid (model_coefficients))
162  {
163  distances.clear ();
164  return;
165  }
166  distances.resize (indices_->size ());
167 
168  const Eigen::Vector3f center (model_coefficients[0], model_coefficients[1], model_coefficients[2]);
169  // Iterate through the 3d points and calculate the distances from them to the sphere
170  for (std::size_t i = 0; i < indices_->size (); ++i)
171  {
172  // Calculate the distance from the point to the sphere as the difference between
173  //dist(point,sphere_origin) and sphere_radius
174  distances[i] = std::abs (((*input_)[(*indices_)[i]].getVector3fMap () - center).norm () - model_coefficients[3]);
175  }
176 }
177 
178 //////////////////////////////////////////////////////////////////////////
179 template <typename PointT> void
181  const Eigen::VectorXf &model_coefficients, const double threshold, Indices &inliers)
182 {
183  // Check if the model is valid given the user constraints
184  if (!isModelValid (model_coefficients))
185  {
186  inliers.clear ();
187  return;
188  }
189 
190  inliers.clear ();
191  error_sqr_dists_.clear ();
192  inliers.reserve (indices_->size ());
193  error_sqr_dists_.reserve (indices_->size ());
194 
195  const float sqr_inner_radius = (model_coefficients[3] <= threshold ? 0.0f : (model_coefficients[3] - threshold) * (model_coefficients[3] - threshold));
196  const float sqr_outer_radius = (model_coefficients[3] + threshold) * (model_coefficients[3] + threshold);
197  const Eigen::Vector3f center (model_coefficients[0], model_coefficients[1], model_coefficients[2]);
198  // Iterate through the 3d points and calculate the distances from them to the sphere
199  for (std::size_t i = 0; i < indices_->size (); ++i)
200  {
201  // To avoid sqrt computation: consider one larger sphere (radius + threshold) and one smaller sphere (radius - threshold).
202  // Valid if point is in larger sphere, but not in smaller sphere.
203  const float sqr_dist = ((*input_)[(*indices_)[i]].getVector3fMap () - center).squaredNorm ();
204  if ((sqr_dist <= sqr_outer_radius) && (sqr_dist >= sqr_inner_radius))
205  {
206  // Returns the indices of the points whose distances are smaller than the threshold
207  inliers.push_back ((*indices_)[i]);
208  // Only compute exact distance if necessary (if point is inlier)
209  error_sqr_dists_.push_back (static_cast<double> (std::abs (std::sqrt (sqr_dist) - model_coefficients[3])));
210  }
211  }
212 }
213 
214 //////////////////////////////////////////////////////////////////////////
215 template <typename PointT> std::size_t
217  const Eigen::VectorXf &model_coefficients, const double threshold) const
218 {
219  // Check if the model is valid given the user constraints
220  if (!isModelValid (model_coefficients))
221  return (0);
222 
223 #if defined (__AVX__) && defined (__AVX2__)
224  return countWithinDistanceAVX (model_coefficients, threshold);
225 #elif defined (__SSE__) && defined (__SSE2__) && defined (__SSE4_1__)
226  return countWithinDistanceSSE (model_coefficients, threshold);
227 #else
228  return countWithinDistanceStandard (model_coefficients, threshold);
229 #endif
230 }
231 
232 //////////////////////////////////////////////////////////////////////////
233 template <typename PointT> std::size_t
235  const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i) const
236 {
237  std::size_t nr_p = 0;
238  const float sqr_inner_radius = (model_coefficients[3] <= threshold ? 0.0f : (model_coefficients[3] - threshold) * (model_coefficients[3] - threshold));
239  const float sqr_outer_radius = (model_coefficients[3] + threshold) * (model_coefficients[3] + threshold);
240  const Eigen::Vector3f center (model_coefficients[0], model_coefficients[1], model_coefficients[2]);
241  // Iterate through the 3d points and calculate the distances from them to the sphere
242  for (; i < indices_->size (); ++i)
243  {
244  // To avoid sqrt computation: consider one larger sphere (radius + threshold) and one smaller sphere (radius - threshold).
245  // Valid if point is in larger sphere, but not in smaller sphere.
246  const float sqr_dist = ((*input_)[(*indices_)[i]].getVector3fMap () - center).squaredNorm ();
247  if ((sqr_dist <= sqr_outer_radius) && (sqr_dist >= sqr_inner_radius))
248  nr_p++;
249  }
250  return (nr_p);
251 }
252 
253 //////////////////////////////////////////////////////////////////////////
254 #if defined (__SSE__) && defined (__SSE2__) && defined (__SSE4_1__)
255 template <typename PointT> std::size_t
257  const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i) const
258 {
259  std::size_t nr_p = 0;
260  const __m128 a_vec = _mm_set1_ps (model_coefficients[0]);
261  const __m128 b_vec = _mm_set1_ps (model_coefficients[1]);
262  const __m128 c_vec = _mm_set1_ps (model_coefficients[2]);
263  // To avoid sqrt computation: consider one larger sphere (radius + threshold) and one smaller sphere (radius - threshold). Valid if point is in larger sphere, but not in smaller sphere.
264  const __m128 sqr_inner_radius = _mm_set1_ps ((model_coefficients[3] <= threshold ? 0.0 : (model_coefficients[3]-threshold)*(model_coefficients[3]-threshold)));
265  const __m128 sqr_outer_radius = _mm_set1_ps ((model_coefficients[3]+threshold)*(model_coefficients[3]+threshold));
266  __m128i res = _mm_set1_epi32(0); // This corresponds to nr_p: 4 32bit integers that, summed together, hold the number of inliers
267  for (; (i + 4) <= indices_->size (); i += 4)
268  {
269  const __m128 sqr_dist = sqr_dist4 (i, a_vec, b_vec, c_vec);
270  const __m128 mask = _mm_and_ps (_mm_cmplt_ps (sqr_inner_radius, sqr_dist), _mm_cmplt_ps (sqr_dist, sqr_outer_radius)); // The mask contains 1 bits if the corresponding points are inliers, else 0 bits
271  res = _mm_add_epi32 (res, _mm_and_si128 (_mm_set1_epi32 (1), _mm_castps_si128 (mask))); // The latter part creates a vector with ones (as 32bit integers) where the points are inliers
272  //const int res = _mm_movemask_ps (mask);
273  //if (res & 1) nr_p++;
274  //if (res & 2) nr_p++;
275  //if (res & 4) nr_p++;
276  //if (res & 8) nr_p++;
277  }
278  nr_p += _mm_extract_epi32 (res, 0);
279  nr_p += _mm_extract_epi32 (res, 1);
280  nr_p += _mm_extract_epi32 (res, 2);
281  nr_p += _mm_extract_epi32 (res, 3);
282 
283  // Process the remaining points (at most 3)
284  nr_p += countWithinDistanceStandard (model_coefficients, threshold, i);
285  return (nr_p);
286 }
287 #endif
288 
289 //////////////////////////////////////////////////////////////////////////
290 #if defined (__AVX__) && defined (__AVX2__)
291 template <typename PointT> std::size_t
293  const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i) const
294 {
295  std::size_t nr_p = 0;
296  const __m256 a_vec = _mm256_set1_ps (model_coefficients[0]);
297  const __m256 b_vec = _mm256_set1_ps (model_coefficients[1]);
298  const __m256 c_vec = _mm256_set1_ps (model_coefficients[2]);
299  // To avoid sqrt computation: consider one larger sphere (radius + threshold) and one smaller sphere (radius - threshold). Valid if point is in larger sphere, but not in smaller sphere.
300  const __m256 sqr_inner_radius = _mm256_set1_ps ((model_coefficients[3] <= threshold ? 0.0 : (model_coefficients[3]-threshold)*(model_coefficients[3]-threshold)));
301  const __m256 sqr_outer_radius = _mm256_set1_ps ((model_coefficients[3]+threshold)*(model_coefficients[3]+threshold));
302  __m256i res = _mm256_set1_epi32(0); // This corresponds to nr_p: 8 32bit integers that, summed together, hold the number of inliers
303  for (; (i + 8) <= indices_->size (); i += 8)
304  {
305  const __m256 sqr_dist = sqr_dist8 (i, a_vec, b_vec, c_vec);
306  const __m256 mask = _mm256_and_ps (_mm256_cmp_ps (sqr_inner_radius, sqr_dist, _CMP_LT_OQ), _mm256_cmp_ps (sqr_dist, sqr_outer_radius, _CMP_LT_OQ)); // The mask contains 1 bits if the corresponding points are inliers, else 0 bits
307  res = _mm256_add_epi32 (res, _mm256_and_si256 (_mm256_set1_epi32 (1), _mm256_castps_si256 (mask))); // The latter part creates a vector with ones (as 32bit integers) where the points are inliers
308  //const int res = _mm256_movemask_ps (mask);
309  //if (res & 1) nr_p++;
310  //if (res & 2) nr_p++;
311  //if (res & 4) nr_p++;
312  //if (res & 8) nr_p++;
313  //if (res & 16) nr_p++;
314  //if (res & 32) nr_p++;
315  //if (res & 64) nr_p++;
316  //if (res & 128) nr_p++;
317  }
318  nr_p += _mm256_extract_epi32 (res, 0);
319  nr_p += _mm256_extract_epi32 (res, 1);
320  nr_p += _mm256_extract_epi32 (res, 2);
321  nr_p += _mm256_extract_epi32 (res, 3);
322  nr_p += _mm256_extract_epi32 (res, 4);
323  nr_p += _mm256_extract_epi32 (res, 5);
324  nr_p += _mm256_extract_epi32 (res, 6);
325  nr_p += _mm256_extract_epi32 (res, 7);
326 
327  // Process the remaining points (at most 7)
328  nr_p += countWithinDistanceStandard (model_coefficients, threshold, i);
329  return (nr_p);
330 }
331 #endif
332 
333 //////////////////////////////////////////////////////////////////////////
334 template <typename PointT> void
336  const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const
337 {
338  optimized_coefficients = model_coefficients;
339 
340  // Needs a set of valid model coefficients
341  if (!isModelValid (model_coefficients))
342  {
343  PCL_ERROR ("[pcl::SampleConsensusModelSphere::optimizeModelCoefficients] Given model is invalid!\n");
344  return;
345  }
346 
347  // Need more than the minimum sample size to make a difference
348  if (inliers.size () <= sample_size_)
349  {
350  PCL_ERROR ("[pcl::SampleConsensusModelSphere::optimizeModelCoefficients] Not enough inliers to refine/optimize the model's coefficients (%lu)! Returning the same coefficients.\n", inliers.size ());
351  return;
352  }
353 
354  OptimizationFunctor functor (this, inliers);
355  Eigen::NumericalDiff<OptimizationFunctor> num_diff (functor);
356  Eigen::LevenbergMarquardt<Eigen::NumericalDiff<OptimizationFunctor>, float> lm (num_diff);
357  int info = lm.minimize (optimized_coefficients);
358 
359  // Compute the L2 norm of the residuals
360  PCL_DEBUG ("[pcl::SampleConsensusModelSphere::optimizeModelCoefficients] LM solver finished with exit code %i, having a residual norm of %g. \nInitial solution: %g %g %g %g \nFinal solution: %g %g %g %g\n",
361  info, lm.fvec.norm (), model_coefficients[0], model_coefficients[1], model_coefficients[2], model_coefficients[3], optimized_coefficients[0], optimized_coefficients[1], optimized_coefficients[2], optimized_coefficients[3]);
362 }
363 
364 //////////////////////////////////////////////////////////////////////////
365 template <typename PointT> void
367  const Indices &, const Eigen::VectorXf &model_coefficients, PointCloud &projected_points, bool) const
368 {
369  // Needs a valid model coefficients
370  if (!isModelValid (model_coefficients))
371  {
372  PCL_ERROR ("[pcl::SampleConsensusModelSphere::projectPoints] Given model is invalid!\n");
373  return;
374  }
375 
376  // Allocate enough space and copy the basics
377  projected_points.resize (input_->size ());
378  projected_points.header = input_->header;
379  projected_points.width = input_->width;
380  projected_points.height = input_->height;
381  projected_points.is_dense = input_->is_dense;
382 
383  PCL_WARN ("[pcl::SampleConsensusModelSphere::projectPoints] Not implemented yet.\n");
384  projected_points.points = input_->points;
385 }
386 
387 //////////////////////////////////////////////////////////////////////////
388 template <typename PointT> bool
390  const std::set<index_t> &indices, const Eigen::VectorXf &model_coefficients, const double threshold) const
391 {
392  // Needs a valid model coefficients
393  if (!isModelValid (model_coefficients))
394  {
395  PCL_ERROR ("[pcl::SampleConsensusModelSphere::doSamplesVerifyModel] Given model is invalid!\n");
396  return (false);
397  }
398 
399  const float sqr_inner_radius = (model_coefficients[3] <= threshold ? 0.0f : (model_coefficients[3] - threshold) * (model_coefficients[3] - threshold));
400  const float sqr_outer_radius = (model_coefficients[3] + threshold) * (model_coefficients[3] + threshold);
401  const Eigen::Vector3f center (model_coefficients[0], model_coefficients[1], model_coefficients[2]);
402  for (const auto &index : indices)
403  {
404  // To avoid sqrt computation: consider one larger sphere (radius + threshold) and one smaller sphere (radius - threshold).
405  // Valid if point is in larger sphere, but not in smaller sphere.
406  const float sqr_dist = ((*input_)[index].getVector3fMap () - center).squaredNorm ();
407  if ((sqr_dist > sqr_outer_radius) || (sqr_dist < sqr_inner_radius))
408  {
409  return (false);
410  }
411  }
412 
413  return (true);
414 }
415 
416 #define PCL_INSTANTIATE_SampleConsensusModelSphere(T) template class PCL_EXPORTS pcl::SampleConsensusModelSphere<T>;
417 
418 #endif // PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_SPHERE_H_
419 
pcl::SampleConsensusModelSphere::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_sphere.hpp:216
pcl::PointCloud::height
std::uint32_t height
The point cloud height (if organized as an image-structure).
Definition: point_cloud.h:393
pcl::SampleConsensusModelSphere::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_sphere.hpp:49
pcl::PointCloud::points
std::vector< PointT, Eigen::aligned_allocator< PointT > > points
The point data.
Definition: point_cloud.h:388
pcl::PointCloud< pcl::PointXYZRGB >
pcl::SampleConsensusModelSphere::getDistancesToModel
void getDistancesToModel(const Eigen::VectorXf &model_coefficients, std::vector< double > &distances) const override
Compute all distances from the cloud data to a given sphere model.
Definition: sac_model_sphere.hpp:157
pcl::SampleConsensusModelSphere::optimizeModelCoefficients
void optimizeModelCoefficients(const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const override
Recompute the sphere coefficients using the given inlier set and return them to the user.
Definition: sac_model_sphere.hpp:335
pcl::PointCloud::width
std::uint32_t width
The point cloud width (if organized as an image-structure).
Definition: point_cloud.h:391
pcl::SampleConsensusModelSphere::countWithinDistanceStandard
std::size_t countWithinDistanceStandard(const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i=0) const
This implementation uses no SIMD instructions.
Definition: sac_model_sphere.hpp:234
pcl::SampleConsensusModelSphere::computeModelCoefficients
bool computeModelCoefficients(const Indices &samples, Eigen::VectorXf &model_coefficients) const override
Check whether the given index samples can form a valid sphere model, compute the model coefficients f...
Definition: sac_model_sphere.hpp:61
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:396
pcl::PointCloud::resize
void resize(std::size_t count)
Resizes the container to contain count elements.
Definition: point_cloud.h:455
pcl::PointCloud::header
pcl::PCLHeader header
The point cloud header.
Definition: point_cloud.h:385
pcl::Indices
IndicesAllocator<> Indices
Type used for indices in PCL.
Definition: types.h:131
pcl::SampleConsensusModelSphere::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 sphere model coefficients.
Definition: sac_model_sphere.hpp:389
pcl::SampleConsensusModelSphere::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 sphere model.
Definition: sac_model_sphere.hpp:366
pcl::SampleConsensusModelSphere::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_sphere.hpp:180
pcl::SampleConsensusModelSphere
SampleConsensusModelSphere defines a model for 3D sphere segmentation.
Definition: sac_model_sphere.h:59