swept_volume_signed_distance.cpp 3.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113
  1. #include "swept_volume_signed_distance.h"
  2. #include "flood_fill.h"
  3. #include "signed_distance.h"
  4. #include "AABB.h"
  5. #include "pseudonormal_test.h"
  6. #include "per_face_normals.h"
  7. #include "per_vertex_normals.h"
  8. #include "per_edge_normals.h"
  9. #include <Eigen/Geometry>
  10. #include <cmath>
  11. IGL_INLINE void igl::swept_volume_signed_distance(
  12. const Eigen::MatrixXd & V,
  13. const Eigen::MatrixXi & F,
  14. const std::function<Eigen::Affine3d(const double t)> & transform,
  15. const size_t & steps,
  16. const Eigen::MatrixXd & GV,
  17. const Eigen::RowVector3i & res,
  18. const double h,
  19. const double isolevel,
  20. const Eigen::VectorXd & S0,
  21. Eigen::VectorXd & S)
  22. {
  23. using namespace std;
  24. using namespace igl;
  25. using namespace Eigen;
  26. S = S0;
  27. const VectorXd t = VectorXd::LinSpaced(steps,0,1);
  28. const bool finite_iso = isfinite(isolevel);
  29. const double extension = (finite_iso ? isolevel : 0) + sqrt(3.0)*h;
  30. Eigen::AlignedBox3d box(
  31. V.colwise().minCoeff().array()-extension,
  32. V.colwise().maxCoeff().array()+extension);
  33. // Precomputation
  34. Eigen::MatrixXd FN,VN,EN;
  35. Eigen::MatrixXi E;
  36. Eigen::VectorXi EMAP;
  37. per_face_normals(V,F,FN);
  38. per_vertex_normals(V,F,PER_VERTEX_NORMALS_WEIGHTING_TYPE_ANGLE,FN,VN);
  39. per_edge_normals(
  40. V,F,PER_EDGE_NORMALS_WEIGHTING_TYPE_UNIFORM,FN,EN,E,EMAP);
  41. AABB<MatrixXd,3> tree;
  42. tree.init(V,F);
  43. for(int ti = 0;ti<t.size();ti++)
  44. {
  45. const Affine3d At = transform(t(ti));
  46. for(int g = 0;g<GV.rows();g++)
  47. {
  48. // Don't bother finding out how deep inside points are.
  49. if(finite_iso && S(g)==S(g) && S(g)<isolevel-sqrt(3.0)*h)
  50. {
  51. continue;
  52. }
  53. const RowVector3d gv =
  54. (GV.row(g) - At.translation().transpose())*At.linear();
  55. // If outside of extended box, then consider it "far away enough"
  56. if(finite_iso && !box.contains(gv.transpose()))
  57. {
  58. continue;
  59. }
  60. RowVector3d c,n;
  61. int i;
  62. double sqrd,s;
  63. //signed_distance_pseudonormal(tree,V,F,FN,VN,EN,EMAP,gv,s,sqrd,i,c,n);
  64. const double min_sqrd =
  65. finite_iso ?
  66. pow(sqrt(3.)*h+isolevel,2) :
  67. numeric_limits<double>::infinity();
  68. sqrd = tree.squared_distance(V,F,gv,min_sqrd,i,c);
  69. if(sqrd<min_sqrd)
  70. {
  71. pseudonormal_test(V,F,FN,VN,EN,EMAP,gv,i,c,s,n);
  72. if(S(g) == S(g))
  73. {
  74. S(g) = min(S(g),s*sqrt(sqrd));
  75. }else
  76. {
  77. S(g) = s*sqrt(sqrd);
  78. }
  79. }
  80. }
  81. }
  82. if(finite_iso)
  83. {
  84. flood_fill(res,S);
  85. }else
  86. {
  87. #ifndef NDEBUG
  88. // Check for nans
  89. for_each(S.data(),S.data()+S.size(),[](const double s){assert(s==s);});
  90. #endif
  91. }
  92. }
  93. IGL_INLINE void igl::swept_volume_signed_distance(
  94. const Eigen::MatrixXd & V,
  95. const Eigen::MatrixXi & F,
  96. const std::function<Eigen::Affine3d(const double t)> & transform,
  97. const size_t & steps,
  98. const Eigen::MatrixXd & GV,
  99. const Eigen::RowVector3i & res,
  100. const double h,
  101. const double isolevel,
  102. Eigen::VectorXd & S)
  103. {
  104. using namespace std;
  105. using namespace igl;
  106. using namespace Eigen;
  107. S = VectorXd::Constant(GV.rows(),1,numeric_limits<double>::quiet_NaN());
  108. return
  109. swept_volume_signed_distance(V,F,transform,steps,GV,res,h,isolevel,S,S);
  110. }