procrustes.h 3.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Stefan Brugger <[email protected]>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #ifndef IGL_PROCRUSTES_H
  9. #define IGL_PROCRUSTES_H
  10. #include "igl_inline.h"
  11. #include <Eigen/Dense>
  12. #include <Eigen/Geometry>
  13. namespace igl
  14. {
  15. /// Solve Procrustes problem in d dimensions. Given two point sets X,Y in R^d
  16. /// find best scale s, orthogonal R and translation t s.t. |s*X*R + t - Y|^2
  17. /// is minimized.
  18. ///
  19. /// @tparam DerivedV point type
  20. /// @tparam Scalar scalar type
  21. /// @tparam DerivedR type of R
  22. /// @tparam DerivedT type of t
  23. /// @param[in] X #V by DIM first list of points
  24. /// @param[in] Y #V by DIM second list of points
  25. /// @param[in] includeScaling if scaling should be allowed
  26. /// @param[in] includeReflections if R is allowed to be a reflection
  27. /// @param[out] scale scaling
  28. /// @param[out] R orthogonal matrix
  29. /// @param[out] t translation
  30. ///
  31. /// #### Example
  32. ///
  33. /// \code{cpp}
  34. /// Eigen::MatrixXd X, Y; (containing 3d points as rows)
  35. /// double scale;
  36. /// Eigen::MatrixXd R;
  37. /// Eigen::VectorXd t;
  38. /// igl::procrustes(X,Y,true,false,scale,R,t);
  39. /// R *= scale;
  40. /// Eigen::MatrixXd Xprime = (X * R).rowwise() + t.transpose();
  41. /// \endcode
  42. ///
  43. template <
  44. typename DerivedX,
  45. typename DerivedY,
  46. typename Scalar,
  47. typename DerivedR,
  48. typename DerivedT>
  49. IGL_INLINE void procrustes(
  50. const Eigen::MatrixBase<DerivedX>& X,
  51. const Eigen::MatrixBase<DerivedY>& Y,
  52. const bool includeScaling,
  53. const bool includeReflections,
  54. Scalar& scale,
  55. Eigen::PlainObjectBase<DerivedR>& R,
  56. Eigen::PlainObjectBase<DerivedT>& t);
  57. /// \overload
  58. /// \brief Same as above but returns Eigen transformation object.
  59. ///
  60. /// @param[out] T transformation that minimizes error
  61. ///
  62. /// #### Example
  63. /// \code{cpp}
  64. /// Eigen::MatrixXd X, Y; (containing 3d points as rows)
  65. /// AffineCompact3d T;
  66. /// igl::procrustes(X,Y,true,false,T);
  67. /// Eigen::MatrixXd Xprime = (X * T.linear()).rowwise() + T.translation().transpose();
  68. /// \endcode
  69. template <
  70. typename DerivedX,
  71. typename DerivedY,
  72. typename Scalar,
  73. int DIM,
  74. int TType>
  75. IGL_INLINE void procrustes(
  76. const Eigen::MatrixBase<DerivedX>& X,
  77. const Eigen::MatrixBase<DerivedY>& Y,
  78. const bool includeScaling,
  79. const bool includeReflections,
  80. Eigen::Transform<Scalar,DIM,TType>& T);
  81. /// \overload
  82. /// @param[out] S S=scale*R, instead of scale and R separately
  83. template <
  84. typename DerivedX,
  85. typename DerivedY,
  86. typename DerivedR,
  87. typename DerivedT>
  88. IGL_INLINE void procrustes(
  89. const Eigen::MatrixBase<DerivedX>& X,
  90. const Eigen::MatrixBase<DerivedY>& Y,
  91. const bool includeScaling,
  92. const bool includeReflections,
  93. Eigen::PlainObjectBase<DerivedR>& S,
  94. Eigen::PlainObjectBase<DerivedT>& t);
  95. /// \overload
  96. /// \brief Convenient wrapper for rigid case (no scaling, no reflections)
  97. template <
  98. typename DerivedX,
  99. typename DerivedY,
  100. typename DerivedR,
  101. typename DerivedT>
  102. IGL_INLINE void procrustes(
  103. const Eigen::MatrixBase<DerivedX>& X,
  104. const Eigen::MatrixBase<DerivedY>& Y,
  105. Eigen::PlainObjectBase<DerivedR>& R,
  106. Eigen::PlainObjectBase<DerivedT>& t);
  107. /// \overload
  108. /// \brief Convenient wrapper for 2D case.
  109. template <
  110. typename DerivedX,
  111. typename DerivedY,
  112. typename Scalar,
  113. typename DerivedT>
  114. IGL_INLINE void procrustes(
  115. const Eigen::MatrixBase<DerivedX>& X,
  116. const Eigen::MatrixBase<DerivedY>& Y,
  117. Eigen::Rotation2D<Scalar>& R,
  118. Eigen::PlainObjectBase<DerivedT>& t);
  119. }
  120. #ifndef IGL_STATIC_LIBRARY
  121. #include "procrustes.cpp"
  122. #endif
  123. #endif