rigid_alignment.cpp 3.2 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2019 Alec Jacobson <[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. #include "rigid_alignment.h"
  9. #include "polar_svd.h"
  10. #include "matlab_format.h"
  11. #include <Eigen/Sparse>
  12. #include <Eigen/Cholesky>
  13. #include <vector>
  14. #include <iostream>
  15. template <
  16. typename DerivedX,
  17. typename DerivedP,
  18. typename DerivedN,
  19. typename DerivedR,
  20. typename Derivedt
  21. >
  22. IGL_INLINE void igl::rigid_alignment(
  23. const Eigen::MatrixBase<DerivedX> & _X,
  24. const Eigen::MatrixBase<DerivedP> & P,
  25. const Eigen::MatrixBase<DerivedN> & N,
  26. Eigen::PlainObjectBase<DerivedR> & R,
  27. Eigen::PlainObjectBase<Derivedt> & t)
  28. {
  29. typedef typename DerivedX::Scalar Scalar;
  30. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> MatrixXS;
  31. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> VectorXS;
  32. typedef Eigen::Matrix<Scalar,3,3> Matrix3S;
  33. const int k = _X.rows();
  34. VectorXS Z = VectorXS::Zero(k,1);
  35. VectorXS I = VectorXS::Ones(k,1);
  36. DerivedX X = _X;
  37. R = DerivedR::Identity(3,3);
  38. t = Derivedt::Zero(1,3);
  39. // See gptoolbox, each iter could be O(1) instead of O(k)
  40. const int max_iters = 5;
  41. for(int iters = 0;iters<max_iters;iters++)
  42. {
  43. MatrixXS A(k*3,6);
  44. A <<
  45. Z, X.col(2),-X.col(1),I,Z,Z,
  46. -X.col(2), Z, X.col(0),Z,I,Z,
  47. X.col(1),-X.col(0), Z,Z,Z,I;
  48. VectorXS B(k*3,1);
  49. B<<
  50. P.col(0)-X.col(0),
  51. P.col(1)-X.col(1),
  52. P.col(2)-X.col(2);
  53. std::vector<Eigen::Triplet<Scalar> > NNIJV;
  54. for(int i = 0;i<k;i++)
  55. {
  56. for(int c = 0;c<3;c++)
  57. {
  58. NNIJV.emplace_back(i,i+k*c,N(i,c));
  59. }
  60. }
  61. Eigen::SparseMatrix<Scalar> NN(k,k*3);
  62. NN.setFromTriplets(NNIJV.begin(),NNIJV.end());
  63. A = (NN * A).eval();
  64. B = (NN * B).eval();
  65. VectorXS u = (A.transpose() * A).ldlt().solve(A.transpose() * B);
  66. Derivedt ti = u.tail(3).transpose();
  67. Matrix3S W;
  68. W<<
  69. 0, u(2),-u(1),
  70. -u(2), 0, u(0),
  71. u(1),-u(0), 0;
  72. // strayed from a perfect rotation. Correct it.
  73. const double x = u.head(3).stableNorm();
  74. DerivedR Ri;
  75. if(x == 0)
  76. {
  77. Ri = DerivedR::Identity(3,3);
  78. }else
  79. {
  80. Ri =
  81. DerivedR::Identity(3,3) +
  82. sin(x)/x*W +
  83. (1.0-cos(x))/(x*x)*W*W;
  84. }
  85. R = (R*Ri).eval();
  86. t = (t*Ri + ti).eval();
  87. X = ((_X*R).rowwise()+t).eval();
  88. }
  89. }
  90. #ifdef IGL_STATIC_LIBRARY
  91. // Explicit template instantiation
  92. template void igl::rigid_alignment<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, 3, 3, 0, 3, 3>, Eigen::Matrix<double, 1, 3, 1, 1, 3> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 3, 3, 0, 3, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, 3, 1, 1, 3> >&);
  93. #endif