signed_distance.cpp 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 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 "signed_distance.h"
  9. #include "per_vertex_normals.h"
  10. #include "per_edge_normals.h"
  11. #include "per_face_normals.h"
  12. #include "get_seconds.h"
  13. #include "point_mesh_squared_distance.h"
  14. IGL_INLINE void igl::signed_distance(
  15. const Eigen::MatrixXd & P,
  16. const Eigen::MatrixXd & V,
  17. const Eigen::MatrixXi & F,
  18. const SignedDistanceType sign_type,
  19. Eigen::VectorXd & S,
  20. Eigen::VectorXi & I,
  21. Eigen::MatrixXd & C,
  22. Eigen::MatrixXd & N)
  23. {
  24. using namespace Eigen;
  25. using namespace std;
  26. assert(V.cols() == 3 && "V should have 3d positions");
  27. assert(P.cols() == 3 && "P should have 3d positions");
  28. // Only unsigned distance is supported for non-triangles
  29. if(sign_type != SIGNED_DISTANCE_TYPE_UNSIGNED)
  30. {
  31. assert(F.cols() == 3 && "F should have triangles");
  32. }
  33. // Prepare distance computation
  34. AABB<MatrixXd,3> tree;
  35. tree.init(V,F);
  36. Eigen::MatrixXd FN,VN,EN;
  37. Eigen::MatrixXi E;
  38. Eigen::VectorXi EMAP;
  39. WindingNumberAABB<Eigen::Vector3d> hier;
  40. switch(sign_type)
  41. {
  42. default:
  43. assert(false && "Unknown SignedDistanceType");
  44. case SIGNED_DISTANCE_TYPE_UNSIGNED:
  45. // do nothing
  46. break;
  47. case SIGNED_DISTANCE_TYPE_DEFAULT:
  48. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  49. hier.set_mesh(V,F);
  50. hier.grow();
  51. break;
  52. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  53. // "Signed Distance Computation Using the Angle Weighted Pseudonormal"
  54. // [Bærentzen & Aanæs 2005]
  55. per_face_normals(V,F,FN);
  56. per_vertex_normals(V,F,PER_VERTEX_NORMALS_WEIGHTING_TYPE_ANGLE,FN,VN);
  57. per_edge_normals(
  58. V,F,PER_EDGE_NORMALS_WEIGHTING_TYPE_UNIFORM,FN,EN,E,EMAP);
  59. N.resize(P.rows(),3);
  60. break;
  61. }
  62. S.resize(P.rows(),1);
  63. I.resize(P.rows(),1);
  64. C.resize(P.rows(),3);
  65. for(int p = 0;p<P.rows();p++)
  66. {
  67. const RowVector3d q = P.row(p);
  68. double s,sqrd;
  69. RowVector3d c;
  70. int i=-1;
  71. switch(sign_type)
  72. {
  73. default:
  74. assert(false && "Unknown SignedDistanceType");
  75. case SIGNED_DISTANCE_TYPE_UNSIGNED:
  76. s = 1.;
  77. sqrd = tree.squared_distance(V,F,q,i,c);
  78. break;
  79. case SIGNED_DISTANCE_TYPE_DEFAULT:
  80. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  81. signed_distance_winding_number(tree,V,F,hier,q,s,sqrd,i,c);
  82. break;
  83. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  84. {
  85. Eigen::RowVector3d n;
  86. signed_distance_pseudonormal(tree,V,F,FN,VN,EN,EMAP,q,s,sqrd,i,c,n);
  87. N.row(p) = n;
  88. break;
  89. }
  90. }
  91. I(p) = i;
  92. S(p) = s*sqrt(sqrd);
  93. C.row(p) = c;
  94. }
  95. }
  96. IGL_INLINE double igl::signed_distance_pseudonormal(
  97. const AABB<Eigen::MatrixXd,3> & tree,
  98. const Eigen::MatrixXd & V,
  99. const Eigen::MatrixXi & F,
  100. const Eigen::MatrixXd & FN,
  101. const Eigen::MatrixXd & VN,
  102. const Eigen::MatrixXd & EN,
  103. const Eigen::VectorXi & EMAP,
  104. const Eigen::RowVector3d & q)
  105. {
  106. double s,sqrd;
  107. Eigen::RowVector3d n,c;
  108. int i = -1;
  109. signed_distance_pseudonormal(tree,V,F,FN,VN,EN,EMAP,q,s,sqrd,i,c,n);
  110. return s*sqrt(sqrd);
  111. }
  112. IGL_INLINE void igl::signed_distance_pseudonormal(
  113. const AABB<Eigen::MatrixXd,3> & tree,
  114. const Eigen::MatrixXd & V,
  115. const Eigen::MatrixXi & F,
  116. const Eigen::MatrixXd & FN,
  117. const Eigen::MatrixXd & VN,
  118. const Eigen::MatrixXd & EN,
  119. const Eigen::VectorXi & EMAP,
  120. const Eigen::RowVector3d & q,
  121. double & s,
  122. double & sqrd,
  123. int & f,
  124. Eigen::RowVector3d & c,
  125. Eigen::RowVector3d & n)
  126. {
  127. using namespace Eigen;
  128. using namespace std;
  129. sqrd = tree.squared_distance(V,F,q,f,c);
  130. const auto & qc = q-c;
  131. RowVector3d b;
  132. AABB<Eigen::MatrixXd,3>::barycentric_coordinates(
  133. c,V.row(F(f,0)),V.row(F(f,1)),V.row(F(f,2)),b);
  134. // Determine which normal to use
  135. const double epsilon = 1e-12;
  136. const int type = (b.array()<=epsilon).cast<int>().sum();
  137. switch(type)
  138. {
  139. case 2:
  140. // Find vertex
  141. for(int x = 0;x<3;x++)
  142. {
  143. if(b(x)>epsilon)
  144. {
  145. n = VN.row(F(f,x));
  146. break;
  147. }
  148. }
  149. break;
  150. case 1:
  151. // Find edge
  152. for(int x = 0;x<3;x++)
  153. {
  154. if(b(x)<=epsilon)
  155. {
  156. n = EN.row(EMAP(F.rows()*x+f));
  157. break;
  158. }
  159. }
  160. break;
  161. default:
  162. assert(false && "all barycentric coords zero.");
  163. case 0:
  164. n = FN.row(f);
  165. break;
  166. }
  167. s = (qc.dot(n) >= 0 ? 1. : -1.);
  168. }
  169. IGL_INLINE double igl::signed_distance_winding_number(
  170. const AABB<Eigen::MatrixXd,3> & tree,
  171. const Eigen::MatrixXd & V,
  172. const Eigen::MatrixXi & F,
  173. const igl::WindingNumberAABB<Eigen::Vector3d> & hier,
  174. const Eigen::RowVector3d & q)
  175. {
  176. double s,sqrd;
  177. Eigen::RowVector3d c;
  178. int i=-1;
  179. signed_distance_winding_number(tree,V,F,hier,q,s,sqrd,i,c);
  180. return s*sqrt(sqrd);
  181. }
  182. IGL_INLINE void igl::signed_distance_winding_number(
  183. const AABB<Eigen::MatrixXd,3> & tree,
  184. const Eigen::MatrixXd & V,
  185. const Eigen::MatrixXi & F,
  186. const igl::WindingNumberAABB<Eigen::Vector3d> & hier,
  187. const Eigen::RowVector3d & q,
  188. double & s,
  189. double & sqrd,
  190. int & i,
  191. Eigen::RowVector3d & c)
  192. {
  193. using namespace Eigen;
  194. using namespace std;
  195. using namespace igl;
  196. sqrd = tree.squared_distance(V,F,q,i,c);
  197. const double w = hier.winding_number(q.transpose());
  198. s = 1.-2.*w;
  199. }
  200. #ifdef IGL_STATIC_LIBRARY
  201. // This template is necessary for the others to compile with clang
  202. // http://stackoverflow.com/questions/27748442/is-clangs-c11-support-reliable
  203. #endif