Point Cloud Library (PCL)  1.14.1-dev
transformation_estimation_lm.h
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2010-2011, 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 #pragma once
42 
43 #include <pcl/registration/transformation_estimation.h>
44 #include <pcl/registration/warp_point_rigid.h>
45 #include <pcl/memory.h>
46 
47 namespace pcl {
48 namespace registration {
49 /** @b TransformationEstimationLM implements Levenberg Marquardt-based
50  * estimation of the transformation aligning the given correspondences.
51  *
52  * \note The class is templated on the source and target point types as well as on the
53  * output scalar of the transformation matrix (i.e., float or double). Default: float.
54  * \author Radu B. Rusu
55  * \ingroup registration
56  */
57 template <typename PointSource, typename PointTarget, typename MatScalar = float>
59 : public TransformationEstimation<PointSource, PointTarget, MatScalar> {
61  using PointCloudSourcePtr = typename PointCloudSource::Ptr;
62  using PointCloudSourceConstPtr = typename PointCloudSource::ConstPtr;
63 
65 
66  using PointIndicesPtr = PointIndices::Ptr;
67  using PointIndicesConstPtr = PointIndices::ConstPtr;
68 
69 public:
70  using Ptr =
71  shared_ptr<TransformationEstimationLM<PointSource, PointTarget, MatScalar>>;
72  using ConstPtr =
73  shared_ptr<const TransformationEstimationLM<PointSource, PointTarget, MatScalar>>;
74 
75  using VectorX = Eigen::Matrix<MatScalar, Eigen::Dynamic, 1>;
76  using Vector4 = Eigen::Matrix<MatScalar, 4, 1>;
77  using Matrix4 =
79 
80  /** \brief Constructor. */
82 
83  /** \brief Copy constructor.
84  * \param[in] src the TransformationEstimationLM object to copy into this
85  */
87  : tmp_src_(src.tmp_src_)
88  , tmp_tgt_(src.tmp_tgt_)
91  , warp_point_(src.warp_point_){};
92 
93  /** \brief Copy operator.
94  * \param[in] src the TransformationEstimationLM object to copy into this
95  */
98  {
99  tmp_src_ = src.tmp_src_;
100  tmp_tgt_ = src.tmp_tgt_;
103  warp_point_ = src.warp_point_;
104  return (*this);
105  }
106 
107  /** \brief Destructor. */
108  ~TransformationEstimationLM() override = default;
109 
110  /** \brief Estimate a rigid rotation transformation between a source and a target
111  * point cloud using LM. \param[in] cloud_src the source point cloud dataset
112  * \param[in] cloud_tgt the target point cloud dataset
113  * \param[out] transformation_matrix the resultant transformation matrix
114  */
115  inline void
117  const pcl::PointCloud<PointTarget>& cloud_tgt,
118  Matrix4& transformation_matrix) const override;
119 
120  /** \brief Estimate a rigid rotation transformation between a source and a target
121  * point cloud using LM. \param[in] cloud_src the source point cloud dataset
122  * \param[in] indices_src the vector of indices describing the points of interest in
123  * \a cloud_src
124  * \param[in] cloud_tgt the target point cloud dataset
125  * \param[out] transformation_matrix the resultant transformation matrix
126  */
127  inline void
129  const pcl::Indices& indices_src,
130  const pcl::PointCloud<PointTarget>& cloud_tgt,
131  Matrix4& transformation_matrix) const override;
132 
133  /** \brief Estimate a rigid rotation transformation between a source and a target
134  * point cloud using LM. \param[in] cloud_src the source point cloud dataset
135  * \param[in] indices_src the vector of indices describing the points of interest in
136  * \a cloud_src
137  * \param[in] cloud_tgt the target point cloud dataset
138  * \param[in] indices_tgt the vector of indices describing the correspondences of the
139  * interest points from \a indices_src
140  * \param[out] transformation_matrix the resultant transformation matrix
141  */
142  inline void
144  const pcl::Indices& indices_src,
145  const pcl::PointCloud<PointTarget>& cloud_tgt,
146  const pcl::Indices& indices_tgt,
147  Matrix4& transformation_matrix) const override;
148 
149  /** \brief Estimate a rigid rotation transformation between a source and a target
150  * point cloud using LM. \param[in] cloud_src the source point cloud dataset
151  * \param[in] cloud_tgt the target point cloud dataset
152  * \param[in] correspondences the vector of correspondences between source and target
153  * point cloud \param[out] transformation_matrix the resultant transformation matrix
154  */
155  inline void
157  const pcl::PointCloud<PointTarget>& cloud_tgt,
158  const pcl::Correspondences& correspondences,
159  Matrix4& transformation_matrix) const override;
160 
161  /** \brief Set the function we use to warp points. Defaults to rigid 6D warp.
162  * \param[in] warp_fcn a shared pointer to an object that warps points
163  */
164  void
167  {
168  warp_point_ = warp_fcn;
169  }
170 
171 protected:
172  /** \brief Compute the distance between a source point and its corresponding target
173  * point \param[in] p_src The source point \param[in] p_tgt The target point \return
174  * The distance between \a p_src and \a p_tgt
175  *
176  * \note Older versions of PCL used this method internally for calculating the
177  * optimization gradient. Since PCL 1.7, a switch has been made to the
178  * computeDistance method using Vector4 types instead. This method is only
179  * kept for API compatibility reasons.
180  */
181  virtual MatScalar
182  computeDistance(const PointSource& p_src, const PointTarget& p_tgt) const
183  {
184  Vector4 s(p_src.x, p_src.y, p_src.z, 0);
185  Vector4 t(p_tgt.x, p_tgt.y, p_tgt.z, 0);
186  return ((s - t).norm());
187  }
188 
189  /** \brief Compute the distance between a source point and its corresponding target
190  * point \param[in] p_src The source point \param[in] p_tgt The target point \return
191  * The distance between \a p_src and \a p_tgt
192  *
193  * \note A different distance function can be defined by creating a subclass of
194  * TransformationEstimationLM and overriding this method.
195  * (See \a TransformationEstimationPointToPlane)
196  */
197  virtual MatScalar
198  computeDistance(const Vector4& p_src, const PointTarget& p_tgt) const
199  {
200  Vector4 t(p_tgt.x, p_tgt.y, p_tgt.z, 0);
201  return ((p_src - t).norm());
202  }
203 
204  /** \brief Temporary pointer to the source dataset. */
205  mutable const PointCloudSource* tmp_src_{nullptr};
206 
207  /** \brief Temporary pointer to the target dataset. */
208  mutable const PointCloudTarget* tmp_tgt_{nullptr};
209 
210  /** \brief Temporary pointer to the source dataset indices. */
211  mutable const pcl::Indices* tmp_idx_src_{nullptr};
212 
213  /** \brief Temporary pointer to the target dataset indices. */
214  mutable const pcl::Indices* tmp_idx_tgt_{nullptr};
215 
216  /** \brief The parameterized function used to warp the source to the target. */
219 
220  /** Base functor all the models that need non linear optimization must
221  * define their own one and implement operator() (const Eigen::VectorXd& x,
222  * Eigen::VectorXd& fvec) or operator() (const Eigen::VectorXf& x, Eigen::VectorXf&
223  * fvec) depending on the chosen _Scalar
224  */
225  template <typename _Scalar, int NX = Eigen::Dynamic, int NY = Eigen::Dynamic>
226  struct Functor {
227  using Scalar = _Scalar;
229  using InputType = Eigen::Matrix<_Scalar, InputsAtCompileTime, 1>;
230  using ValueType = Eigen::Matrix<_Scalar, ValuesAtCompileTime, 1>;
231  using JacobianType =
232  Eigen::Matrix<_Scalar, ValuesAtCompileTime, InputsAtCompileTime>;
233 
234  /** \brief Empty Constructor. */
236 
237  /** \brief Constructor
238  * \param[in] m_data_points number of data points to evaluate.
239  */
240  Functor(int m_data_points) : m_data_points_(m_data_points) {}
241 
242  /** \brief Destructor. */
243  virtual ~Functor() = default;
244 
245  /** \brief Get the number of values. */
246  int
247  values() const
248  {
249  return (m_data_points_);
250  }
251 
252  protected:
254  };
255 
256  struct OptimizationFunctor : public Functor<MatScalar> {
258 
259  /** Functor constructor
260  * \param[in] m_data_points the number of data points to evaluate
261  * \param[in,out] estimator pointer to the estimator object
262  */
263  OptimizationFunctor(int m_data_points, const TransformationEstimationLM* estimator)
264  : Functor<MatScalar>(m_data_points), estimator_(estimator)
265  {}
266 
267  /** Copy constructor
268  * \param[in] src the optimization functor to copy into this
269  */
271  : Functor<MatScalar>(src.m_data_points_), estimator_()
272  {
273  *this = src;
274  }
275 
276  /** Copy operator
277  * \param[in] src the optimization functor to copy into this
278  */
279  inline OptimizationFunctor&
281  {
283  estimator_ = src.estimator_;
284  return (*this);
285  }
286 
287  /** \brief Destructor. */
288  ~OptimizationFunctor() override = default;
289 
290  /** Fill fvec from x. For the current state vector x fill the f values
291  * \param[in] x state vector
292  * \param[out] fvec f values vector
293  */
294  int
295  operator()(const VectorX& x, VectorX& fvec) const;
296 
298  };
299 
300  struct OptimizationFunctorWithIndices : public Functor<MatScalar> {
302 
303  /** Functor constructor
304  * \param[in] m_data_points the number of data points to evaluate
305  * \param[in,out] estimator pointer to the estimator object
306  */
308  const TransformationEstimationLM* estimator)
309  : Functor<MatScalar>(m_data_points), estimator_(estimator)
310  {}
311 
312  /** Copy constructor
313  * \param[in] src the optimization functor to copy into this
314  */
316  : Functor<MatScalar>(src.m_data_points_), estimator_()
317  {
318  *this = src;
319  }
320 
321  /** Copy operator
322  * \param[in] src the optimization functor to copy into this
323  */
326  {
328  estimator_ = src.estimator_;
329  return (*this);
330  }
331 
332  /** \brief Destructor. */
333  ~OptimizationFunctorWithIndices() override = default;
334 
335  /** Fill fvec from x. For the current state vector x fill the f values
336  * \param[in] x state vector
337  * \param[out] fvec f values vector
338  */
339  int
340  operator()(const VectorX& x, VectorX& fvec) const;
341 
343  };
344 
345 public:
347 };
348 } // namespace registration
349 } // namespace pcl
350 
351 #include <pcl/registration/impl/transformation_estimation_lm.hpp>
shared_ptr< PointCloud< PointSource > > Ptr
Definition: point_cloud.h:413
shared_ptr< const PointCloud< PointSource > > ConstPtr
Definition: point_cloud.h:414
TransformationEstimation represents the base class for methods for transformation estimation based on...
TransformationEstimationLM implements Levenberg Marquardt-based estimation of the transformation alig...
pcl::registration::WarpPointRigid< PointSource, PointTarget, MatScalar >::Ptr warp_point_
The parameterized function used to warp the source to the target.
TransformationEstimationLM(const TransformationEstimationLM &src)
Copy constructor.
shared_ptr< const TransformationEstimationLM< PointSource, PointTarget, MatScalar > > ConstPtr
virtual MatScalar computeDistance(const PointSource &p_src, const PointTarget &p_tgt) const
Compute the distance between a source point and its corresponding target point.
shared_ptr< TransformationEstimationLM< PointSource, PointTarget, MatScalar > > Ptr
void estimateRigidTransformation(const pcl::PointCloud< PointSource > &cloud_src, const pcl::PointCloud< PointTarget > &cloud_tgt, Matrix4 &transformation_matrix) const override
Estimate a rigid rotation transformation between a source and a target point cloud using LM.
virtual MatScalar computeDistance(const Vector4 &p_src, const PointTarget &p_tgt) const
Compute the distance between a source point and its corresponding target point.
Eigen::Matrix< MatScalar, Eigen::Dynamic, 1 > VectorX
const pcl::Indices * tmp_idx_tgt_
Temporary pointer to the target dataset indices.
const pcl::Indices * tmp_idx_src_
Temporary pointer to the source dataset indices.
const PointCloudSource * tmp_src_
Temporary pointer to the source dataset.
~TransformationEstimationLM() override=default
Destructor.
typename TransformationEstimation< PointSource, PointTarget, MatScalar >::Matrix4 Matrix4
TransformationEstimationLM & operator=(const TransformationEstimationLM &src)
Copy operator.
const PointCloudTarget * tmp_tgt_
Temporary pointer to the target dataset.
void setWarpFunction(const typename WarpPointRigid< PointSource, PointTarget, MatScalar >::Ptr &warp_fcn)
Set the function we use to warp points.
shared_ptr< WarpPointRigid< PointSourceT, PointTargetT, Scalar > > Ptr
#define PCL_MAKE_ALIGNED_OPERATOR_NEW
Macro to signal a class requires a custom allocator.
Definition: memory.h:63
Defines functions, macros and traits for allocating and using memory.
std::vector< pcl::Correspondence, Eigen::aligned_allocator< pcl::Correspondence > > Correspondences
IndicesAllocator<> Indices
Type used for indices in PCL.
Definition: types.h:133
shared_ptr< ::pcl::PointIndices > Ptr
Definition: PointIndices.h:13
shared_ptr< const ::pcl::PointIndices > ConstPtr
Definition: PointIndices.h:14
Base functor all the models that need non linear optimization must define their own one and implement...
Eigen::Matrix< _Scalar, InputsAtCompileTime, 1 > InputType
Eigen::Matrix< _Scalar, ValuesAtCompileTime, 1 > ValueType
Eigen::Matrix< _Scalar, ValuesAtCompileTime, InputsAtCompileTime > JacobianType
OptimizationFunctor & operator=(const OptimizationFunctor &src)
Copy operator.
OptimizationFunctor(int m_data_points, const TransformationEstimationLM *estimator)
Functor constructor.
int operator()(const VectorX &x, VectorX &fvec) const
Fill fvec from x.
const TransformationEstimationLM< PointSource, PointTarget, MatScalar > * estimator_
OptimizationFunctorWithIndices(const OptimizationFunctorWithIndices &src)
Copy constructor.
OptimizationFunctorWithIndices(int m_data_points, const TransformationEstimationLM *estimator)
Functor constructor.
OptimizationFunctorWithIndices & operator=(const OptimizationFunctorWithIndices &src)
Copy operator.
const TransformationEstimationLM< PointSource, PointTarget, MatScalar > * estimator_