swept_volume_signed_distance.cpp 3.5 KB

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