peel_outer_hull_layers.cpp 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126
  1. #include "peel_outer_hull_layers.h"
  2. #include "per_face_normals.h"
  3. #include "outer_hull.h"
  4. #include <vector>
  5. #include <iostream>
  6. //#define IGL_PEEL_OUTER_HULL_LAYERS_DEBUG
  7. #ifdef IGL_PEEL_OUTER_HULL_LAYERS_DEBUG
  8. #include "writePLY.h"
  9. #include "STR.h"
  10. #endif
  11. using namespace std;
  12. template <
  13. typename DerivedV,
  14. typename DerivedF,
  15. typename DerivedN,
  16. typename Derivedodd,
  17. typename Derivedflip>
  18. IGL_INLINE size_t igl::peel_outer_hull_layers(
  19. const Eigen::PlainObjectBase<DerivedV > & V,
  20. const Eigen::PlainObjectBase<DerivedF > & F,
  21. const Eigen::PlainObjectBase<DerivedN > & N,
  22. Eigen::PlainObjectBase<Derivedodd > & odd,
  23. Eigen::PlainObjectBase<Derivedflip > & flip)
  24. {
  25. using namespace Eigen;
  26. using namespace std;
  27. typedef typename DerivedF::Index Index;
  28. typedef Matrix<typename DerivedF::Scalar,Dynamic,DerivedF::ColsAtCompileTime> MatrixXF;
  29. typedef Matrix<typename DerivedN::Scalar,Dynamic,DerivedN::ColsAtCompileTime> MatrixXN;
  30. typedef Matrix<Index,Dynamic,1> MatrixXI;
  31. typedef Matrix<typename Derivedflip::Scalar,Dynamic,Derivedflip::ColsAtCompileTime> MatrixXflip;
  32. const Index m = F.rows();
  33. #ifdef IGL_PEEL_OUTER_HULL_LAYERS_DEBUG
  34. cout<<"peel outer hull layers..."<<endl;
  35. #endif
  36. #ifdef IGL_PEEL_OUTER_HULL_LAYERS_DEBUG
  37. cout<<"resize output ..."<<endl;
  38. #endif
  39. // keep track of iteration parity and whether flipped in hull
  40. MatrixXF Fr = F;
  41. MatrixXN Nr = N;
  42. odd.resize(m,1);
  43. flip.resize(m,1);
  44. // Keep track of index map
  45. MatrixXI IM = MatrixXI::LinSpaced(m,0,m-1);
  46. // This is O(n * layers)
  47. bool odd_iter = true;
  48. MatrixXI P(m,1);
  49. Index iter = 0;
  50. while(Fr.size() > 0)
  51. {
  52. assert(Fr.rows() == IM.rows());
  53. // Compute outer hull of current Fr
  54. MatrixXF Fo;
  55. MatrixXI Jo;
  56. MatrixXflip flipr;
  57. #ifdef IGL_PEEL_OUTER_HULL_LAYERS_DEBUG
  58. cout<<"calling outer hull..."<<endl;
  59. writePLY(STR("outer-hull-input-"<<iter<<".ply"),V,Fr);
  60. #endif
  61. outer_hull(V,Fr,Nr,Fo,Jo,flipr);
  62. #ifdef IGL_PEEL_OUTER_HULL_LAYERS_DEBUG
  63. writePLY(STR("outer-hull-output-"<<iter<<".ply"),V,Fo);
  64. cout<<"reindex, flip..."<<endl;
  65. #endif
  66. assert(Fo.rows() == Jo.rows());
  67. // all faces in Fo of Fr
  68. vector<bool> in_outer(Fr.rows(),false);
  69. for(Index g = 0;g<Jo.rows();g++)
  70. {
  71. odd(IM(Jo(g))) = odd_iter;
  72. P(IM(Jo(g))) = iter;
  73. in_outer[Jo(g)] = true;
  74. flip(IM(Jo(g))) = flipr(Jo(g));
  75. }
  76. // Fr = Fr - Fo
  77. // update IM
  78. MatrixXF prev_Fr = Fr;
  79. MatrixXN prev_Nr = Nr;
  80. MatrixXI prev_IM = IM;
  81. Fr.resize(prev_Fr.rows() - Fo.rows(),F.cols());
  82. Nr.resize(Fr.rows(),3);
  83. IM.resize(Fr.rows());
  84. {
  85. Index g = 0;
  86. for(Index f = 0;f<prev_Fr.rows();f++)
  87. {
  88. if(!in_outer[f])
  89. {
  90. Fr.row(g) = prev_Fr.row(f);
  91. Nr.row(g) = prev_Nr.row(f);
  92. IM(g) = prev_IM(f);
  93. g++;
  94. }
  95. }
  96. }
  97. odd_iter = !odd_iter;
  98. iter++;
  99. }
  100. return iter;
  101. }
  102. template <
  103. typename DerivedV,
  104. typename DerivedF,
  105. typename Derivedodd,
  106. typename Derivedflip>
  107. IGL_INLINE size_t igl::peel_outer_hull_layers(
  108. const Eigen::PlainObjectBase<DerivedV > & V,
  109. const Eigen::PlainObjectBase<DerivedF > & F,
  110. Eigen::PlainObjectBase<Derivedodd > & odd,
  111. Eigen::PlainObjectBase<Derivedflip > & flip)
  112. {
  113. Eigen::Matrix<typename DerivedV::Scalar,DerivedF::RowsAtCompileTime,3> N;
  114. per_face_normals(V,F,N);
  115. return peel_outer_hull_layers(V,F,N,odd,flip);
  116. }
  117. #ifdef IGL_STATIC_LIBRARY
  118. // Explicit template specialization
  119. template size_t igl::peel_outer_hull_layers<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<bool, -1, 1, 0, -1, 1>, Eigen::Matrix<bool, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 1, 0, -1, 1> >&);
  120. template size_t igl::peel_outer_hull_layers<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<bool, -1, 1, 0, -1, 1>, Eigen::Matrix<bool, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 1, 0, -1, 1> >&);
  121. #endif