Point Cloud Library (PCL)  1.9.1
gicp.h
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2010, 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_GICP_H_
42 #define PCL_GICP_H_
43 
44 #include <pcl/registration/icp.h>
45 #include <pcl/registration/bfgs.h>
46 
47 namespace pcl
48 {
49  /** \brief GeneralizedIterativeClosestPoint is an ICP variant that implements the
50  * generalized iterative closest point algorithm as described by Alex Segal et al. in
51  * http://www.robots.ox.ac.uk/~avsegal/resources/papers/Generalized_ICP.pdf
52  * The approach is based on using anistropic cost functions to optimize the alignment
53  * after closest point assignments have been made.
54  * The original code uses GSL and ANN while in ours we use an eigen mapped BFGS and
55  * FLANN.
56  * \author Nizar Sallem
57  * \ingroup registration
58  */
59  template <typename PointSource, typename PointTarget>
60  class GeneralizedIterativeClosestPoint : public IterativeClosestPoint<PointSource, PointTarget>
61  {
62  public:
81 
85 
89 
92 
93  typedef std::vector< Eigen::Matrix3d, Eigen::aligned_allocator<Eigen::Matrix3d> > MatricesVector;
94  typedef boost::shared_ptr< MatricesVector > MatricesVectorPtr;
95  typedef boost::shared_ptr< const MatricesVector > MatricesVectorConstPtr;
96 
99 
100  typedef boost::shared_ptr< GeneralizedIterativeClosestPoint<PointSource, PointTarget> > Ptr;
101  typedef boost::shared_ptr< const GeneralizedIterativeClosestPoint<PointSource, PointTarget> > ConstPtr;
102 
103 
104  typedef Eigen::Matrix<double, 6, 1> Vector6d;
105 
106  /** \brief Empty constructor. */
108  : k_correspondences_(20)
109  , gicp_epsilon_(0.001)
110  , rotation_epsilon_(2e-3)
111  , mahalanobis_(0)
113  {
115  reg_name_ = "GeneralizedIterativeClosestPoint";
116  max_iterations_ = 200;
121  this, _1, _2, _3, _4, _5);
122  }
123 
124  /** \brief Provide a pointer to the input dataset
125  * \param cloud the const boost shared pointer to a PointCloud message
126  */
127  inline void
129  {
130 
131  if (cloud->points.empty ())
132  {
133  PCL_ERROR ("[pcl::%s::setInputSource] Invalid or empty point cloud dataset given!\n", getClassName ().c_str ());
134  return;
135  }
136  PointCloudSource input = *cloud;
137  // Set all the point.data[3] values to 1 to aid the rigid transformation
138  for (size_t i = 0; i < input.size (); ++i)
139  input[i].data[3] = 1.0;
140 
142  input_covariances_.reset ();
143  }
144 
145  /** \brief Provide a pointer to the covariances of the input source (if computed externally!).
146  * If not set, GeneralizedIterativeClosestPoint will compute the covariances itself.
147  * Make sure to set the covariances AFTER setting the input source point cloud (setting the input source point cloud will reset the covariances).
148  * \param[in] target the input point cloud target
149  */
150  inline void
152  {
153  input_covariances_ = covariances;
154  }
155 
156  /** \brief Provide a pointer to the input target (e.g., the point cloud that we want to align the input source to)
157  * \param[in] target the input point cloud target
158  */
159  inline void
161  {
163  target_covariances_.reset ();
164  }
165 
166  /** \brief Provide a pointer to the covariances of the input target (if computed externally!).
167  * If not set, GeneralizedIterativeClosestPoint will compute the covariances itself.
168  * Make sure to set the covariances AFTER setting the input source point cloud (setting the input source point cloud will reset the covariances).
169  * \param[in] target the input point cloud target
170  */
171  inline void
173  {
174  target_covariances_ = covariances;
175  }
176 
177  /** \brief Estimate a rigid rotation transformation between a source and a target point cloud using an iterative
178  * non-linear Levenberg-Marquardt approach.
179  * \param[in] cloud_src the source point cloud dataset
180  * \param[in] indices_src the vector of indices describing the points of interest in \a cloud_src
181  * \param[in] cloud_tgt the target point cloud dataset
182  * \param[in] indices_tgt the vector of indices describing the correspondences of the interest points from \a indices_src
183  * \param[out] transformation_matrix the resultant transformation matrix
184  */
185  void
187  const std::vector<int> &indices_src,
188  const PointCloudTarget &cloud_tgt,
189  const std::vector<int> &indices_tgt,
190  Eigen::Matrix4f &transformation_matrix);
191 
192  /** \brief \return Mahalanobis distance matrix for the given point index */
193  inline const Eigen::Matrix3d& mahalanobis(size_t index) const
194  {
195  assert(index < mahalanobis_.size());
196  return mahalanobis_[index];
197  }
198 
199  /** \brief Computes rotation matrix derivative.
200  * rotation matrix is obtainded from rotation angles x[3], x[4] and x[5]
201  * \return d/d_rx, d/d_ry and d/d_rz respectively in g[3], g[4] and g[5]
202  * param x array representing 3D transformation
203  * param R rotation matrix
204  * param g gradient vector
205  */
206  void
207  computeRDerivative(const Vector6d &x, const Eigen::Matrix3d &R, Vector6d &g) const;
208 
209  /** \brief Set the rotation epsilon (maximum allowable difference between two
210  * consecutive rotations) in order for an optimization to be considered as having
211  * converged to the final solution.
212  * \param epsilon the rotation epsilon
213  */
214  inline void
215  setRotationEpsilon (double epsilon) { rotation_epsilon_ = epsilon; }
216 
217  /** \brief Get the rotation epsilon (maximum allowable difference between two
218  * consecutive rotations) as set by the user.
219  */
220  inline double
222 
223  /** \brief Set the number of neighbors used when selecting a point neighbourhood
224  * to compute covariances.
225  * A higher value will bring more accurate covariance matrix but will make
226  * covariances computation slower.
227  * \param k the number of neighbors to use when computing covariances
228  */
229  void
231 
232  /** \brief Get the number of neighbors used when computing covariances as set by
233  * the user
234  */
235  int
237 
238  /** set maximum number of iterations at the optimization step
239  * \param[in] max maximum number of iterations for the optimizer
240  */
241  void
243 
244  ///\return maximum number of iterations at the optimization step
245  int
247 
248  protected:
249 
250  /** \brief The number of neighbors used for covariances computation.
251  * default: 20
252  */
254 
255  /** \brief The epsilon constant for gicp paper; this is NOT the convergence
256  * tolerance
257  * default: 0.001
258  */
260 
261  /** The epsilon constant for rotation error. (In GICP the transformation epsilon
262  * is split in rotation part and translation part).
263  * default: 2e-3
264  */
266 
267  /** \brief base transformation */
268  Eigen::Matrix4f base_transformation_;
269 
270  /** \brief Temporary pointer to the source dataset. */
272 
273  /** \brief Temporary pointer to the target dataset. */
275 
276  /** \brief Temporary pointer to the source dataset indices. */
277  const std::vector<int> *tmp_idx_src_;
278 
279  /** \brief Temporary pointer to the target dataset indices. */
280  const std::vector<int> *tmp_idx_tgt_;
281 
282 
283  /** \brief Input cloud points covariances. */
285 
286  /** \brief Target cloud points covariances. */
288 
289  /** \brief Mahalanobis matrices holder. */
290  std::vector<Eigen::Matrix3d> mahalanobis_;
291 
292  /** \brief maximum number of optimizations */
294 
295  /** \brief compute points covariances matrices according to the K nearest
296  * neighbors. K is set via setCorrespondenceRandomness() method.
297  * \param cloud pointer to point cloud
298  * \param tree KD tree performer for nearest neighbors search
299  * \param[out] cloud_covariances covariances matrices for each point in the cloud
300  */
301  template<typename PointT>
303  const typename pcl::search::KdTree<PointT>::Ptr tree,
304  MatricesVector& cloud_covariances);
305 
306  /** \return trace of mat1^t . mat2
307  * \param mat1 matrix of dimension nxm
308  * \param mat2 matrix of dimension nxp
309  */
310  inline double
311  matricesInnerProd(const Eigen::MatrixXd& mat1, const Eigen::MatrixXd& mat2) const
312  {
313  double r = 0.;
314  size_t n = mat1.rows();
315  // tr(mat1^t.mat2)
316  for(size_t i = 0; i < n; i++)
317  for(size_t j = 0; j < n; j++)
318  r += mat1 (j, i) * mat2 (i,j);
319  return r;
320  }
321 
322  /** \brief Rigid transformation computation method with initial guess.
323  * \param output the transformed input point cloud dataset using the rigid transformation found
324  * \param guess the initial guess of the transformation to compute
325  */
326  void
327  computeTransformation (PointCloudSource &output, const Eigen::Matrix4f &guess);
328 
329  /** \brief Search for the closest nearest neighbor of a given point.
330  * \param query the point to search a nearest neighbour for
331  * \param index vector of size 1 to store the index of the nearest neighbour found
332  * \param distance vector of size 1 to store the distance to nearest neighbour found
333  */
334  inline bool
335  searchForNeighbors (const PointSource &query, std::vector<int>& index, std::vector<float>& distance)
336  {
337  int k = tree_->nearestKSearch (query, 1, index, distance);
338  if (k == 0)
339  return (false);
340  return (true);
341  }
342 
343  /// \brief compute transformation matrix from transformation matrix
344  void applyState(Eigen::Matrix4f &t, const Vector6d& x) const;
345 
346  /// \brief optimization functor structure
348  {
350  : BFGSDummyFunctor<double,6> (), gicp_(gicp) {}
351  double operator() (const Vector6d& x);
352  void df(const Vector6d &x, Vector6d &df);
353  void fdf(const Vector6d &x, double &f, Vector6d &df);
354 
356  };
357 
358  boost::function<void(const pcl::PointCloud<PointSource> &cloud_src,
359  const std::vector<int> &src_indices,
360  const pcl::PointCloud<PointTarget> &cloud_tgt,
361  const std::vector<int> &tgt_indices,
362  Eigen::Matrix4f &transformation_matrix)> rigid_transformation_estimation_;
363  };
364 }
365 
366 #include <pcl/registration/impl/gicp.hpp>
367 
368 #endif //#ifndef PCL_GICP_H_
pcl::GeneralizedIterativeClosestPoint::computeRDerivative
void computeRDerivative(const Vector6d &x, const Eigen::Matrix3d &R, Vector6d &g) const
Computes rotation matrix derivative.
Definition: gicp.hpp:127
pcl::GeneralizedIterativeClosestPoint::target_covariances_
MatricesVectorPtr target_covariances_
Target cloud points covariances.
Definition: gicp.h:287
pcl
This file defines compatibility wrappers for low level I/O functions.
Definition: convolution.h:45
pcl::IterativeClosestPoint::setInputSource
virtual void setInputSource(const PointCloudSourceConstPtr &cloud)
Provide a pointer to the input source (e.g., the point cloud that we want to align to the target)
Definition: icp.h:178
pcl::PointCloud< PointSource >::Ptr
boost::shared_ptr< PointCloud< PointSource > > Ptr
Definition: point_cloud.h:428
pcl::GeneralizedIterativeClosestPoint::base_transformation_
Eigen::Matrix4f base_transformation_
base transformation
Definition: gicp.h:268
pcl::GeneralizedIterativeClosestPoint::setRotationEpsilon
void setRotationEpsilon(double epsilon)
Set the rotation epsilon (maximum allowable difference between two consecutive rotations) in order fo...
Definition: gicp.h:215
pcl::Registration
Registration represents the base registration class for general purpose, ICP-like methods.
Definition: registration.h:62
pcl::geometry::distance
float distance(const PointT &p1, const PointT &p2)
Definition: geometry.h:61
pcl::GeneralizedIterativeClosestPoint::gicp_epsilon_
double gicp_epsilon_
The epsilon constant for gicp paper; this is NOT the convergence tolerance default: 0....
Definition: gicp.h:259
pcl::GeneralizedIterativeClosestPoint::estimateRigidTransformationBFGS
void estimateRigidTransformationBFGS(const PointCloudSource &cloud_src, const std::vector< int > &indices_src, const PointCloudTarget &cloud_tgt, const std::vector< int > &indices_tgt, Eigen::Matrix4f &transformation_matrix)
Estimate a rigid rotation transformation between a source and a target point cloud using an iterative...
Definition: gicp.hpp:182
pcl::GeneralizedIterativeClosestPoint::InputKdTreePtr
Registration< PointSource, PointTarget >::KdTreePtr InputKdTreePtr
Definition: gicp.h:98
pcl::GeneralizedIterativeClosestPoint::applyState
void applyState(Eigen::Matrix4f &t, const Vector6d &x) const
compute transformation matrix from transformation matrix
Definition: gicp.hpp:468
pcl::PointCloud::size
size_t size() const
Definition: point_cloud.h:448
pcl::GeneralizedIterativeClosestPoint::PointIndicesPtr
PointIndices::Ptr PointIndicesPtr
Definition: gicp.h:90
pcl::Registration< PointSource, PointTarget, float >::corr_dist_threshold_
double corr_dist_threshold_
The maximum distance threshold between two correspondent points in source <-> target.
Definition: registration.h:532
pcl::IterativeClosestPoint
IterativeClosestPoint provides a base implementation of the Iterative Closest Point algorithm.
Definition: icp.h:94
pcl::GeneralizedIterativeClosestPoint::MatricesVectorPtr
boost::shared_ptr< MatricesVector > MatricesVectorPtr
Definition: gicp.h:94
pcl::Registration< PointSource, PointTarget, float >::tree_
KdTreePtr tree_
A pointer to the spatial search object.
Definition: registration.h:485
pcl::GeneralizedIterativeClosestPoint::setSourceCovariances
void setSourceCovariances(const MatricesVectorPtr &covariances)
Provide a pointer to the covariances of the input source (if computed externally!).
Definition: gicp.h:151
pcl::Registration< PointSource, PointTarget, float >::getClassName
const std::string & getClassName() const
Abstract class get name method.
Definition: registration.h:422
pcl::Registration< PointSource, PointTarget, float >::max_iterations_
int max_iterations_
The maximum number of iterations the internal optimization should run for.
Definition: registration.h:496
pcl::Registration< PointSource, PointTarget, float >::min_number_correspondences_
int min_number_correspondences_
The minimum number of correspondences that the algorithm needs before attempting to estimate the tran...
Definition: registration.h:546
pcl::IterativeClosestPoint::setInputTarget
virtual void setInputTarget(const PointCloudTargetConstPtr &cloud)
Provide a pointer to the input target (e.g., the point cloud that we want to align to the target)
Definition: icp.h:213
pcl::PointCloud< PointSource >
pcl::GeneralizedIterativeClosestPoint::max_inner_iterations_
int max_inner_iterations_
maximum number of optimizations
Definition: gicp.h:293
pcl::GeneralizedIterativeClosestPoint::PointCloudSource
pcl::PointCloud< PointSource > PointCloudSource
Definition: gicp.h:82
pcl::PointIndices::Ptr
boost::shared_ptr< ::pcl::PointIndices > Ptr
Definition: PointIndices.h:22
pcl::GeneralizedIterativeClosestPoint::setInputTarget
void setInputTarget(const PointCloudTargetConstPtr &target)
Provide a pointer to the input target (e.g., the point cloud that we want to align the input source t...
Definition: gicp.h:160
pcl::GeneralizedIterativeClosestPoint::setCorrespondenceRandomness
void setCorrespondenceRandomness(int k)
Set the number of neighbors used when selecting a point neighbourhood to compute covariances.
Definition: gicp.h:230
pcl::GeneralizedIterativeClosestPoint::tmp_src_
const PointCloudSource * tmp_src_
Temporary pointer to the source dataset.
Definition: gicp.h:271
pcl::GeneralizedIterativeClosestPoint::searchForNeighbors
bool searchForNeighbors(const PointSource &query, std::vector< int > &index, std::vector< float > &distance)
Search for the closest nearest neighbor of a given point.
Definition: gicp.h:335
pcl::GeneralizedIterativeClosestPoint::PointIndicesConstPtr
PointIndices::ConstPtr PointIndicesConstPtr
Definition: gicp.h:91
pcl::search::KdTree::Ptr
boost::shared_ptr< KdTree< PointT, Tree > > Ptr
Definition: kdtree.h:79
pcl::Registration< PointSource, PointTarget, float >::transformation_epsilon_
double transformation_epsilon_
The maximum difference between two consecutive transformations in order to consider convergence (user...
Definition: registration.h:516
pcl::GeneralizedIterativeClosestPoint::tmp_tgt_
const PointCloudTarget * tmp_tgt_
Temporary pointer to the target dataset.
Definition: gicp.h:274
pcl::GeneralizedIterativeClosestPoint::tmp_idx_tgt_
const std::vector< int > * tmp_idx_tgt_
Temporary pointer to the target dataset indices.
Definition: gicp.h:280
pcl::GeneralizedIterativeClosestPoint::OptimizationFunctorWithIndices
optimization functor structure
Definition: gicp.h:347
pcl::GeneralizedIterativeClosestPoint::computeTransformation
void computeTransformation(PointCloudSource &output, const Eigen::Matrix4f &guess)
Rigid transformation computation method with initial guess.
Definition: gicp.hpp:345
pcl::GeneralizedIterativeClosestPoint::ConstPtr
boost::shared_ptr< const GeneralizedIterativeClosestPoint< PointSource, PointTarget > > ConstPtr
Definition: gicp.h:101
pcl::GeneralizedIterativeClosestPoint::Ptr
boost::shared_ptr< GeneralizedIterativeClosestPoint< PointSource, PointTarget > > Ptr
Definition: gicp.h:100
pcl::GeneralizedIterativeClosestPoint::Vector6d
Eigen::Matrix< double, 6, 1 > Vector6d
Definition: gicp.h:104
pcl::GeneralizedIterativeClosestPoint::OptimizationFunctorWithIndices::fdf
void fdf(const Vector6d &x, double &f, Vector6d &df)
Definition: gicp.hpp:308
pcl::GeneralizedIterativeClosestPoint::setInputSource
void setInputSource(const PointCloudSourceConstPtr &cloud)
Provide a pointer to the input dataset.
Definition: gicp.h:128
pcl::GeneralizedIterativeClosestPoint::OptimizationFunctorWithIndices::OptimizationFunctorWithIndices
OptimizationFunctorWithIndices(const GeneralizedIterativeClosestPoint *gicp)
Definition: gicp.h:349
pcl::GeneralizedIterativeClosestPoint::computeCovariances
void computeCovariances(typename pcl::PointCloud< PointT >::ConstPtr cloud, const typename pcl::search::KdTree< PointT >::Ptr tree, MatricesVector &cloud_covariances)
compute points covariances matrices according to the K nearest neighbors.
Definition: gicp.hpp:49
pcl::GeneralizedIterativeClosestPoint::setTargetCovariances
void setTargetCovariances(const MatricesVectorPtr &covariances)
Provide a pointer to the covariances of the input target (if computed externally!).
Definition: gicp.h:172
pcl::GeneralizedIterativeClosestPoint::MatricesVectorConstPtr
boost::shared_ptr< const MatricesVector > MatricesVectorConstPtr
Definition: gicp.h:95
pcl::GeneralizedIterativeClosestPoint::InputKdTree
Registration< PointSource, PointTarget >::KdTree InputKdTree
Definition: gicp.h:97
BFGSDummyFunctor
Definition: bfgs.h:74
pcl::GeneralizedIterativeClosestPoint::matricesInnerProd
double matricesInnerProd(const Eigen::MatrixXd &mat1, const Eigen::MatrixXd &mat2) const
Definition: gicp.h:311
pcl::PointIndices::ConstPtr
boost::shared_ptr< ::pcl::PointIndices const > ConstPtr
Definition: PointIndices.h:23
pcl::GeneralizedIterativeClosestPoint::PointCloudSourceConstPtr
PointCloudSource::ConstPtr PointCloudSourceConstPtr
Definition: gicp.h:84
pcl::GeneralizedIterativeClosestPoint::rotation_epsilon_
double rotation_epsilon_
The epsilon constant for rotation error.
Definition: gicp.h:265
pcl::GeneralizedIterativeClosestPoint::PointCloudTargetConstPtr
PointCloudTarget::ConstPtr PointCloudTargetConstPtr
Definition: gicp.h:88
pcl::PointCloud< PointSource >::ConstPtr
boost::shared_ptr< const PointCloud< PointSource > > ConstPtr
Definition: point_cloud.h:429
pcl::GeneralizedIterativeClosestPoint::getCorrespondenceRandomness
int getCorrespondenceRandomness()
Get the number of neighbors used when computing covariances as set by the user.
Definition: gicp.h:236
pcl::GeneralizedIterativeClosestPoint::getRotationEpsilon
double getRotationEpsilon()
Get the rotation epsilon (maximum allowable difference between two consecutive rotations) as set by t...
Definition: gicp.h:221
pcl::GeneralizedIterativeClosestPoint::rigid_transformation_estimation_
boost::function< void(const pcl::PointCloud< PointSource > &cloud_src, const std::vector< int > &src_indices, const pcl::PointCloud< PointTarget > &cloud_tgt, const std::vector< int > &tgt_indices, Eigen::Matrix4f &transformation_matrix)> rigid_transformation_estimation_
Definition: gicp.h:362
pcl::GeneralizedIterativeClosestPoint::mahalanobis
const Eigen::Matrix3d & mahalanobis(size_t index) const
Definition: gicp.h:193
pcl::GeneralizedIterativeClosestPoint::PointCloudTargetPtr
PointCloudTarget::Ptr PointCloudTargetPtr
Definition: gicp.h:87
pcl::GeneralizedIterativeClosestPoint::OptimizationFunctorWithIndices::df
void df(const Vector6d &x, Vector6d &df)
Definition: gicp.hpp:272
pcl::GeneralizedIterativeClosestPoint::OptimizationFunctorWithIndices::operator()
double operator()(const Vector6d &x)
Definition: gicp.hpp:247
pcl::GeneralizedIterativeClosestPoint::PointCloudSourcePtr
PointCloudSource::Ptr PointCloudSourcePtr
Definition: gicp.h:83
pcl::GeneralizedIterativeClosestPoint::tmp_idx_src_
const std::vector< int > * tmp_idx_src_
Temporary pointer to the source dataset indices.
Definition: gicp.h:277
pcl::GeneralizedIterativeClosestPoint::OptimizationFunctorWithIndices::gicp_
const GeneralizedIterativeClosestPoint * gicp_
Definition: gicp.h:355
pcl::GeneralizedIterativeClosestPoint::input_covariances_
MatricesVectorPtr input_covariances_
Input cloud points covariances.
Definition: gicp.h:284
pcl::GeneralizedIterativeClosestPoint::getMaximumOptimizerIterations
int getMaximumOptimizerIterations()
Definition: gicp.h:246
pcl::GeneralizedIterativeClosestPoint::setMaximumOptimizerIterations
void setMaximumOptimizerIterations(int max)
set maximum number of iterations at the optimization step
Definition: gicp.h:242
pcl::GeneralizedIterativeClosestPoint::mahalanobis_
std::vector< Eigen::Matrix3d > mahalanobis_
Mahalanobis matrices holder.
Definition: gicp.h:290
pcl::GeneralizedIterativeClosestPoint
GeneralizedIterativeClosestPoint is an ICP variant that implements the generalized iterative closest ...
Definition: gicp.h:60
pcl::GeneralizedIterativeClosestPoint::MatricesVector
std::vector< Eigen::Matrix3d, Eigen::aligned_allocator< Eigen::Matrix3d > > MatricesVector
Definition: gicp.h:93
pcl::GeneralizedIterativeClosestPoint::PointCloudTarget
pcl::PointCloud< PointTarget > PointCloudTarget
Definition: gicp.h:86
pcl::GeneralizedIterativeClosestPoint::k_correspondences_
int k_correspondences_
The number of neighbors used for covariances computation.
Definition: gicp.h:253
pcl::Registration< PointSource, PointTarget, float >::reg_name_
std::string reg_name_
The registration method name.
Definition: registration.h:482
pcl::GeneralizedIterativeClosestPoint::GeneralizedIterativeClosestPoint
GeneralizedIterativeClosestPoint()
Empty constructor.
Definition: gicp.h:107