signed_distance.cpp 8.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326
  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 "get_seconds.h"
  10. #include "per_edge_normals.h"
  11. #include "per_face_normals.h"
  12. #include "per_vertex_normals.h"
  13. #include "point_mesh_squared_distance.h"
  14. #include "pseudonormal_test.h"
  15. IGL_INLINE void igl::signed_distance(
  16. const Eigen::MatrixXd & P,
  17. const Eigen::MatrixXd & V,
  18. const Eigen::MatrixXi & F,
  19. const SignedDistanceType sign_type,
  20. Eigen::VectorXd & S,
  21. Eigen::VectorXi & I,
  22. Eigen::MatrixXd & C,
  23. Eigen::MatrixXd & N)
  24. {
  25. using namespace Eigen;
  26. using namespace std;
  27. const int dim = V.cols();
  28. assert((V.cols() == 3||V.cols() == 2) && "V should have 3d or 2d positions");
  29. assert((P.cols() == 3||V.cols() == 2) && "P should have 3d or 2d positions");
  30. assert(V.cols() == P.cols() && "V should have same dimension as P");
  31. // Only unsigned distance is supported for non-triangles
  32. if(sign_type != SIGNED_DISTANCE_TYPE_UNSIGNED)
  33. {
  34. assert(F.cols() == dim && "F should have co-dimension 0 simplices");
  35. }
  36. // Prepare distance computation
  37. AABB<MatrixXd,3> tree3;
  38. AABB<MatrixXd,2> tree2;
  39. switch(dim)
  40. {
  41. default:
  42. case 3:
  43. tree3.init(V,F);
  44. break;
  45. case 2:
  46. tree2.init(V,F);
  47. break;
  48. }
  49. Eigen::MatrixXd FN,VN,EN;
  50. Eigen::MatrixXi E;
  51. Eigen::VectorXi EMAP;
  52. WindingNumberAABB<Eigen::Vector3d> hier3;
  53. switch(sign_type)
  54. {
  55. default:
  56. assert(false && "Unknown SignedDistanceType");
  57. case SIGNED_DISTANCE_TYPE_UNSIGNED:
  58. // do nothing
  59. break;
  60. case SIGNED_DISTANCE_TYPE_DEFAULT:
  61. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  62. switch(dim)
  63. {
  64. default:
  65. case 3:
  66. hier3.set_mesh(V,F);
  67. hier3.grow();
  68. break;
  69. case 2:
  70. // no precomp, no hierarchy
  71. break;
  72. }
  73. break;
  74. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  75. switch(dim)
  76. {
  77. default:
  78. case 3:
  79. // "Signed Distance Computation Using the Angle Weighted Pseudonormal"
  80. // [Bærentzen & Aanæs 2005]
  81. per_face_normals(V,F,FN);
  82. per_vertex_normals(V,F,PER_VERTEX_NORMALS_WEIGHTING_TYPE_ANGLE,FN,VN);
  83. per_edge_normals(
  84. V,F,PER_EDGE_NORMALS_WEIGHTING_TYPE_UNIFORM,FN,EN,E,EMAP);
  85. break;
  86. case 2:
  87. FN.resize(F.rows(),2);
  88. VN = MatrixXd::Zero(V.rows(),2);
  89. for(int e = 0;e<F.rows();e++)
  90. {
  91. // rotate edge vector
  92. FN(e,0) = -(V(F(e,1),1)-V(F(e,0),1));
  93. FN(e,1) = (V(F(e,1),0)-V(F(e,0),0));
  94. FN.row(e).normalize();
  95. // add to vertex normal
  96. VN.row(F(e,1)) += FN.row(e);
  97. VN.row(F(e,0)) += FN.row(e);
  98. }
  99. // normalize to average
  100. VN.rowwise().normalize();
  101. break;
  102. }
  103. N.resize(P.rows(),dim);
  104. break;
  105. }
  106. S.resize(P.rows(),1);
  107. I.resize(P.rows(),1);
  108. C.resize(P.rows(),dim);
  109. for(int p = 0;p<P.rows();p++)
  110. {
  111. const RowVectorXd q = P.row(p);
  112. const RowVector3d q3 = q;
  113. const RowVector2d q2 = q;
  114. double s,sqrd;
  115. RowVectorXd c;
  116. RowVector3d c3;
  117. RowVector2d c2;
  118. int i=-1;
  119. switch(sign_type)
  120. {
  121. default:
  122. assert(false && "Unknown SignedDistanceType");
  123. case SIGNED_DISTANCE_TYPE_UNSIGNED:
  124. s = 1.;
  125. sqrd = dim==3?
  126. tree3.squared_distance(V,F,q3,i,c3):
  127. tree2.squared_distance(V,F,q2,i,c2);
  128. break;
  129. case SIGNED_DISTANCE_TYPE_DEFAULT:
  130. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  131. dim==3 ?
  132. signed_distance_winding_number(tree3,V,F,hier3,q3,s,sqrd,i,c3):
  133. signed_distance_winding_number(tree2,V,F,q2,s,sqrd,i,c2);
  134. break;
  135. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  136. {
  137. RowVector3d n3;
  138. RowVector2d n2;
  139. dim==3 ?
  140. signed_distance_pseudonormal(tree3,V,F,FN,VN,EN,EMAP,q3,s,sqrd,i,c3,n3):
  141. signed_distance_pseudonormal(tree2,V,F,FN,VN,q2,s,sqrd,i,c2,n2);
  142. Eigen::RowVectorXd n;
  143. (dim==3 ? n = n3 : n = n3);
  144. N.row(p) = n;
  145. break;
  146. }
  147. }
  148. I(p) = i;
  149. S(p) = s*sqrt(sqrd);
  150. C.row(p) = (dim==3 ? c=c3 : c=c2);
  151. }
  152. }
  153. IGL_INLINE double igl::signed_distance_pseudonormal(
  154. const AABB<Eigen::MatrixXd,3> & tree,
  155. const Eigen::MatrixXd & V,
  156. const Eigen::MatrixXi & F,
  157. const Eigen::MatrixXd & FN,
  158. const Eigen::MatrixXd & VN,
  159. const Eigen::MatrixXd & EN,
  160. const Eigen::VectorXi & EMAP,
  161. const Eigen::RowVector3d & q)
  162. {
  163. double s,sqrd;
  164. Eigen::RowVector3d n,c;
  165. int i = -1;
  166. signed_distance_pseudonormal(tree,V,F,FN,VN,EN,EMAP,q,s,sqrd,i,c,n);
  167. return s*sqrt(sqrd);
  168. }
  169. IGL_INLINE void igl::signed_distance_pseudonormal(
  170. const Eigen::MatrixXd & P,
  171. const Eigen::MatrixXd & V,
  172. const Eigen::MatrixXi & F,
  173. const AABB<Eigen::MatrixXd,3> & tree,
  174. const Eigen::MatrixXd & FN,
  175. const Eigen::MatrixXd & VN,
  176. const Eigen::MatrixXd & EN,
  177. const Eigen::VectorXi & EMAP,
  178. Eigen::VectorXd & S,
  179. Eigen::VectorXi & I,
  180. Eigen::MatrixXd & C,
  181. Eigen::MatrixXd & N)
  182. {
  183. using namespace Eigen;
  184. const size_t np = P.rows();
  185. S.resize(np,1);
  186. I.resize(np,1);
  187. N.resize(np,3);
  188. C.resize(np,3);
  189. # pragma omp parallel for if(np>1000)
  190. for(size_t p = 0;p<np;p++)
  191. {
  192. double s,sqrd;
  193. RowVector3d n,c;
  194. int i = -1;
  195. RowVector3d q = P.row(p);
  196. signed_distance_pseudonormal(tree,V,F,FN,VN,EN,EMAP,q,s,sqrd,i,c,n);
  197. S(p) = s*sqrt(sqrd);
  198. I(p) = i;
  199. N.row(p) = n;
  200. C.row(p) = c;
  201. }
  202. // igl::AABB<MatrixXd,3> tree_P;
  203. // MatrixXi J = VectorXi::LinSpaced(P.rows(),0,P.rows()-1);
  204. // tree_P.init(P,J);
  205. // tree.squared_distance(V,F,tree_P,P,J,S,I,C);
  206. //# pragma omp parallel for if(np>1000)
  207. // for(size_t p = 0;p<np;p++)
  208. // {
  209. // RowVector3d c = C.row(p);
  210. // RowVector3d q = P.row(p);
  211. // const int f = I(p);
  212. // double s;
  213. // RowVector3d n;
  214. // pseudonormal_test(V,F,FN,VN,EN,EMAP,q,f,c,s,n);
  215. // N.row(p) = n;
  216. // S(p) = s*sqrt(S(p));
  217. // }
  218. }
  219. IGL_INLINE void igl::signed_distance_pseudonormal(
  220. const AABB<Eigen::MatrixXd,3> & tree,
  221. const Eigen::MatrixXd & V,
  222. const Eigen::MatrixXi & F,
  223. const Eigen::MatrixXd & FN,
  224. const Eigen::MatrixXd & VN,
  225. const Eigen::MatrixXd & EN,
  226. const Eigen::VectorXi & EMAP,
  227. const Eigen::RowVector3d & q,
  228. double & s,
  229. double & sqrd,
  230. int & f,
  231. Eigen::RowVector3d & c,
  232. Eigen::RowVector3d & n)
  233. {
  234. using namespace Eigen;
  235. using namespace std;
  236. sqrd = tree.squared_distance(V,F,q,f,c);
  237. pseudonormal_test(V,F,FN,VN,EN,EMAP,q,f,c,s,n);
  238. }
  239. IGL_INLINE void igl::signed_distance_pseudonormal(
  240. const AABB<Eigen::MatrixXd,2> & tree,
  241. const Eigen::MatrixXd & V,
  242. const Eigen::MatrixXi & F,
  243. const Eigen::MatrixXd & FN,
  244. const Eigen::MatrixXd & VN,
  245. const Eigen::RowVector2d & q,
  246. double & s,
  247. double & sqrd,
  248. int & f,
  249. Eigen::RowVector2d & c,
  250. Eigen::RowVector2d & n)
  251. {
  252. using namespace Eigen;
  253. using namespace std;
  254. sqrd = tree.squared_distance(V,F,q,f,c);
  255. pseudonormal_test(V,F,FN,VN,q,f,c,s,n);
  256. }
  257. IGL_INLINE double igl::signed_distance_winding_number(
  258. const AABB<Eigen::MatrixXd,3> & tree,
  259. const Eigen::MatrixXd & V,
  260. const Eigen::MatrixXi & F,
  261. const igl::WindingNumberAABB<Eigen::Vector3d> & hier,
  262. const Eigen::RowVector3d & q)
  263. {
  264. double s,sqrd;
  265. Eigen::RowVector3d c;
  266. int i=-1;
  267. signed_distance_winding_number(tree,V,F,hier,q,s,sqrd,i,c);
  268. return s*sqrt(sqrd);
  269. }
  270. IGL_INLINE void igl::signed_distance_winding_number(
  271. const AABB<Eigen::MatrixXd,3> & tree,
  272. const Eigen::MatrixXd & V,
  273. const Eigen::MatrixXi & F,
  274. const igl::WindingNumberAABB<Eigen::Matrix<double,3,1> > & hier,
  275. const Eigen::Matrix<double,1,3> & q,
  276. double & s,
  277. double & sqrd,
  278. int & i,
  279. Eigen::Matrix<double,1,3> & c)
  280. {
  281. using namespace Eigen;
  282. using namespace std;
  283. using namespace igl;
  284. sqrd = tree.squared_distance(V,F,q,i,c);
  285. const double w = hier.winding_number(q.transpose());
  286. s = 1.-2.*w;
  287. }
  288. IGL_INLINE void igl::signed_distance_winding_number(
  289. const AABB<Eigen::MatrixXd,2> & tree,
  290. const Eigen::MatrixXd & V,
  291. const Eigen::MatrixXi & F,
  292. const Eigen::Matrix<double,1,2> & q,
  293. double & s,
  294. double & sqrd,
  295. int & i,
  296. Eigen::Matrix<double,1,2> & c)
  297. {
  298. using namespace Eigen;
  299. using namespace std;
  300. using namespace igl;
  301. sqrd = tree.squared_distance(V,F,q,i,c);
  302. double w;
  303. winding_number_2(V.data(), V.rows(), F.data(), F.rows(), q.data(), 1, &w);
  304. s = 1.-2.*w;
  305. }
  306. #ifdef IGL_STATIC_LIBRARY
  307. // This template is necessary for the others to compile with clang
  308. // http://stackoverflow.com/questions/27748442/is-clangs-c11-support-reliable
  309. #endif