iterative_closest_point.cpp 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117
  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 "iterative_closest_point.h"
  9. #include "AABB.h"
  10. #include "per_face_normals.h"
  11. #include "random_points_on_mesh.h"
  12. #include "placeholders.h"
  13. #include "rigid_alignment.h"
  14. #include <cassert>
  15. #include <iostream>
  16. template <
  17. typename DerivedVX,
  18. typename DerivedFX,
  19. typename DerivedVY,
  20. typename DerivedFY,
  21. typename DerivedR,
  22. typename Derivedt
  23. >
  24. IGL_INLINE void igl::iterative_closest_point(
  25. const Eigen::MatrixBase<DerivedVX> & VX,
  26. const Eigen::MatrixBase<DerivedFX> & FX,
  27. const Eigen::MatrixBase<DerivedVY> & VY,
  28. const Eigen::MatrixBase<DerivedFY> & FY,
  29. const int num_samples,
  30. const int max_iters,
  31. Eigen::PlainObjectBase<DerivedR> & R,
  32. Eigen::PlainObjectBase<Derivedt> & t)
  33. {
  34. assert(VX.cols() == 3 && "X should be a mesh in 3D");
  35. assert(VY.cols() == 3 && "Y should be a mesh in 3D");
  36. typedef typename DerivedVX::Scalar Scalar;
  37. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> MatrixXS;
  38. // Precompute BVH on Y
  39. AABB<DerivedVY,3> Ytree;
  40. Ytree.init(VY,FY);
  41. MatrixXS NY;
  42. per_face_normals(VY,FY,NY);
  43. return iterative_closest_point(
  44. VX,FX,VY,FY,Ytree,NY,num_samples,max_iters,R,t);
  45. }
  46. template <
  47. typename DerivedVX,
  48. typename DerivedFX,
  49. typename DerivedVY,
  50. typename DerivedFY,
  51. typename DerivedNY,
  52. typename DerivedR,
  53. typename Derivedt
  54. >
  55. IGL_INLINE void igl::iterative_closest_point(
  56. const Eigen::MatrixBase<DerivedVX> & VX,
  57. const Eigen::MatrixBase<DerivedFX> & FX,
  58. const Eigen::MatrixBase<DerivedVY> & VY,
  59. const Eigen::MatrixBase<DerivedFY> & FY,
  60. const igl::AABB<DerivedVY,3> & Ytree,
  61. const Eigen::MatrixBase<DerivedNY> & NY,
  62. const int num_samples,
  63. const int max_iters,
  64. Eigen::PlainObjectBase<DerivedR> & R,
  65. Eigen::PlainObjectBase<Derivedt> & t)
  66. {
  67. typedef typename DerivedVX::Scalar Scalar;
  68. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> MatrixXS;
  69. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> VectorXS;
  70. typedef Eigen::Matrix<Scalar,3,3> Matrix3S;
  71. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,3> MatrixX3S;
  72. typedef Eigen::Matrix<Scalar,1,3> RowVector3S;
  73. R.setIdentity(3,3);
  74. t.setConstant(1,3,0);
  75. for(int iter = 0;iter<max_iters;iter++)
  76. {
  77. // Sample random points on X
  78. MatrixXS X;
  79. {
  80. Eigen::VectorXi XI;
  81. MatrixX3S B;
  82. MatrixXS VXRT = (VX*R).rowwise()+t;
  83. random_points_on_mesh(num_samples,VXRT,FX,B,XI,X);
  84. }
  85. // Compute closest point
  86. Eigen::VectorXi I;
  87. MatrixXS P;
  88. {
  89. VectorXS sqrD;
  90. Ytree.squared_distance(VY,FY,X,sqrD,I,P);
  91. }
  92. // Use better normals?
  93. MatrixXS N = NY(I,igl::placeholders::all);
  94. //MatrixXS N = (X - P).rowwise().normalized();
  95. // fit rotation,translation
  96. Matrix3S Rup;
  97. RowVector3S tup;
  98. // Note: Should try out Szymon Rusinkiewicz's new symmetric icp
  99. rigid_alignment(X,P,N,Rup,tup);
  100. // update running rigid transformation
  101. R = (R*Rup).eval();
  102. t = (t*Rup + tup).eval();
  103. // Better stopping condition?
  104. }
  105. }
  106. #ifdef IGL_STATIC_LIBRARY
  107. // Explicit template instantiation
  108. template void igl::iterative_closest_point<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -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<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, int, int, Eigen::PlainObjectBase<Eigen::Matrix<double, 3, 3, 0, 3, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, 3, 1, 1, 3> >&);
  109. #endif