signed_distance.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292
  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. template <typename Kernel>
  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. assert(V.cols() == 3 && "V should have 3d positions");
  28. assert(P.cols() == 3 && "P should have 3d positions");
  29. assert(F.cols() == 3 && "F should have triangles");
  30. typedef typename Kernel::FT FT;
  31. typedef typename Kernel::Point_3 Point_3;
  32. typedef typename CGAL::Triangle_3<Kernel> Triangle_3;
  33. typedef typename std::vector<Triangle_3>::iterator Iterator;
  34. typedef typename CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  35. typedef typename CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  36. typedef typename CGAL::AABB_tree<AABB_triangle_traits> Tree;
  37. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  38. // Prepare distance computation
  39. Tree tree;
  40. vector<Triangle_3 > T;
  41. point_mesh_squared_distance_precompute(V,F,tree,T);
  42. Eigen::MatrixXd FN,VN,EN;
  43. Eigen::MatrixXi E;
  44. Eigen::VectorXi EMAP;
  45. WindingNumberAABB<Eigen::Vector3d> hier;
  46. switch(sign_type)
  47. {
  48. default:
  49. assert(false && "Unknown SignedDistanceType");
  50. case SIGNED_DISTANCE_TYPE_DEFAULT:
  51. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  52. hier.set_mesh(V,F);
  53. hier.grow();
  54. break;
  55. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  56. // "Signed Distance Computation Using the Angle Weighted Pseudonormal"
  57. // [Bærentzen & Aanæs 2005]
  58. per_face_normals(V,F,FN);
  59. per_vertex_normals(V,F,PER_VERTEX_NORMALS_WEIGHTING_TYPE_ANGLE,FN,VN);
  60. per_edge_normals(
  61. V,F,PER_EDGE_NORMALS_WEIGHTING_TYPE_UNIFORM,FN,EN,E,EMAP);
  62. N.resize(P.rows(),3);
  63. break;
  64. }
  65. S.resize(P.rows(),1);
  66. I.resize(P.rows(),1);
  67. C.resize(P.rows(),3);
  68. for(int p = 0;p<P.rows();p++)
  69. {
  70. const Point_3 q(P(p,0),P(p,1),P(p,2));
  71. double s,sqrd;
  72. Point_and_primitive_id pp;
  73. switch(sign_type)
  74. {
  75. default:
  76. assert(false && "Unknown SignedDistanceType");
  77. case SIGNED_DISTANCE_TYPE_DEFAULT:
  78. case SIGNED_DISTANCE_TYPE_WINDING_NUMBER:
  79. signed_distance_winding_number(tree,hier,q,s,sqrd,pp);
  80. break;
  81. case SIGNED_DISTANCE_TYPE_PSEUDONORMAL:
  82. {
  83. Vector3d n(0,0,0);
  84. signed_distance_pseudonormal(tree,T,F,FN,VN,EN,EMAP,q,s,sqrd,pp,n);
  85. N.row(p) = n;
  86. break;
  87. }
  88. }
  89. I(p) = pp.second - T.begin();
  90. S(p) = s*sqrt(sqrd);
  91. C(p,0) = pp.first.x();
  92. C(p,1) = pp.first.y();
  93. C(p,2) = pp.first.z();
  94. }
  95. }
  96. template <typename Kernel>
  97. IGL_INLINE typename Kernel::FT igl::signed_distance_pseudonormal(
  98. const CGAL::AABB_tree<
  99. CGAL::AABB_traits<Kernel,
  100. CGAL::AABB_triangle_primitive<Kernel,
  101. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  102. >
  103. >
  104. > & tree,
  105. const std::vector<CGAL::Triangle_3<Kernel> > & T,
  106. const Eigen::MatrixXi & F,
  107. const Eigen::MatrixXd & FN,
  108. const Eigen::MatrixXd & VN,
  109. const Eigen::MatrixXd & EN,
  110. const Eigen::VectorXi & EMAP,
  111. const typename Kernel::Point_3 & q)
  112. {
  113. typename CGAL::AABB_tree<
  114. CGAL::AABB_traits<Kernel,
  115. CGAL::AABB_triangle_primitive<Kernel,
  116. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator>
  117. >
  118. >::Point_and_primitive_id pp;
  119. typename Kernel::FT s,sqrd;
  120. Eigen::Vector3d n(0,0,0);
  121. signed_distance_pseudonormal<Kernel>(tree,T,F,FN,VN,EN,EMAP,q,s,sqrd,pp,n);
  122. return s*sqrt(sqrd);
  123. }
  124. template <typename Kernel>
  125. IGL_INLINE void igl::signed_distance_pseudonormal(
  126. const CGAL::AABB_tree<
  127. CGAL::AABB_traits<Kernel,
  128. CGAL::AABB_triangle_primitive<Kernel,
  129. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  130. >
  131. >
  132. > & tree,
  133. const std::vector<CGAL::Triangle_3<Kernel> > & T,
  134. const Eigen::MatrixXi & F,
  135. const Eigen::MatrixXd & FN,
  136. const Eigen::MatrixXd & VN,
  137. const Eigen::MatrixXd & EN,
  138. const Eigen::VectorXi & EMAP,
  139. const typename Kernel::Point_3 & q,
  140. typename Kernel::FT & s,
  141. typename Kernel::FT & sqrd,
  142. typename CGAL::AABB_tree<
  143. CGAL::AABB_traits<Kernel,
  144. CGAL::AABB_triangle_primitive<Kernel,
  145. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator>
  146. >
  147. >::Point_and_primitive_id & pp,
  148. Eigen::Vector3d & n)
  149. {
  150. using namespace Eigen;
  151. using namespace std;
  152. typedef typename Kernel::FT FT;
  153. typedef typename Kernel::Point_3 Point_3;
  154. typedef typename CGAL::Triangle_3<Kernel> Triangle_3;
  155. typedef typename std::vector<Triangle_3>::iterator Iterator;
  156. typedef typename CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  157. typedef typename CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  158. typedef typename CGAL::AABB_tree<AABB_triangle_traits> Tree;
  159. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  160. pp = tree.closest_point_and_primitive(q);
  161. Point_3 & p = pp.first;
  162. const auto & qp = q-p;
  163. sqrd = qp.squared_length();
  164. Vector3d v(qp.x(),qp.y(),qp.z());
  165. const int f = pp.second - T.begin();
  166. const Triangle_3 & t = *pp.second;
  167. // barycentric coordinates
  168. const auto & area = [&p,&t](const int i, const int j)->FT
  169. {
  170. return sqrt(Triangle_3(p,t.vertex(i),t.vertex(j)).squared_area());
  171. };
  172. Vector3d b(area(1,2),area(2,0),area(0,1));
  173. b /= b.sum();
  174. // Determine which normal to use
  175. const double epsilon = 1e-12;
  176. const int type = (b.array()<=epsilon).cast<int>().sum();
  177. switch(type)
  178. {
  179. case 2:
  180. // Find vertex
  181. for(int c = 0;c<3;c++)
  182. {
  183. if(b(c)>epsilon)
  184. {
  185. n = VN.row(F(f,c));
  186. break;
  187. }
  188. }
  189. break;
  190. case 1:
  191. // Find edge
  192. for(int c = 0;c<3;c++)
  193. {
  194. if(b(c)<=epsilon)
  195. {
  196. n = EN.row(EMAP(F.rows()*c+f));
  197. break;
  198. }
  199. }
  200. break;
  201. default:
  202. assert(false && "all barycentric coords zero.");
  203. case 0:
  204. n = FN.row(f);
  205. break;
  206. }
  207. s = (v.dot(n) >= 0 ? 1. : -1.);
  208. }
  209. template <typename Kernel>
  210. IGL_INLINE typename Kernel::FT igl::signed_distance_winding_number(
  211. const CGAL::AABB_tree<
  212. CGAL::AABB_traits<Kernel,
  213. CGAL::AABB_triangle_primitive<Kernel,
  214. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  215. >
  216. >
  217. > & tree,
  218. const igl::WindingNumberAABB<Eigen::Vector3d> & hier,
  219. const typename Kernel::Point_3 & q)
  220. {
  221. typename CGAL::AABB_tree<
  222. CGAL::AABB_traits<Kernel,
  223. CGAL::AABB_triangle_primitive<Kernel,
  224. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator>
  225. >
  226. >::Point_and_primitive_id pp;
  227. typename Kernel::FT s,sqrd;
  228. signed_distance_winding_number<Kernel>(tree,hier,q,s,sqrd,pp);
  229. return s*sqrt(sqrd);
  230. }
  231. template <typename Kernel>
  232. IGL_INLINE void igl::signed_distance_winding_number(
  233. const CGAL::AABB_tree<
  234. CGAL::AABB_traits<Kernel,
  235. CGAL::AABB_triangle_primitive<Kernel,
  236. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator
  237. >
  238. >
  239. > & tree,
  240. const igl::WindingNumberAABB<Eigen::Vector3d> & hier,
  241. const typename Kernel::Point_3 & q,
  242. typename Kernel::FT & s,
  243. typename Kernel::FT & sqrd,
  244. typename CGAL::AABB_tree<
  245. CGAL::AABB_traits<Kernel,
  246. CGAL::AABB_triangle_primitive<Kernel,
  247. typename std::vector<CGAL::Triangle_3<Kernel> >::iterator>
  248. >
  249. >::Point_and_primitive_id & pp)
  250. {
  251. typedef typename Kernel::FT FT;
  252. typedef typename Kernel::Point_3 Point_3;
  253. typedef typename CGAL::Triangle_3<Kernel> Triangle_3;
  254. typedef typename std::vector<Triangle_3>::iterator Iterator;
  255. typedef typename CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  256. typedef typename CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  257. typedef typename CGAL::AABB_tree<AABB_triangle_traits> Tree;
  258. typedef typename Tree::Point_and_primitive_id Point_and_primitive_id;
  259. using namespace Eigen;
  260. using namespace std;
  261. using namespace igl;
  262. pp = tree.closest_point_and_primitive(q);
  263. Point_3 & p = pp.first;
  264. const auto & qp = q-p;
  265. const Vector3d eq(q.x(),q.y(),q.z());
  266. sqrd = qp.squared_length();
  267. const double w = hier.winding_number(eq);
  268. s = 1.-2.*w;
  269. }
  270. #ifdef IGL_STATIC_LIBRARY
  271. // This template is necessary for the others to compile with clang
  272. // http://stackoverflow.com/questions/27748442/is-clangs-c11-support-reliable
  273. template void igl::signed_distance<CGAL::Epick>( const Eigen::MatrixXd & , const Eigen::MatrixXd & , const Eigen::MatrixXi & , const SignedDistanceType , Eigen::VectorXd & , Eigen::VectorXi &, Eigen::MatrixXd & , Eigen::MatrixXd & );
  274. template CGAL::Epick::FT igl::signed_distance_winding_number<CGAL::Epick>(CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Epick, CGAL::AABB_triangle_primitive<CGAL::Epick, std::vector<CGAL::Triangle_3<CGAL::Epick>, std::allocator<CGAL::Triangle_3<CGAL::Epick> > >::iterator, CGAL::Boolean_tag<false> > > > const&, igl::WindingNumberAABB<Eigen::Matrix<double, 3, 1, 0, 3, 1> > const&, CGAL::Epick::Point_3 const&);
  275. template CGAL::Epick::FT igl::signed_distance_pseudonormal<CGAL::Epick>(CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Epick, CGAL::AABB_triangle_primitive<CGAL::Epick, std::vector<CGAL::Triangle_3<CGAL::Epick>, std::allocator<CGAL::Triangle_3<CGAL::Epick> > >::iterator, CGAL::Boolean_tag<false> > > > const&, std::vector<CGAL::Triangle_3<CGAL::Epick>, std::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, CGAL::Epick::Point_3 const&);
  276. template void igl::signed_distance_pseudonormal<CGAL::Simple_cartesian<double> >(CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> > > > const&, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > > const&, Eigen::Matrix<int, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<double, -1, -1, 0, -1, -1> const&, Eigen::Matrix<int, -1, 1, 0, -1, 1> const&, CGAL::Simple_cartesian<double>::Point_3 const&, CGAL::Simple_cartesian<double>::FT&, CGAL::Simple_cartesian<double>::FT&, CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> > > >::Point_and_primitive_id&, Eigen::Matrix<double, 3, 1, 0, 3, 1>&);
  277. template void igl::signed_distance_winding_number<CGAL::Simple_cartesian<double> >(CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> > > > const&, igl::WindingNumberAABB<Eigen::Matrix<double, 3, 1, 0, 3, 1> > const&, CGAL::Simple_cartesian<double>::Point_3 const&, CGAL::Simple_cartesian<double>::FT&, CGAL::Simple_cartesian<double>::FT&, CGAL::AABB_tree<CGAL::AABB_traits<CGAL::Simple_cartesian<double>, CGAL::AABB_triangle_primitive<CGAL::Simple_cartesian<double>, std::vector<CGAL::Triangle_3<CGAL::Simple_cartesian<double> >, std::allocator<CGAL::Triangle_3<CGAL::Simple_cartesian<double> > > >::iterator, CGAL::Boolean_tag<false> > > >::Point_and_primitive_id&);
  278. #endif