Point Cloud Library (PCL)  1.11.1-dev
eigen.h
1 /*
2  * Software License Agreement (BSD License)
3  *
4  * Point Cloud Library (PCL) - www.pointclouds.org
5  * Copyright (c) 2010-2012, Willow Garage, Inc.
6  * Copyright (C) 2010 Gael Guennebaud <gael.guennebaud@inria.fr>
7  * Copyright (C) 2009 Hauke Heibel <hauke.heibel@gmail.com>
8  * Copyright (c) 2012-, Open Perception, Inc.
9  *
10  * All rights reserved.
11  *
12  * Redistribution and use in source and binary forms, with or without
13  * modification, are permitted provided that the following conditions
14  * are met:
15  *
16  * * Redistributions of source code must retain the above copyright
17  * notice, this list of conditions and the following disclaimer.
18  * * Redistributions in binary form must reproduce the above
19  * copyright notice, this list of conditions and the following
20  * disclaimer in the documentation and/or other materials provided
21  * with the distribution.
22  * * Neither the name of the copyright holder(s) nor the names of its
23  * contributors may be used to endorse or promote products derived
24  * from this software without specific prior written permission.
25  *
26  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
27  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
28  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
29  * FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
30  * COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
31  * INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
32  * BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
33  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
34  * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
35  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
36  * ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
37  * POSSIBILITY OF SUCH DAMAGE.
38  *
39  * $Id$
40  *
41  */
42 
43 #pragma once
44 
45 #ifndef NOMINMAX
46 #define NOMINMAX
47 #endif
48 
49 #if defined __GNUC__
50 # pragma GCC system_header
51 #elif defined __SUNPRO_CC
52 # pragma disable_warn
53 #endif
54 
55 #include <pcl/ModelCoefficients.h>
56 
57 #include <Eigen/StdVector>
58 #include <Eigen/Geometry>
59 #include <Eigen/LU>
60 
61 namespace pcl
62 {
63  /** \brief Compute the roots of a quadratic polynom x^2 + b*x + c = 0
64  * \param[in] b linear parameter
65  * \param[in] c constant parameter
66  * \param[out] roots solutions of x^2 + b*x + c = 0
67  */
68  template <typename Scalar, typename Roots> void
69  computeRoots2 (const Scalar &b, const Scalar &c, Roots &roots);
70 
71  /** \brief computes the roots of the characteristic polynomial of the input matrix m, which are the eigenvalues
72  * \param[in] m input matrix
73  * \param[out] roots roots of the characteristic polynomial of the input matrix m, which are the eigenvalues
74  */
75  template <typename Matrix, typename Roots> void
76  computeRoots (const Matrix &m, Roots &roots);
77 
78  /** \brief determine the smallest eigenvalue and its corresponding eigenvector
79  * \param[in] mat input matrix that needs to be symmetric and positive semi definite
80  * \param[out] eigenvalue the smallest eigenvalue of the input matrix
81  * \param[out] eigenvector the corresponding eigenvector to the smallest eigenvalue of the input matrix
82  * \ingroup common
83  */
84  template <typename Matrix, typename Vector> void
85  eigen22 (const Matrix &mat, typename Matrix::Scalar &eigenvalue, Vector &eigenvector);
86 
87  /** \brief determine the smallest eigenvalue and its corresponding eigenvector
88  * \param[in] mat input matrix that needs to be symmetric and positive semi definite
89  * \param[out] eigenvectors the corresponding eigenvector to the smallest eigenvalue of the input matrix
90  * \param[out] eigenvalues the smallest eigenvalue of the input matrix
91  * \ingroup common
92  */
93  template <typename Matrix, typename Vector> void
94  eigen22 (const Matrix &mat, Matrix &eigenvectors, Vector &eigenvalues);
95 
96  /** \brief determines the corresponding eigenvector to the given eigenvalue of the symmetric positive semi definite input matrix
97  * \param[in] mat symmetric positive semi definite input matrix
98  * \param[in] eigenvalue the eigenvalue which corresponding eigenvector is to be computed
99  * \param[out] eigenvector the corresponding eigenvector for the input eigenvalue
100  * \ingroup common
101  */
102  template <typename Matrix, typename Vector> void
103  computeCorrespondingEigenVector (const Matrix &mat, const typename Matrix::Scalar &eigenvalue, Vector &eigenvector);
104 
105  /** \brief determines the eigenvector and eigenvalue of the smallest eigenvalue of the symmetric positive semi definite input matrix
106  * \param[in] mat symmetric positive semi definite input matrix
107  * \param[out] eigenvalue smallest eigenvalue of the input matrix
108  * \param[out] eigenvector the corresponding eigenvector for the input eigenvalue
109  * \note if the smallest eigenvalue is not unique, this function may return any eigenvector that is consistent to the eigenvalue.
110  * \ingroup common
111  */
112  template <typename Matrix, typename Vector> void
113  eigen33 (const Matrix &mat, typename Matrix::Scalar &eigenvalue, Vector &eigenvector);
114 
115  /** \brief determines the eigenvalues of the symmetric positive semi definite input matrix
116  * \param[in] mat symmetric positive semi definite input matrix
117  * \param[out] evals resulting eigenvalues in ascending order
118  * \ingroup common
119  */
120  template <typename Matrix, typename Vector> void
121  eigen33 (const Matrix &mat, Vector &evals);
122 
123  /** \brief determines the eigenvalues and corresponding eigenvectors of the symmetric positive semi definite input matrix
124  * \param[in] mat symmetric positive semi definite input matrix
125  * \param[out] evecs corresponding eigenvectors in correct order according to eigenvalues
126  * \param[out] evals resulting eigenvalues in ascending order
127  * \ingroup common
128  */
129  template <typename Matrix, typename Vector> void
130  eigen33 (const Matrix &mat, Matrix &evecs, Vector &evals);
131 
132  /** \brief Calculate the inverse of a 2x2 matrix
133  * \param[in] matrix matrix to be inverted
134  * \param[out] inverse the resultant inverted matrix
135  * \note only the upper triangular part is taken into account => non symmetric matrices will give wrong results
136  * \return determinant of the original matrix => if 0 no inverse exists => result is invalid
137  * \ingroup common
138  */
139  template <typename Matrix> typename Matrix::Scalar
140  invert2x2 (const Matrix &matrix, Matrix &inverse);
141 
142  /** \brief Calculate the inverse of a 3x3 symmetric matrix.
143  * \param[in] matrix matrix to be inverted
144  * \param[out] inverse the resultant inverted matrix
145  * \note only the upper triangular part is taken into account => non symmetric matrices will give wrong results
146  * \return determinant of the original matrix => if 0 no inverse exists => result is invalid
147  * \ingroup common
148  */
149  template <typename Matrix> typename Matrix::Scalar
150  invert3x3SymMatrix (const Matrix &matrix, Matrix &inverse);
151 
152  /** \brief Calculate the inverse of a general 3x3 matrix.
153  * \param[in] matrix matrix to be inverted
154  * \param[out] inverse the resultant inverted matrix
155  * \return determinant of the original matrix => if 0 no inverse exists => result is invalid
156  * \ingroup common
157  */
158  template <typename Matrix> typename Matrix::Scalar
159  invert3x3Matrix (const Matrix &matrix, Matrix &inverse);
160 
161  /** \brief Calculate the determinant of a 3x3 matrix.
162  * \param[in] matrix matrix
163  * \return determinant of the matrix
164  * \ingroup common
165  */
166  template <typename Matrix> typename Matrix::Scalar
167  determinant3x3Matrix (const Matrix &matrix);
168 
169  /** \brief Get the unique 3D rotation that will rotate \a z_axis into (0,0,1) and \a y_direction into a vector
170  * with x=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
171  * \param[in] z_axis the z-axis
172  * \param[in] y_direction the y direction
173  * \param[out] transformation the resultant 3D rotation
174  * \ingroup common
175  */
176  inline void
177  getTransFromUnitVectorsZY (const Eigen::Vector3f& z_axis,
178  const Eigen::Vector3f& y_direction,
179  Eigen::Affine3f& transformation);
180 
181  /** \brief Get the unique 3D rotation that will rotate \a z_axis into (0,0,1) and \a y_direction into a vector
182  * with x=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
183  * \param[in] z_axis the z-axis
184  * \param[in] y_direction the y direction
185  * \return the resultant 3D rotation
186  * \ingroup common
187  */
188  inline Eigen::Affine3f
189  getTransFromUnitVectorsZY (const Eigen::Vector3f& z_axis,
190  const Eigen::Vector3f& y_direction);
191 
192  /** \brief Get the unique 3D rotation that will rotate \a x_axis into (1,0,0) and \a y_direction into a vector
193  * with z=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
194  * \param[in] x_axis the x-axis
195  * \param[in] y_direction the y direction
196  * \param[out] transformation the resultant 3D rotation
197  * \ingroup common
198  */
199  inline void
200  getTransFromUnitVectorsXY (const Eigen::Vector3f& x_axis,
201  const Eigen::Vector3f& y_direction,
202  Eigen::Affine3f& transformation);
203 
204  /** \brief Get the unique 3D rotation that will rotate \a x_axis into (1,0,0) and \a y_direction into a vector
205  * with z=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
206  * \param[in] x_axis the x-axis
207  * \param[in] y_direction the y direction
208  * \return the resulting 3D rotation
209  * \ingroup common
210  */
211  inline Eigen::Affine3f
212  getTransFromUnitVectorsXY (const Eigen::Vector3f& x_axis,
213  const Eigen::Vector3f& y_direction);
214 
215  /** \brief Get the unique 3D rotation that will rotate \a z_axis into (0,0,1) and \a y_direction into a vector
216  * with x=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
217  * \param[in] y_direction the y direction
218  * \param[in] z_axis the z-axis
219  * \param[out] transformation the resultant 3D rotation
220  * \ingroup common
221  */
222  inline void
223  getTransformationFromTwoUnitVectors (const Eigen::Vector3f& y_direction,
224  const Eigen::Vector3f& z_axis,
225  Eigen::Affine3f& transformation);
226 
227  /** \brief Get the unique 3D rotation that will rotate \a z_axis into (0,0,1) and \a y_direction into a vector
228  * with x=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
229  * \param[in] y_direction the y direction
230  * \param[in] z_axis the z-axis
231  * \return transformation the resultant 3D rotation
232  * \ingroup common
233  */
234  inline Eigen::Affine3f
235  getTransformationFromTwoUnitVectors (const Eigen::Vector3f& y_direction,
236  const Eigen::Vector3f& z_axis);
237 
238  /** \brief Get the transformation that will translate \a origin to (0,0,0) and rotate \a z_axis into (0,0,1)
239  * and \a y_direction into a vector with x=0 (or into (0,1,0) should \a y_direction be orthogonal to \a z_axis)
240  * \param[in] y_direction the y direction
241  * \param[in] z_axis the z-axis
242  * \param[in] origin the origin
243  * \param[in] transformation the resultant transformation matrix
244  * \ingroup common
245  */
246  inline void
247  getTransformationFromTwoUnitVectorsAndOrigin (const Eigen::Vector3f& y_direction,
248  const Eigen::Vector3f& z_axis,
249  const Eigen::Vector3f& origin,
250  Eigen::Affine3f& transformation);
251 
252  /** \brief Extract the Euler angles (XYZ-convention) from the given transformation
253  * \param[in] t the input transformation matrix
254  * \param[in] roll the resulting roll angle
255  * \param[in] pitch the resulting pitch angle
256  * \param[in] yaw the resulting yaw angle
257  * \ingroup common
258  */
259  template <typename Scalar> void
260  getEulerAngles (const Eigen::Transform<Scalar, 3, Eigen::Affine> &t, Scalar &roll, Scalar &pitch, Scalar &yaw);
261 
262  inline void
263  getEulerAngles (const Eigen::Affine3f &t, float &roll, float &pitch, float &yaw)
264  {
265  getEulerAngles<float> (t, roll, pitch, yaw);
266  }
267 
268  inline void
269  getEulerAngles (const Eigen::Affine3d &t, double &roll, double &pitch, double &yaw)
270  {
271  getEulerAngles<double> (t, roll, pitch, yaw);
272  }
273 
274  /** Extract x,y,z and the Euler angles (XYZ-convention) from the given transformation
275  * \param[in] t the input transformation matrix
276  * \param[out] x the resulting x translation
277  * \param[out] y the resulting y translation
278  * \param[out] z the resulting z translation
279  * \param[out] roll the resulting roll angle
280  * \param[out] pitch the resulting pitch angle
281  * \param[out] yaw the resulting yaw angle
282  * \ingroup common
283  */
284  template <typename Scalar> void
285  getTranslationAndEulerAngles (const Eigen::Transform<Scalar, 3, Eigen::Affine> &t,
286  Scalar &x, Scalar &y, Scalar &z,
287  Scalar &roll, Scalar &pitch, Scalar &yaw);
288 
289  inline void
290  getTranslationAndEulerAngles (const Eigen::Affine3f &t,
291  float &x, float &y, float &z,
292  float &roll, float &pitch, float &yaw)
293  {
294  getTranslationAndEulerAngles<float> (t, x, y, z, roll, pitch, yaw);
295  }
296 
297  inline void
298  getTranslationAndEulerAngles (const Eigen::Affine3d &t,
299  double &x, double &y, double &z,
300  double &roll, double &pitch, double &yaw)
301  {
302  getTranslationAndEulerAngles<double> (t, x, y, z, roll, pitch, yaw);
303  }
304 
305  /** \brief Create a transformation from the given translation and Euler angles (XYZ-convention)
306  * \param[in] x the input x translation
307  * \param[in] y the input y translation
308  * \param[in] z the input z translation
309  * \param[in] roll the input roll angle
310  * \param[in] pitch the input pitch angle
311  * \param[in] yaw the input yaw angle
312  * \param[out] t the resulting transformation matrix
313  * \ingroup common
314  */
315  template <typename Scalar> void
316  getTransformation (Scalar x, Scalar y, Scalar z, Scalar roll, Scalar pitch, Scalar yaw,
317  Eigen::Transform<Scalar, 3, Eigen::Affine> &t);
318 
319  inline void
320  getTransformation (float x, float y, float z, float roll, float pitch, float yaw,
321  Eigen::Affine3f &t)
322  {
323  return (getTransformation<float> (x, y, z, roll, pitch, yaw, t));
324  }
325 
326  inline void
327  getTransformation (double x, double y, double z, double roll, double pitch, double yaw,
328  Eigen::Affine3d &t)
329  {
330  return (getTransformation<double> (x, y, z, roll, pitch, yaw, t));
331  }
332 
333  /** \brief Create a transformation from the given translation and Euler angles (XYZ-convention)
334  * \param[in] x the input x translation
335  * \param[in] y the input y translation
336  * \param[in] z the input z translation
337  * \param[in] roll the input roll angle
338  * \param[in] pitch the input pitch angle
339  * \param[in] yaw the input yaw angle
340  * \return the resulting transformation matrix
341  * \ingroup common
342  */
343  inline Eigen::Affine3f
344  getTransformation (float x, float y, float z, float roll, float pitch, float yaw)
345  {
346  Eigen::Affine3f t;
347  getTransformation<float> (x, y, z, roll, pitch, yaw, t);
348  return (t);
349  }
350 
351  /** \brief Write a matrix to an output stream
352  * \param[in] matrix the matrix to output
353  * \param[out] file the output stream
354  * \ingroup common
355  */
356  template <typename Derived> void
357  saveBinary (const Eigen::MatrixBase<Derived>& matrix, std::ostream& file);
358 
359  /** \brief Read a matrix from an input stream
360  * \param[out] matrix the resulting matrix, read from the input stream
361  * \param[in,out] file the input stream
362  * \ingroup common
363  */
364  template <typename Derived> void
365  loadBinary (Eigen::MatrixBase<Derived> const& matrix, std::istream& file);
366 
367 // PCL_EIGEN_SIZE_MIN_PREFER_DYNAMIC gives the min between compile-time sizes. 0 has absolute priority, followed by 1,
368 // followed by Dynamic, followed by other finite values. The reason for giving Dynamic the priority over
369 // finite values is that min(3, Dynamic) should be Dynamic, since that could be anything between 0 and 3.
370 #define PCL_EIGEN_SIZE_MIN_PREFER_DYNAMIC(a,b) ((int (a) == 0 || int (b) == 0) ? 0 \
371  : (int (a) == 1 || int (b) == 1) ? 1 \
372  : (int (a) == Eigen::Dynamic || int (b) == Eigen::Dynamic) ? Eigen::Dynamic \
373  : (int (a) <= int (b)) ? int (a) : int (b))
374 
375  /** \brief Returns the transformation between two point sets.
376  * The algorithm is based on:
377  * "Least-squares estimation of transformation parameters between two point patterns",
378  * Shinji Umeyama, PAMI 1991, DOI: 10.1109/34.88573
379  *
380  * It estimates parameters \f$ c, \mathbf{R}, \f$ and \f$ \mathbf{t} \f$ such that
381  * \f{align*}
382  * \frac{1}{n} \sum_{i=1}^n \vert\vert y_i - (c\mathbf{R}x_i + \mathbf{t}) \vert\vert_2^2
383  * \f}
384  * is minimized.
385  *
386  * The algorithm is based on the analysis of the covariance matrix
387  * \f$ \Sigma_{\mathbf{x}\mathbf{y}} \in \mathbb{R}^{d \times d} \f$
388  * of the input point sets \f$ \mathbf{x} \f$ and \f$ \mathbf{y} \f$ where
389  * \f$d\f$ is corresponding to the dimension (which is typically small).
390  * The analysis is involving the SVD having a complexity of \f$O(d^3)\f$
391  * though the actual computational effort lies in the covariance
392  * matrix computation which has an asymptotic lower bound of \f$O(dm)\f$ when
393  * the input point sets have dimension \f$d \times m\f$.
394  *
395  * \param[in] src Source points \f$ \mathbf{x} = \left( x_1, \hdots, x_n \right) \f$
396  * \param[in] dst Destination points \f$ \mathbf{y} = \left( y_1, \hdots, y_n \right) \f$.
397  * \param[in] with_scaling Sets \f$ c=1 \f$ when <code>false</code> is passed. (default: false)
398  * \return The homogeneous transformation
399  * \f{align*}
400  * T = \begin{bmatrix} c\mathbf{R} & \mathbf{t} \\ \mathbf{0} & 1 \end{bmatrix}
401  * \f}
402  * minimizing the resudiual above. This transformation is always returned as an
403  * Eigen::Matrix.
404  */
405  template <typename Derived, typename OtherDerived>
406  typename Eigen::internal::umeyama_transform_matrix_type<Derived, OtherDerived>::type
407  umeyama (const Eigen::MatrixBase<Derived>& src, const Eigen::MatrixBase<OtherDerived>& dst, bool with_scaling = false);
408 
409 /** \brief Transform a point using an affine matrix
410  * \param[in] point_in the vector to be transformed
411  * \param[out] point_out the transformed vector
412  * \param[in] transformation the transformation matrix
413  *
414  * \note Can be used with \c point_in = \c point_out
415  */
416  template<typename Scalar> inline void
417  transformPoint (const Eigen::Matrix<Scalar, 3, 1> &point_in,
418  Eigen::Matrix<Scalar, 3, 1> &point_out,
419  const Eigen::Transform<Scalar, 3, Eigen::Affine> &transformation)
420  {
421  Eigen::Matrix<Scalar, 4, 1> point;
422  point << point_in, 1.0;
423  point_out = (transformation * point).template head<3> ();
424  }
425 
426  inline void
427  transformPoint (const Eigen::Vector3f &point_in,
428  Eigen::Vector3f &point_out,
429  const Eigen::Affine3f &transformation)
430  {
431  transformPoint<float> (point_in, point_out, transformation);
432  }
433 
434  inline void
435  transformPoint (const Eigen::Vector3d &point_in,
436  Eigen::Vector3d &point_out,
437  const Eigen::Affine3d &transformation)
438  {
439  transformPoint<double> (point_in, point_out, transformation);
440  }
441 
442 /** \brief Transform a vector using an affine matrix
443  * \param[in] vector_in the vector to be transformed
444  * \param[out] vector_out the transformed vector
445  * \param[in] transformation the transformation matrix
446  *
447  * \note Can be used with \c vector_in = \c vector_out
448  */
449  template <typename Scalar> inline void
450  transformVector (const Eigen::Matrix<Scalar, 3, 1> &vector_in,
451  Eigen::Matrix<Scalar, 3, 1> &vector_out,
452  const Eigen::Transform<Scalar, 3, Eigen::Affine> &transformation)
453  {
454  vector_out = transformation.linear () * vector_in;
455  }
456 
457  inline void
458  transformVector (const Eigen::Vector3f &vector_in,
459  Eigen::Vector3f &vector_out,
460  const Eigen::Affine3f &transformation)
461  {
462  transformVector<float> (vector_in, vector_out, transformation);
463  }
464 
465  inline void
466  transformVector (const Eigen::Vector3d &vector_in,
467  Eigen::Vector3d &vector_out,
468  const Eigen::Affine3d &transformation)
469  {
470  transformVector<double> (vector_in, vector_out, transformation);
471  }
472 
473 /** \brief Transform a line using an affine matrix
474  * \param[in] line_in the line to be transformed
475  * \param[out] line_out the transformed line
476  * \param[in] transformation the transformation matrix
477  *
478  * Lines must be filled in this form:\n
479  * line[0-2] = Origin coordinates of the vector\n
480  * line[3-5] = Direction vector
481  *
482  * \note Can be used with \c line_in = \c line_out
483  */
484  template <typename Scalar> bool
485  transformLine (const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &line_in,
486  Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &line_out,
487  const Eigen::Transform<Scalar, 3, Eigen::Affine> &transformation);
488 
489  inline bool
490  transformLine (const Eigen::VectorXf &line_in,
491  Eigen::VectorXf &line_out,
492  const Eigen::Affine3f &transformation)
493  {
494  return (transformLine<float> (line_in, line_out, transformation));
495  }
496 
497  inline bool
498  transformLine (const Eigen::VectorXd &line_in,
499  Eigen::VectorXd &line_out,
500  const Eigen::Affine3d &transformation)
501  {
502  return (transformLine<double> (line_in, line_out, transformation));
503  }
504 
505 /** \brief Transform plane vectors using an affine matrix
506  * \param[in] plane_in the plane coefficients to be transformed
507  * \param[out] plane_out the transformed plane coefficients to fill
508  * \param[in] transformation the transformation matrix
509  *
510  * The plane vectors are filled in the form ax+by+cz+d=0
511  * Can be used with non Hessian form planes coefficients
512  * Can be used with \c plane_in = \c plane_out
513  */
514  template <typename Scalar> void
515  transformPlane (const Eigen::Matrix<Scalar, 4, 1> &plane_in,
516  Eigen::Matrix<Scalar, 4, 1> &plane_out,
517  const Eigen::Transform<Scalar, 3, Eigen::Affine> &transformation);
518 
519  inline void
520  transformPlane (const Eigen::Matrix<double, 4, 1> &plane_in,
521  Eigen::Matrix<double, 4, 1> &plane_out,
522  const Eigen::Transform<double, 3, Eigen::Affine> &transformation)
523  {
524  transformPlane<double> (plane_in, plane_out, transformation);
525  }
526 
527  inline void
528  transformPlane (const Eigen::Matrix<float, 4, 1> &plane_in,
529  Eigen::Matrix<float, 4, 1> &plane_out,
530  const Eigen::Transform<float, 3, Eigen::Affine> &transformation)
531  {
532  transformPlane<float> (plane_in, plane_out, transformation);
533  }
534 
535 /** \brief Transform plane vectors using an affine matrix
536  * \param[in] plane_in the plane coefficients to be transformed
537  * \param[out] plane_out the transformed plane coefficients to fill
538  * \param[in] transformation the transformation matrix
539  *
540  * The plane vectors are filled in the form ax+by+cz+d=0
541  * Can be used with non Hessian form planes coefficients
542  * Can be used with \c plane_in = \c plane_out
543  * \warning ModelCoefficients stores floats only !
544  */
545  template<typename Scalar> void
547  pcl::ModelCoefficients::Ptr plane_out,
548  const Eigen::Transform<Scalar, 3, Eigen::Affine> &transformation);
549 
550  inline void
552  pcl::ModelCoefficients::Ptr plane_out,
553  const Eigen::Transform<double, 3, Eigen::Affine> &transformation)
554  {
555  transformPlane<double> (plane_in, plane_out, transformation);
556  }
557 
558  inline void
560  pcl::ModelCoefficients::Ptr plane_out,
561  const Eigen::Transform<float, 3, Eigen::Affine> &transformation)
562  {
563  transformPlane<float> (plane_in, plane_out, transformation);
564  }
565 
566 /** \brief Check coordinate system integrity
567  * \param[in] line_x the first axis
568  * \param[in] line_y the second axis
569  * \param[in] norm_limit the limit to ignore norm rounding errors
570  * \param[in] dot_limit the limit to ignore dot product rounding errors
571  * \return True if the coordinate system is consistent, false otherwise.
572  *
573  * Lines must be filled in this form:\n
574  * line[0-2] = Origin coordinates of the vector\n
575  * line[3-5] = Direction vector
576  *
577  * Can be used like this :\n
578  * line_x = X axis and line_y = Y axis\n
579  * line_x = Z axis and line_y = X axis\n
580  * line_x = Y axis and line_y = Z axis\n
581  * Because X^Y = Z, Z^X = Y and Y^Z = X.
582  * Do NOT invert line order !
583  *
584  * Determine whether a coordinate system is consistent or not by checking :\n
585  * Line origins: They must be the same for the 2 lines\n
586  * Norm: The 2 lines must be normalized\n
587  * Dot products: Must be 0 or perpendicular vectors
588  */
589  template<typename Scalar> bool
590  checkCoordinateSystem (const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &line_x,
591  const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &line_y,
592  const Scalar norm_limit = 1e-3,
593  const Scalar dot_limit = 1e-3);
594 
595  inline bool
596  checkCoordinateSystem (const Eigen::Matrix<double, Eigen::Dynamic, 1> &line_x,
597  const Eigen::Matrix<double, Eigen::Dynamic, 1> &line_y,
598  const double norm_limit = 1e-3,
599  const double dot_limit = 1e-3)
600  {
601  return (checkCoordinateSystem<double> (line_x, line_y, norm_limit, dot_limit));
602  }
603 
604  inline bool
605  checkCoordinateSystem (const Eigen::Matrix<float, Eigen::Dynamic, 1> &line_x,
606  const Eigen::Matrix<float, Eigen::Dynamic, 1> &line_y,
607  const float norm_limit = 1e-3,
608  const float dot_limit = 1e-3)
609  {
610  return (checkCoordinateSystem<float> (line_x, line_y, norm_limit, dot_limit));
611  }
612 
613 /** \brief Check coordinate system integrity
614  * \param[in] origin the origin of the coordinate system
615  * \param[in] x_direction the first axis
616  * \param[in] y_direction the second axis
617  * \param[in] norm_limit the limit to ignore norm rounding errors
618  * \param[in] dot_limit the limit to ignore dot product rounding errors
619  * \return True if the coordinate system is consistent, false otherwise.
620  *
621  * Read the other variant for more information
622  */
623  template <typename Scalar> inline bool
624  checkCoordinateSystem (const Eigen::Matrix<Scalar, 3, 1> &origin,
625  const Eigen::Matrix<Scalar, 3, 1> &x_direction,
626  const Eigen::Matrix<Scalar, 3, 1> &y_direction,
627  const Scalar norm_limit = 1e-3,
628  const Scalar dot_limit = 1e-3)
629  {
630  Eigen::Matrix<Scalar, Eigen::Dynamic, 1> line_x;
631  Eigen::Matrix<Scalar, Eigen::Dynamic, 1> line_y;
632  line_x << origin, x_direction;
633  line_y << origin, y_direction;
634  return (checkCoordinateSystem<Scalar> (line_x, line_y, norm_limit, dot_limit));
635  }
636 
637  inline bool
638  checkCoordinateSystem (const Eigen::Matrix<double, 3, 1> &origin,
639  const Eigen::Matrix<double, 3, 1> &x_direction,
640  const Eigen::Matrix<double, 3, 1> &y_direction,
641  const double norm_limit = 1e-3,
642  const double dot_limit = 1e-3)
643  {
644  Eigen::Matrix<double, Eigen::Dynamic, 1> line_x;
645  Eigen::Matrix<double, Eigen::Dynamic, 1> line_y;
646  line_x.resize (6);
647  line_y.resize (6);
648  line_x << origin, x_direction;
649  line_y << origin, y_direction;
650  return (checkCoordinateSystem<double> (line_x, line_y, norm_limit, dot_limit));
651  }
652 
653  inline bool
654  checkCoordinateSystem (const Eigen::Matrix<float, 3, 1> &origin,
655  const Eigen::Matrix<float, 3, 1> &x_direction,
656  const Eigen::Matrix<float, 3, 1> &y_direction,
657  const float norm_limit = 1e-3,
658  const float dot_limit = 1e-3)
659  {
660  Eigen::Matrix<float, Eigen::Dynamic, 1> line_x;
661  Eigen::Matrix<float, Eigen::Dynamic, 1> line_y;
662  line_x.resize (6);
663  line_y.resize (6);
664  line_x << origin, x_direction;
665  line_y << origin, y_direction;
666  return (checkCoordinateSystem<float> (line_x, line_y, norm_limit, dot_limit));
667  }
668 
669 /** \brief Compute the transformation between two coordinate systems
670  * \param[in] from_line_x X axis from the origin coordinate system
671  * \param[in] from_line_y Y axis from the origin coordinate system
672  * \param[in] to_line_x X axis from the destination coordinate system
673  * \param[in] to_line_y Y axis from the destination coordinate system
674  * \param[out] transformation the transformation matrix to fill
675  * \return true if transformation was filled, false otherwise.
676  *
677  * Line must be filled in this form:\n
678  * line[0-2] = Coordinate system origin coordinates \n
679  * line[3-5] = Direction vector (norm doesn't matter)
680  */
681  template <typename Scalar> bool
682  transformBetween2CoordinateSystems (const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> from_line_x,
683  const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> from_line_y,
684  const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> to_line_x,
685  const Eigen::Matrix<Scalar, Eigen::Dynamic, 1> to_line_y,
686  Eigen::Transform<Scalar, 3, Eigen::Affine> &transformation);
687 
688  inline bool
689  transformBetween2CoordinateSystems (const Eigen::Matrix<double, Eigen::Dynamic, 1> from_line_x,
690  const Eigen::Matrix<double, Eigen::Dynamic, 1> from_line_y,
691  const Eigen::Matrix<double, Eigen::Dynamic, 1> to_line_x,
692  const Eigen::Matrix<double, Eigen::Dynamic, 1> to_line_y,
693  Eigen::Transform<double, 3, Eigen::Affine> &transformation)
694  {
695  return (transformBetween2CoordinateSystems<double> (from_line_x, from_line_y, to_line_x, to_line_y, transformation));
696  }
697 
698  inline bool
699  transformBetween2CoordinateSystems (const Eigen::Matrix<float, Eigen::Dynamic, 1> from_line_x,
700  const Eigen::Matrix<float, Eigen::Dynamic, 1> from_line_y,
701  const Eigen::Matrix<float, Eigen::Dynamic, 1> to_line_x,
702  const Eigen::Matrix<float, Eigen::Dynamic, 1> to_line_y,
703  Eigen::Transform<float, 3, Eigen::Affine> &transformation)
704  {
705  return (transformBetween2CoordinateSystems<float> (from_line_x, from_line_y, to_line_x, to_line_y, transformation));
706  }
707 
708 }
709 
710 #include <pcl/common/impl/eigen.hpp>
711 
712 #if defined __SUNPRO_CC
713 # pragma enable_warn
714 #endif
pcl::saveBinary
void saveBinary(const Eigen::MatrixBase< Derived > &matrix, std::ostream &file)
Write a matrix to an output stream.
Definition: eigen.hpp:623
pcl
Definition: convolution.h:46
pcl::transformVector
void transformVector(const Eigen::Matrix< Scalar, 3, 1 > &vector_in, Eigen::Matrix< Scalar, 3, 1 > &vector_out, const Eigen::Transform< Scalar, 3, Eigen::Affine > &transformation)
Transform a vector using an affine matrix.
Definition: eigen.h:450
pcl::getTransFromUnitVectorsXY
void getTransFromUnitVectorsXY(const Eigen::Vector3f &x_axis, const Eigen::Vector3f &y_direction, Eigen::Affine3f &transformation)
Get the unique 3D rotation that will rotate x_axis into (1,0,0) and y_direction into a vector with z=...
Definition: eigen.hpp:528
pcl::getTransformationFromTwoUnitVectors
void getTransformationFromTwoUnitVectors(const Eigen::Vector3f &y_direction, const Eigen::Vector3f &z_axis, Eigen::Affine3f &transformation)
Get the unique 3D rotation that will rotate z_axis into (0,0,1) and y_direction into a vector with x=...
Definition: eigen.hpp:554
pcl::getTransFromUnitVectorsZY
void getTransFromUnitVectorsZY(const Eigen::Vector3f &z_axis, const Eigen::Vector3f &y_direction, Eigen::Affine3f &transformation)
Get the unique 3D rotation that will rotate z_axis into (0,0,1) and y_direction into a vector with x=...
Definition: eigen.hpp:502
pcl::getTransformation
void getTransformation(Scalar x, Scalar y, Scalar z, Scalar roll, Scalar pitch, Scalar yaw, Eigen::Transform< Scalar, 3, Eigen::Affine > &t)
Create a transformation from the given translation and Euler angles (XYZ-convention)
Definition: eigen.hpp:608
pcl::checkCoordinateSystem
bool checkCoordinateSystem(const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > &line_x, const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > &line_y, const Scalar norm_limit=1e-3, const Scalar dot_limit=1e-3)
Check coordinate system integrity.
Definition: eigen.hpp:788
pcl::umeyama
Eigen::internal::umeyama_transform_matrix_type< Derived, OtherDerived >::type umeyama(const Eigen::MatrixBase< Derived > &src, const Eigen::MatrixBase< OtherDerived > &dst, bool with_scaling=false)
Returns the transformation between two point sets.
Definition: eigen.hpp:660
pcl::ModelCoefficients::ConstPtr
shared_ptr< const ::pcl::ModelCoefficients > ConstPtr
Definition: ModelCoefficients.h:23
pcl::invert2x2
Matrix::Scalar invert2x2(const Matrix &matrix, Matrix &inverse)
Calculate the inverse of a 2x2 matrix.
Definition: eigen.hpp:405
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::computeRoots
void computeRoots(const Matrix &m, Roots &roots)
computes the roots of the characteristic polynomial of the input matrix m, which are the eigenvalues
Definition: eigen.hpp:68
pcl::eigen22
void eigen22(const Matrix &mat, typename Matrix::Scalar &eigenvalue, Vector &eigenvector)
determine the smallest eigenvalue and its corresponding eigenvector
Definition: eigen.hpp:133
pcl::invert3x3Matrix
Matrix::Scalar invert3x3Matrix(const Matrix &matrix, Matrix &inverse)
Calculate the inverse of a general 3x3 matrix.
Definition: eigen.hpp:459
pcl::getTranslationAndEulerAngles
void getTranslationAndEulerAngles(const Eigen::Transform< Scalar, 3, Eigen::Affine > &t, Scalar &x, Scalar &y, Scalar &z, Scalar &roll, Scalar &pitch, Scalar &yaw)
Extract x,y,z and the Euler angles (XYZ-convention) from the given transformation.
Definition: eigen.hpp:594
pcl::getTransformationFromTwoUnitVectorsAndOrigin
void getTransformationFromTwoUnitVectorsAndOrigin(const Eigen::Vector3f &y_direction, const Eigen::Vector3f &z_axis, const Eigen::Vector3f &origin, Eigen::Affine3f &transformation)
Get the transformation that will translate origin to (0,0,0) and rotate z_axis into (0,...
Definition: eigen.hpp:573
pcl::invert3x3SymMatrix
Matrix::Scalar invert3x3SymMatrix(const Matrix &matrix, Matrix &inverse)
Calculate the inverse of a 3x3 symmetric matrix.
Definition: eigen.hpp:424
pcl::getEulerAngles
void getEulerAngles(const Eigen::Transform< Scalar, 3, Eigen::Affine > &t, Scalar &roll, Scalar &pitch, Scalar &yaw)
Extract the Euler angles (XYZ-convention) from the given transformation.
Definition: eigen.hpp:585
pcl::transformBetween2CoordinateSystems
bool transformBetween2CoordinateSystems(const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > from_line_x, const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > from_line_y, const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > to_line_x, const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > to_line_y, Eigen::Transform< Scalar, 3, Eigen::Affine > &transformation)
Compute the transformation between two coordinate systems.
Definition: eigen.hpp:854
pcl::transformPlane
void transformPlane(const Eigen::Matrix< Scalar, 4, 1 > &plane_in, Eigen::Matrix< Scalar, 4, 1 > &plane_out, const Eigen::Transform< Scalar, 3, Eigen::Affine > &transformation)
Transform plane vectors using an affine matrix.
Definition: eigen.hpp:758
pcl::transformLine
bool transformLine(const Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > &line_in, Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > &line_out, const Eigen::Transform< Scalar, 3, Eigen::Affine > &transformation)
Transform a line using an affine matrix.
Definition: eigen.hpp:736
pcl::ModelCoefficients::Ptr
shared_ptr< ::pcl::ModelCoefficients > Ptr
Definition: ModelCoefficients.h:22
pcl::computeCorrespondingEigenVector
void computeCorrespondingEigenVector(const Matrix &mat, const typename Matrix::Scalar &eigenvalue, Vector &eigenvector)
determines the corresponding eigenvector to the given eigenvalue of the symmetric positive semi defin...
Definition: eigen.hpp:226
pcl::loadBinary
void loadBinary(Eigen::MatrixBase< Derived > const &matrix, std::istream &file)
Read a matrix from an input stream.
Definition: eigen.hpp:638
pcl::computeRoots2
void computeRoots2(const Scalar &b, const Scalar &c, Roots &roots)
Compute the roots of a quadratic polynom x^2 + b*x + c = 0.
Definition: eigen.hpp:53
pcl::determinant3x3Matrix
Matrix::Scalar determinant3x3Matrix(const Matrix &matrix)
Calculate the determinant of a 3x3 matrix.
Definition: eigen.hpp:492
pcl::transformPoint
void transformPoint(const Eigen::Matrix< Scalar, 3, 1 > &point_in, Eigen::Matrix< Scalar, 3, 1 > &point_out, const Eigen::Transform< Scalar, 3, Eigen::Affine > &transformation)
Transform a point using an affine matrix.
Definition: eigen.h:417