Point Cloud Library (PCL)  1.11.1-dev
sac_model_plane.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_PLANE_H_
42 #define PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_PLANE_H_
43 
44 #include <pcl/sample_consensus/sac_model_plane.h>
45 #include <pcl/common/centroid.h>
46 #include <pcl/common/eigen.h>
47 #include <pcl/common/concatenate.h>
48 
49 //////////////////////////////////////////////////////////////////////////
50 template <typename PointT> bool
52 {
53  if (samples.size () != sample_size_)
54  {
55  PCL_ERROR ("[pcl::SampleConsensusModelPlane::isSampleGood] Wrong number of samples (is %lu, should be %lu)!\n", samples.size (), sample_size_);
56  return (false);
57  }
58  // Get the values at the two points
59  pcl::Array4fMapConst p0 = (*input_)[samples[0]].getArray4fMap ();
60  pcl::Array4fMapConst p1 = (*input_)[samples[1]].getArray4fMap ();
61  pcl::Array4fMapConst p2 = (*input_)[samples[2]].getArray4fMap ();
62 
63  Eigen::Array4f dy1dy2 = (p1-p0) / (p2-p0);
64 
65  return ( (dy1dy2[0] != dy1dy2[1]) || (dy1dy2[2] != dy1dy2[1]) );
66 }
67 
68 //////////////////////////////////////////////////////////////////////////
69 template <typename PointT> bool
71  const Indices &samples, Eigen::VectorXf &model_coefficients) const
72 {
73  // Need 3 samples
74  if (samples.size () != sample_size_)
75  {
76  PCL_ERROR ("[pcl::SampleConsensusModelPlane::computeModelCoefficients] Invalid set of samples given (%lu)!\n", samples.size ());
77  return (false);
78  }
79 
80  pcl::Array4fMapConst p0 = (*input_)[samples[0]].getArray4fMap ();
81  pcl::Array4fMapConst p1 = (*input_)[samples[1]].getArray4fMap ();
82  pcl::Array4fMapConst p2 = (*input_)[samples[2]].getArray4fMap ();
83 
84  // Compute the segment values (in 3d) between p1 and p0
85  Eigen::Array4f p1p0 = p1 - p0;
86  // Compute the segment values (in 3d) between p2 and p0
87  Eigen::Array4f p2p0 = p2 - p0;
88 
89  // Avoid some crashes by checking for collinearity here
90  Eigen::Array4f dy1dy2 = p1p0 / p2p0;
91  if ( (dy1dy2[0] == dy1dy2[1]) && (dy1dy2[2] == dy1dy2[1]) ) // Check for collinearity
92  {
93  return (false);
94  }
95 
96  // Compute the plane coefficients from the 3 given points in a straightforward manner
97  // calculate the plane normal n = (p2-p1) x (p3-p1) = cross (p2-p1, p3-p1)
98  model_coefficients.resize (model_size_);
99  model_coefficients[0] = p1p0[1] * p2p0[2] - p1p0[2] * p2p0[1];
100  model_coefficients[1] = p1p0[2] * p2p0[0] - p1p0[0] * p2p0[2];
101  model_coefficients[2] = p1p0[0] * p2p0[1] - p1p0[1] * p2p0[0];
102  model_coefficients[3] = 0.0f;
103 
104  // Normalize
105  model_coefficients.normalize ();
106 
107  // ... + d = 0
108  model_coefficients[3] = -1.0f * (model_coefficients.template head<4>().dot (p0.matrix ()));
109 
110  PCL_DEBUG ("[pcl::SampleConsensusModelPlane::computeModelCoefficients] Model is (%g,%g,%g,%g).\n",
111  model_coefficients[0], model_coefficients[1], model_coefficients[2], model_coefficients[3]);
112  return (true);
113 }
114 
115 #define AT(POS) ((*input_)[(*indices_)[(POS)]])
116 
117 #ifdef __AVX__
118 // This function computes the distances of 8 points to the plane
119 template <typename PointT> inline __m256 pcl::SampleConsensusModelPlane<PointT>::dist8 (const std::size_t i, const __m256 &a_vec, const __m256 &b_vec, const __m256 &c_vec, const __m256 &d_vec, const __m256 &abs_help) const
120 {
121  // The andnot-function realizes an abs-operation: the sign bit is removed
122  return _mm256_andnot_ps (abs_help,
123  _mm256_add_ps (_mm256_add_ps (_mm256_mul_ps (a_vec, _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)),
124  _mm256_mul_ps (b_vec, _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))),
125  _mm256_add_ps (_mm256_mul_ps (c_vec, _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)),
126  d_vec))); // TODO this could be replaced by three fmadd-instructions (if available), but the speed gain would probably be minimal
127 }
128 #endif // ifdef __AVX__
129 
130 #ifdef __SSE__
131 // This function computes the distances of 4 points to the plane
132 template <typename PointT> inline __m128 pcl::SampleConsensusModelPlane<PointT>::dist4 (const std::size_t i, const __m128 &a_vec, const __m128 &b_vec, const __m128 &c_vec, const __m128 &d_vec, const __m128 &abs_help) const
133 {
134  // The andnot-function realizes an abs-operation: the sign bit is removed
135  return _mm_andnot_ps (abs_help,
136  _mm_add_ps (_mm_add_ps (_mm_mul_ps (a_vec, _mm_set_ps (AT(i ).x, AT(i+1).x, AT(i+2).x, AT(i+3).x)),
137  _mm_mul_ps (b_vec, _mm_set_ps (AT(i ).y, AT(i+1).y, AT(i+2).y, AT(i+3).y))),
138  _mm_add_ps (_mm_mul_ps (c_vec, _mm_set_ps (AT(i ).z, AT(i+1).z, AT(i+2).z, AT(i+3).z)),
139  d_vec))); // TODO this could be replaced by three fmadd-instructions (if available), but the speed gain would probably be minimal
140 }
141 #endif // ifdef __SSE__
142 
143 #undef AT
144 
145 //////////////////////////////////////////////////////////////////////////
146 template <typename PointT> void
148  const Eigen::VectorXf &model_coefficients, std::vector<double> &distances) const
149 {
150  // Needs a valid set of model coefficients
151  if (!isModelValid (model_coefficients))
152  {
153  PCL_ERROR ("[pcl::SampleConsensusModelPlane::getDistancesToModel] Given model is invalid!\n");
154  return;
155  }
156 
157  distances.resize (indices_->size ());
158 
159  // Iterate through the 3d points and calculate the distances from them to the plane
160  for (std::size_t i = 0; i < indices_->size (); ++i)
161  {
162  // Calculate the distance from the point to the plane normal as the dot product
163  // D = (P-A).N/|N|
164  /*distances[i] = std::abs (model_coefficients[0] * (*input_)[(*indices_)[i]].x +
165  model_coefficients[1] * (*input_)[(*indices_)[i]].y +
166  model_coefficients[2] * (*input_)[(*indices_)[i]].z +
167  model_coefficients[3]);*/
168  Eigen::Vector4f pt ((*input_)[(*indices_)[i]].x,
169  (*input_)[(*indices_)[i]].y,
170  (*input_)[(*indices_)[i]].z,
171  1.0f);
172  distances[i] = std::abs (model_coefficients.dot (pt));
173  }
174 }
175 
176 //////////////////////////////////////////////////////////////////////////
177 template <typename PointT> void
179  const Eigen::VectorXf &model_coefficients, const double threshold, Indices &inliers)
180 {
181  // Needs a valid set of model coefficients
182  if (!isModelValid (model_coefficients))
183  {
184  PCL_ERROR ("[pcl::SampleConsensusModelPlane::selectWithinDistance] Given model is invalid!\n");
185  return;
186  }
187 
188  inliers.clear ();
189  error_sqr_dists_.clear ();
190  inliers.reserve (indices_->size ());
191  error_sqr_dists_.reserve (indices_->size ());
192 
193  // Iterate through the 3d points and calculate the distances from them to the plane
194  for (std::size_t i = 0; i < indices_->size (); ++i)
195  {
196  // Calculate the distance from the point to the plane normal as the dot product
197  // D = (P-A).N/|N|
198  Eigen::Vector4f pt ((*input_)[(*indices_)[i]].x,
199  (*input_)[(*indices_)[i]].y,
200  (*input_)[(*indices_)[i]].z,
201  1.0f);
202 
203  float distance = std::abs (model_coefficients.dot (pt));
204 
205  if (distance < threshold)
206  {
207  // Returns the indices of the points whose distances are smaller than the threshold
208  inliers.push_back ((*indices_)[i]);
209  error_sqr_dists_.push_back (static_cast<double> (distance));
210  }
211  }
212 }
213 
214 //////////////////////////////////////////////////////////////////////////
215 template <typename PointT> std::size_t
217  const Eigen::VectorXf &model_coefficients, const double threshold) const
218 {
219  // Needs a valid set of model coefficients
220  if (!isModelValid (model_coefficients))
221  {
222  PCL_ERROR ("[pcl::SampleConsensusModelPlane::countWithinDistance] Given model is invalid!\n");
223  return (0);
224  }
225 #if defined (__AVX__) && defined (__AVX2__)
226  return countWithinDistanceAVX (model_coefficients, threshold);
227 #elif defined (__SSE__) && defined (__SSE2__) && defined (__SSE4_1__)
228  return countWithinDistanceSSE (model_coefficients, threshold);
229 #else
230  return countWithinDistanceStandard (model_coefficients, threshold);
231 #endif
232 }
233 
234 //////////////////////////////////////////////////////////////////////////
235 template <typename PointT> std::size_t
237  const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i) const
238 {
239  std::size_t nr_p = 0;
240  // Iterate through the 3d points and calculate the distances from them to the plane
241  for (; i < indices_->size (); ++i)
242  {
243  // Calculate the distance from the point to the plane normal as the dot product
244  // D = (P-A).N/|N|
245  Eigen::Vector4f pt ((*input_)[(*indices_)[i]].x,
246  (*input_)[(*indices_)[i]].y,
247  (*input_)[(*indices_)[i]].z,
248  1.0f);
249  if (std::abs (model_coefficients.dot (pt)) < threshold)
250  {
251  nr_p++;
252  }
253  }
254  return (nr_p);
255 }
256 
257 //////////////////////////////////////////////////////////////////////////
258 #if defined (__SSE__) && defined (__SSE2__) && defined (__SSE4_1__)
259 template <typename PointT> std::size_t
261  const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i) const
262 {
263  std::size_t nr_p = 0;
264  const __m128 a_vec = _mm_set1_ps (model_coefficients[0]);
265  const __m128 b_vec = _mm_set1_ps (model_coefficients[1]);
266  const __m128 c_vec = _mm_set1_ps (model_coefficients[2]);
267  const __m128 d_vec = _mm_set1_ps (model_coefficients[3]);
268  const __m128 threshold_vec = _mm_set1_ps (threshold);
269  const __m128 abs_help = _mm_set1_ps (-0.0F); // -0.0F (negative zero) means that all bits are 0, only the sign bit is 1
270  __m128i res = _mm_set1_epi32(0); // This corresponds to nr_p: 4 32bit integers that, summed together, hold the number of inliers
271  for (; (i + 4) <= indices_->size (); i += 4)
272  {
273  const __m128 mask = _mm_cmplt_ps (dist4 (i, a_vec, b_vec, c_vec, d_vec, abs_help), threshold_vec); // The mask contains 1 bits if the corresponding points are inliers, else 0 bits
274  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
275  //const int res = _mm_movemask_ps (mask);
276  //if (res & 1) nr_p++;
277  //if (res & 2) nr_p++;
278  //if (res & 4) nr_p++;
279  //if (res & 8) nr_p++;
280  }
281  nr_p += _mm_extract_epi32 (res, 0);
282  nr_p += _mm_extract_epi32 (res, 1);
283  nr_p += _mm_extract_epi32 (res, 2);
284  nr_p += _mm_extract_epi32 (res, 3);
285 
286  // Process the remaining points (at most 3)
287  nr_p += countWithinDistanceStandard(model_coefficients, threshold, i);
288  return (nr_p);
289 }
290 #endif
291 
292 //////////////////////////////////////////////////////////////////////////
293 #if defined (__AVX__) && defined (__AVX2__)
294 template <typename PointT> std::size_t
296  const Eigen::VectorXf &model_coefficients, const double threshold, std::size_t i) const
297 {
298  std::size_t nr_p = 0;
299  const __m256 a_vec = _mm256_set1_ps (model_coefficients[0]);
300  const __m256 b_vec = _mm256_set1_ps (model_coefficients[1]);
301  const __m256 c_vec = _mm256_set1_ps (model_coefficients[2]);
302  const __m256 d_vec = _mm256_set1_ps (model_coefficients[3]);
303  const __m256 threshold_vec = _mm256_set1_ps (threshold);
304  const __m256 abs_help = _mm256_set1_ps (-0.0F); // -0.0F (negative zero) means that all bits are 0, only the sign bit is 1
305  __m256i res = _mm256_set1_epi32(0); // This corresponds to nr_p: 8 32bit integers that, summed together, hold the number of inliers
306  for (; (i + 8) <= indices_->size (); i += 8)
307  {
308  const __m256 mask = _mm256_cmp_ps (dist8 (i, a_vec, b_vec, c_vec, d_vec, abs_help), threshold_vec, _CMP_LT_OQ); // The mask contains 1 bits if the corresponding points are inliers, else 0 bits
309  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
310  //const int res = _mm256_movemask_ps (mask);
311  //if (res & 1) nr_p++;
312  //if (res & 2) nr_p++;
313  //if (res & 4) nr_p++;
314  //if (res & 8) nr_p++;
315  //if (res & 16) nr_p++;
316  //if (res & 32) nr_p++;
317  //if (res & 64) nr_p++;
318  //if (res & 128) nr_p++;
319  }
320  nr_p += _mm256_extract_epi32 (res, 0);
321  nr_p += _mm256_extract_epi32 (res, 1);
322  nr_p += _mm256_extract_epi32 (res, 2);
323  nr_p += _mm256_extract_epi32 (res, 3);
324  nr_p += _mm256_extract_epi32 (res, 4);
325  nr_p += _mm256_extract_epi32 (res, 5);
326  nr_p += _mm256_extract_epi32 (res, 6);
327  nr_p += _mm256_extract_epi32 (res, 7);
328 
329  // Process the remaining points (at most 7)
330  nr_p += countWithinDistanceStandard(model_coefficients, threshold, i);
331  return (nr_p);
332 }
333 #endif
334 
335 //////////////////////////////////////////////////////////////////////////
336 template <typename PointT> void
338  const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const
339 {
340  // Needs a valid set of model coefficients
341  if (!isModelValid (model_coefficients))
342  {
343  PCL_ERROR ("[pcl::SampleConsensusModelPlane::optimizeModelCoefficients] Given model is invalid!\n");
344  optimized_coefficients = model_coefficients;
345  return;
346  }
347 
348  // Need more than the minimum sample size to make a difference
349  if (inliers.size () <= sample_size_)
350  {
351  PCL_ERROR ("[pcl::SampleConsensusModelPlane::optimizeModelCoefficients] Not enough inliers found to optimize model coefficients (%lu)! Returning the same coefficients.\n", inliers.size ());
352  optimized_coefficients = model_coefficients;
353  return;
354  }
355 
356  Eigen::Vector4f plane_parameters;
357 
358  // Use Least-Squares to fit the plane through all the given sample points and find out its coefficients
359  EIGEN_ALIGN16 Eigen::Matrix3f covariance_matrix;
360  Eigen::Vector4f xyz_centroid;
361 
362  computeMeanAndCovarianceMatrix (*input_, inliers, covariance_matrix, xyz_centroid);
363 
364  // Compute the model coefficients
365  EIGEN_ALIGN16 Eigen::Vector3f::Scalar eigen_value;
366  EIGEN_ALIGN16 Eigen::Vector3f eigen_vector;
367  pcl::eigen33 (covariance_matrix, eigen_value, eigen_vector);
368 
369  // Hessian form (D = nc . p_plane (centroid here) + p)
370  optimized_coefficients.resize (model_size_);
371  optimized_coefficients[0] = eigen_vector [0];
372  optimized_coefficients[1] = eigen_vector [1];
373  optimized_coefficients[2] = eigen_vector [2];
374  optimized_coefficients[3] = 0.0f;
375  optimized_coefficients[3] = -1.0f * optimized_coefficients.dot (xyz_centroid);
376 
377  // Make sure it results in a valid model
378  if (!isModelValid (optimized_coefficients))
379  {
380  optimized_coefficients = model_coefficients;
381  }
382 }
383 
384 //////////////////////////////////////////////////////////////////////////
385 template <typename PointT> void
387  const Indices &inliers, const Eigen::VectorXf &model_coefficients, PointCloud &projected_points, bool copy_data_fields) const
388 {
389  // Needs a valid set of model coefficients
390  if (!isModelValid (model_coefficients))
391  {
392  PCL_ERROR ("[pcl::SampleConsensusModelPlane::projectPoints] Given model is invalid!\n");
393  return;
394  }
395 
396  projected_points.header = input_->header;
397  projected_points.is_dense = input_->is_dense;
398 
399  Eigen::Vector4f mc (model_coefficients[0], model_coefficients[1], model_coefficients[2], 0.0f);
400 
401  // normalize the vector perpendicular to the plane...
402  mc.normalize ();
403  // ... and store the resulting normal as a local copy of the model coefficients
404  Eigen::Vector4f tmp_mc = model_coefficients;
405  tmp_mc[0] = mc[0];
406  tmp_mc[1] = mc[1];
407  tmp_mc[2] = mc[2];
408 
409  // Copy all the data fields from the input cloud to the projected one?
410  if (copy_data_fields)
411  {
412  // Allocate enough space and copy the basics
413  projected_points.resize (input_->size ());
414  projected_points.width = input_->width;
415  projected_points.height = input_->height;
416 
417  using FieldList = typename pcl::traits::fieldList<PointT>::type;
418  // Iterate over each point
419  for (std::size_t i = 0; i < input_->size (); ++i)
420  // Iterate over each dimension
421  pcl::for_each_type <FieldList> (NdConcatenateFunctor <PointT, PointT> ((*input_)[i], projected_points[i]));
422 
423  // Iterate through the 3d points and calculate the distances from them to the plane
424  for (const auto &inlier : inliers)
425  {
426  // Calculate the distance from the point to the plane
427  Eigen::Vector4f p ((*input_)[inlier].x,
428  (*input_)[inlier].y,
429  (*input_)[inlier].z,
430  1);
431  // use normalized coefficients to calculate the scalar projection
432  float distance_to_plane = tmp_mc.dot (p);
433 
434  pcl::Vector4fMap pp = projected_points[inlier].getVector4fMap ();
435  pp.matrix () = p - mc * distance_to_plane; // mc[3] = 0, therefore the 3rd coordinate is safe
436  }
437  }
438  else
439  {
440  // Allocate enough space and copy the basics
441  projected_points.resize (inliers.size ());
442  projected_points.width = inliers.size ();
443  projected_points.height = 1;
444 
445  using FieldList = typename pcl::traits::fieldList<PointT>::type;
446  // Iterate over each point
447  for (std::size_t i = 0; i < inliers.size (); ++i)
448  {
449  // Iterate over each dimension
450  pcl::for_each_type <FieldList> (NdConcatenateFunctor <PointT, PointT> ((*input_)[inliers[i]], projected_points[i]));
451  }
452 
453  // Iterate through the 3d points and calculate the distances from them to the plane
454  for (std::size_t i = 0; i < inliers.size (); ++i)
455  {
456  // Calculate the distance from the point to the plane
457  Eigen::Vector4f p ((*input_)[inliers[i]].x,
458  (*input_)[inliers[i]].y,
459  (*input_)[inliers[i]].z,
460  1.0f);
461  // use normalized coefficients to calculate the scalar projection
462  float distance_to_plane = tmp_mc.dot (p);
463 
464  pcl::Vector4fMap pp = projected_points[i].getVector4fMap ();
465  pp.matrix () = p - mc * distance_to_plane; // mc[3] = 0, therefore the 3rd coordinate is safe
466  }
467  }
468 }
469 
470 //////////////////////////////////////////////////////////////////////////
471 template <typename PointT> bool
473  const std::set<index_t> &indices, const Eigen::VectorXf &model_coefficients, const double threshold) const
474 {
475  // Needs a valid set of model coefficients
476  if (!isModelValid (model_coefficients))
477  {
478  PCL_ERROR ("[pcl::SampleConsensusModelPlane::doSamplesVerifyModel] Given model is invalid!\n");
479  return (false);
480  }
481 
482  for (const auto &index : indices)
483  {
484  Eigen::Vector4f pt ((*input_)[index].x,
485  (*input_)[index].y,
486  (*input_)[index].z,
487  1.0f);
488  if (std::abs (model_coefficients.dot (pt)) > threshold)
489  {
490  return (false);
491  }
492  }
493 
494  return (true);
495 }
496 
497 #define PCL_INSTANTIATE_SampleConsensusModelPlane(T) template class PCL_EXPORTS pcl::SampleConsensusModelPlane<T>;
498 
499 #endif // PCL_SAMPLE_CONSENSUS_IMPL_SAC_MODEL_PLANE_H_
500 
pcl::computeMeanAndCovarianceMatrix
unsigned int computeMeanAndCovarianceMatrix(const pcl::PointCloud< PointT > &cloud, Eigen::Matrix< Scalar, 3, 3 > &covariance_matrix, Eigen::Matrix< Scalar, 4, 1 > &centroid)
Compute the normalized 3x3 covariance matrix and the centroid of a given set of points in a single lo...
Definition: centroid.hpp:485
pcl::PointCloud::height
std::uint32_t height
The point cloud height (if organized as an image-structure).
Definition: point_cloud.h:394
pcl::SampleConsensusModelPlane::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_plane.hpp:178
pcl::SampleConsensusModelPlane::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_plane.hpp:216
pcl::geometry::distance
float distance(const PointT &p1, const PointT &p2)
Definition: geometry.h:60
pcl::SampleConsensusModelPlane
SampleConsensusModelPlane defines a model for 3D plane segmentation.
Definition: sac_model_plane.h:142
pcl::NdConcatenateFunctor
Helper functor structure for concatenate.
Definition: concatenate.h:49
pcl::PointCloud
PointCloud represents the base class in PCL for storing collections of 3D points.
Definition: distances.h:55
pcl::eigen33
void eigen33(const Matrix &mat, typename Matrix::Scalar &eigenvalue, Vector &eigenvector)
determines the eigenvector and eigenvalue of the smallest eigenvalue of the symmetric positive semi d...
Definition: eigen.hpp:296
pcl::SampleConsensusModelPlane::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_plane.hpp:236
pcl::SampleConsensusModelPlane::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 plane model.
Definition: sac_model_plane.hpp:386
pcl::PointCloud::width
std::uint32_t width
The point cloud width (if organized as an image-structure).
Definition: point_cloud.h:392
pcl::SampleConsensusModelPlane::optimizeModelCoefficients
void optimizeModelCoefficients(const Indices &inliers, const Eigen::VectorXf &model_coefficients, Eigen::VectorXf &optimized_coefficients) const override
Recompute the plane coefficients using the given inlier set and return them to the user.
Definition: sac_model_plane.hpp:337
pcl::SampleConsensusModelPlane::computeModelCoefficients
bool computeModelCoefficients(const Indices &samples, Eigen::VectorXf &model_coefficients) const override
Check whether the given index samples can form a valid plane model, compute the model coefficients fr...
Definition: sac_model_plane.hpp:70
pcl::Array4fMapConst
const Eigen::Map< const Eigen::Array4f, Eigen::Aligned > Array4fMapConst
Definition: point_types.hpp:181
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::SampleConsensusModelPlane::getDistancesToModel
void getDistancesToModel(const Eigen::VectorXf &model_coefficients, std::vector< double > &distances) const override
Compute all distances from the cloud data to a given plane model.
Definition: sac_model_plane.hpp:147
pcl::SampleConsensusModelPlane::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 plane model coefficients.
Definition: sac_model_plane.hpp:472
pcl::Vector4fMap
Eigen::Map< Eigen::Vector4f, Eigen::Aligned > Vector4fMap
Definition: point_types.hpp:184
centroid.h