points_inside_component.cpp 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Qingnan Zhou <qnzhou@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 "points_inside_component.h"
  9. #include <cassert>
  10. #include <list>
  11. #include <limits>
  12. #include <vector>
  13. #include <CGAL/AABB_tree.h>
  14. #include <CGAL/AABB_traits.h>
  15. #include <CGAL/AABB_triangle_primitive.h>
  16. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  17. #include "order_facets_around_edge.h"
  18. #include "assign_scalar.h"
  19. namespace igl {
  20. namespace cgal {
  21. namespace points_inside_component_helper {
  22. typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
  23. typedef Kernel::Ray_3 Ray_3;
  24. typedef Kernel::Point_3 Point_3;
  25. typedef Kernel::Vector_3 Vector_3;
  26. typedef Kernel::Triangle_3 Triangle;
  27. typedef Kernel::Plane_3 Plane_3;
  28. typedef std::vector<Triangle>::iterator Iterator;
  29. typedef CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  30. typedef CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  31. typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
  32. enum ElementType { VERTEX, EDGE, FACE };
  33. template<typename DerivedV, typename DerivedF, typename DerivedI>
  34. ElementType determine_element_type(
  35. const Eigen::PlainObjectBase<DerivedV>& V,
  36. const Eigen::PlainObjectBase<DerivedF>& F,
  37. const Eigen::PlainObjectBase<DerivedI>& I,
  38. const size_t fid, const Point_3& p,
  39. size_t& element_index) {
  40. const Eigen::Vector3i f = F.row(I(fid, 0));
  41. const Point_3 p0(V(f[0], 0), V(f[0], 1), V(f[0], 2));
  42. const Point_3 p1(V(f[1], 0), V(f[1], 1), V(f[1], 2));
  43. const Point_3 p2(V(f[2], 0), V(f[2], 1), V(f[2], 2));
  44. if (p == p0) { element_index = 0; return VERTEX; }
  45. if (p == p1) { element_index = 1; return VERTEX; }
  46. if (p == p2) { element_index = 2; return VERTEX; }
  47. if (CGAL::collinear(p0, p1, p)) { element_index = 2; return EDGE; }
  48. if (CGAL::collinear(p1, p2, p)) { element_index = 0; return EDGE; }
  49. if (CGAL::collinear(p2, p0, p)) { element_index = 1; return EDGE; }
  50. element_index = 0;
  51. return FACE;
  52. }
  53. template<typename DerivedF, typename DerivedI>
  54. void extract_adj_faces(
  55. const Eigen::PlainObjectBase<DerivedF>& F,
  56. const Eigen::PlainObjectBase<DerivedI>& I,
  57. const size_t s, const size_t d,
  58. std::vector<int>& adj_faces) {
  59. const size_t num_faces = I.rows();
  60. for (size_t i=0; i<num_faces; i++) {
  61. Eigen::Vector3i f = F.row(I(i, 0));
  62. if (((size_t)f[0] == s && (size_t)f[1] == d) ||
  63. ((size_t)f[1] == s && (size_t)f[2] == d) ||
  64. ((size_t)f[2] == s && (size_t)f[0] == d)) {
  65. adj_faces.push_back((I(i, 0)+1) * -1);
  66. continue;
  67. }
  68. if (((size_t)f[0] == d && (size_t)f[1] == s) ||
  69. ((size_t)f[1] == d && (size_t)f[2] == s) ||
  70. ((size_t)f[2] == d && (size_t)f[0] == s)) {
  71. adj_faces.push_back(I(i, 0)+1);
  72. continue;
  73. }
  74. }
  75. }
  76. template<typename DerivedF, typename DerivedI>
  77. void extract_adj_vertices(
  78. const Eigen::PlainObjectBase<DerivedF>& F,
  79. const Eigen::PlainObjectBase<DerivedI>& I,
  80. const size_t v, std::vector<int>& adj_vertices) {
  81. std::set<size_t> unique_adj_vertices;
  82. const size_t num_faces = I.rows();
  83. for (size_t i=0; i<num_faces; i++) {
  84. Eigen::Vector3i f = F.row(I(i, 0));
  85. if ((size_t)f[0] == v) {
  86. unique_adj_vertices.insert(f[1]);
  87. unique_adj_vertices.insert(f[2]);
  88. } else if ((size_t)f[1] == v) {
  89. unique_adj_vertices.insert(f[0]);
  90. unique_adj_vertices.insert(f[2]);
  91. } else if ((size_t)f[2] == v) {
  92. unique_adj_vertices.insert(f[0]);
  93. unique_adj_vertices.insert(f[1]);
  94. }
  95. }
  96. adj_vertices.resize(unique_adj_vertices.size());
  97. std::copy(unique_adj_vertices.begin(),
  98. unique_adj_vertices.end(),
  99. adj_vertices.begin());
  100. }
  101. template<typename DerivedV, typename DerivedF, typename DerivedI>
  102. bool determine_point_edge_orientation(
  103. const Eigen::PlainObjectBase<DerivedV>& V,
  104. const Eigen::PlainObjectBase<DerivedF>& F,
  105. const Eigen::PlainObjectBase<DerivedI>& I,
  106. const Point_3& query, size_t s, size_t d) {
  107. // Algorithm:
  108. //
  109. // Order the adj faces around the edge (s,d) clockwise using
  110. // query point as pivot. (i.e. The first face of the ordering
  111. // is directly after the pivot point, and the last face is
  112. // directly before the pivot.)
  113. //
  114. // The point is outside if the first and last faces of the
  115. // ordering forms a convex angle. This check can be done
  116. // without any construction by looking at the orientation of the
  117. // faces. The angle is convex iff the first face contains (s,d)
  118. // as an edge and the last face contains (d,s) as an edge.
  119. //
  120. // The point is inside if the first and last faces of the
  121. // ordering forms a concave angle. That is the first face
  122. // contains (d,s) as an edge and the last face contains (s,d) as
  123. // an edge.
  124. //
  125. // In the special case of duplicated faces. I.e. multiple faces
  126. // sharing the same 3 corners, but not necessarily the same
  127. // orientation. The ordering will always rank faces containing
  128. // edge (s,d) before faces containing edge (d,s).
  129. //
  130. // Therefore, if there are any duplicates of the first faces,
  131. // the ordering will always choose the one with edge (s,d) if
  132. // possible. The same for the last face.
  133. //
  134. // In the very degenerated case where the first and last face
  135. // are duplicates, but with different orientations, it is
  136. // equally valid to think the angle formed by them is either 0
  137. // or 360 degrees. By default, 0 degree is used, and thus the
  138. // query point is outside.
  139. std::vector<int> adj_faces;
  140. extract_adj_faces(F, I, s, d, adj_faces);
  141. const size_t num_adj_faces = adj_faces.size();
  142. assert(num_adj_faces > 0);
  143. DerivedV pivot_point(1, 3);
  144. igl::cgal::assign_scalar(query.x(), pivot_point(0, 0));
  145. igl::cgal::assign_scalar(query.y(), pivot_point(0, 1));
  146. igl::cgal::assign_scalar(query.z(), pivot_point(0, 2));
  147. Eigen::VectorXi order;
  148. order_facets_around_edge(V, F, s, d,
  149. adj_faces, pivot_point, order);
  150. assert((size_t)order.size() == num_adj_faces);
  151. if (adj_faces[order[0]] > 0 &&
  152. adj_faces[order[num_adj_faces-1] < 0]) {
  153. return true;
  154. } else if (adj_faces[order[0]] < 0 &&
  155. adj_faces[order[num_adj_faces-1] > 0]) {
  156. return false;
  157. } else {
  158. throw "The input mesh does not represent a valid volume";
  159. }
  160. throw "The input mesh does not represent a valid volume";
  161. return false;
  162. }
  163. template<typename DerivedV, typename DerivedF, typename DerivedI>
  164. bool determine_point_vertex_orientation(
  165. const Eigen::PlainObjectBase<DerivedV>& V,
  166. const Eigen::PlainObjectBase<DerivedF>& F,
  167. const Eigen::PlainObjectBase<DerivedI>& I,
  168. const Point_3& query, size_t s) {
  169. std::vector<int> adj_vertices;
  170. extract_adj_vertices(F, I, s, adj_vertices);
  171. const size_t num_adj_vertices = adj_vertices.size();
  172. std::vector<Point_3> adj_points;
  173. for (size_t i=0; i<num_adj_vertices; i++) {
  174. const size_t vi = adj_vertices[i];
  175. adj_points.emplace_back(V(vi,0), V(vi,1), V(vi,2));
  176. }
  177. // A plane is on the exterior if all adj_points lies on or to
  178. // one side of the plane.
  179. auto is_on_exterior = [&](const Plane_3& separator) {
  180. size_t positive=0;
  181. size_t negative=0;
  182. size_t coplanar=0;
  183. for (const auto& point : adj_points) {
  184. switch(separator.oriented_side(point)) {
  185. case CGAL::ON_POSITIVE_SIDE:
  186. positive++;
  187. break;
  188. case CGAL::ON_NEGATIVE_SIDE:
  189. negative++;
  190. break;
  191. case CGAL::ON_ORIENTED_BOUNDARY:
  192. coplanar++;
  193. break;
  194. default:
  195. throw "Unknown plane-point orientation";
  196. }
  197. }
  198. auto query_orientation = separator.oriented_side(query);
  199. bool r =
  200. (positive == 0 && query_orientation == CGAL::POSITIVE)
  201. ||
  202. (negative == 0 && query_orientation == CGAL::NEGATIVE);
  203. return r;
  204. };
  205. size_t d = std::numeric_limits<size_t>::max();
  206. Point_3 p(V(s,0), V(s,1), V(s,2));
  207. for (size_t i=0; i<num_adj_vertices; i++) {
  208. const size_t vi = adj_vertices[i];
  209. for (size_t j=i+1; j<num_adj_vertices; j++) {
  210. Plane_3 separator(p, adj_points[i], adj_points[j]);
  211. if (separator.is_degenerate()) {
  212. throw "Input mesh contains degenerated faces";
  213. }
  214. if (is_on_exterior(separator)) {
  215. d = vi;
  216. assert(!CGAL::collinear(p, adj_points[i], query));
  217. break;
  218. }
  219. }
  220. if (d < (size_t)V.rows()) break;
  221. }
  222. if (d > (size_t)V.rows()) {
  223. // All adj faces are coplanar, use the first edge.
  224. d = adj_vertices[0];
  225. }
  226. return determine_point_edge_orientation(V, F, I, query, s, d);
  227. }
  228. template<typename DerivedV, typename DerivedF, typename DerivedI>
  229. bool determine_point_face_orientation(
  230. const Eigen::PlainObjectBase<DerivedV>& V,
  231. const Eigen::PlainObjectBase<DerivedF>& F,
  232. const Eigen::PlainObjectBase<DerivedI>& I,
  233. const Point_3& query, size_t fid) {
  234. // Algorithm: A point is on the inside of a face if the
  235. // tetrahedron formed by them is negatively oriented.
  236. Eigen::Vector3i f = F.row(I(fid, 0));
  237. const Point_3 v0(V(f[0], 0), V(f[0], 1), V(f[0], 2));
  238. const Point_3 v1(V(f[1], 0), V(f[1], 1), V(f[1], 2));
  239. const Point_3 v2(V(f[2], 0), V(f[2], 1), V(f[2], 2));
  240. auto result = CGAL::orientation(v0, v1, v2, query);
  241. if (result == CGAL::COPLANAR) {
  242. throw "Cannot determine inside/outside because query point lies exactly on the input surface.";
  243. }
  244. return result == CGAL::NEGATIVE;
  245. }
  246. }
  247. }
  248. }
  249. template<typename DerivedV, typename DerivedF, typename DerivedI,
  250. typename DerivedP, typename DerivedB>
  251. IGL_INLINE void igl::cgal::points_inside_component(
  252. const Eigen::PlainObjectBase<DerivedV>& V,
  253. const Eigen::PlainObjectBase<DerivedF>& F,
  254. const Eigen::PlainObjectBase<DerivedI>& I,
  255. const Eigen::PlainObjectBase<DerivedP>& P,
  256. Eigen::PlainObjectBase<DerivedB>& inside) {
  257. using namespace igl::cgal::points_inside_component_helper;
  258. if (F.rows() <= 0 || I.rows() <= 0) {
  259. throw "Inside check cannot be done on empty facet component.";
  260. }
  261. const size_t num_faces = I.rows();
  262. std::vector<Triangle> triangles;
  263. for (size_t i=0; i<num_faces; i++) {
  264. const Eigen::Vector3i f = F.row(I(i, 0));
  265. triangles.emplace_back(
  266. Point_3(V(f[0], 0), V(f[0], 1), V(f[0], 2)),
  267. Point_3(V(f[1], 0), V(f[1], 1), V(f[1], 2)),
  268. Point_3(V(f[2], 0), V(f[2], 1), V(f[2], 2)));
  269. if (triangles.back().is_degenerate()) {
  270. throw "Input facet components contains degenerated triangles";
  271. }
  272. }
  273. Tree tree(triangles.begin(), triangles.end());
  274. tree.accelerate_distance_queries();
  275. const size_t num_queries = P.rows();
  276. inside.resize(num_queries, 1);
  277. for (size_t i=0; i<num_queries; i++) {
  278. const Point_3 query(P(i,0), P(i,1), P(i,2));
  279. auto projection = tree.closest_point_and_primitive(query);
  280. auto closest_point = projection.first;
  281. size_t fid = projection.second - triangles.begin();
  282. size_t element_index;
  283. switch (determine_element_type(
  284. V, F, I, fid, closest_point, element_index)) {
  285. case VERTEX:
  286. {
  287. const size_t s = F(I(fid, 0), element_index);
  288. inside(i,0) = determine_point_vertex_orientation(
  289. V, F, I, query, s);
  290. }
  291. break;
  292. case EDGE:
  293. {
  294. const size_t s = F(I(fid, 0), (element_index+1)%3);
  295. const size_t d = F(I(fid, 0), (element_index+2)%3);
  296. inside(i,0) = determine_point_edge_orientation(
  297. V, F, I, query, s, d);
  298. }
  299. break;
  300. case FACE:
  301. inside(i,0) = determine_point_face_orientation(V, F, I, query, fid);
  302. break;
  303. default:
  304. throw "Unknow closest element type!";
  305. }
  306. }
  307. }
  308. template<typename DerivedV, typename DerivedF, typename DerivedP,
  309. typename DerivedB>
  310. IGL_INLINE void igl::cgal::points_inside_component(
  311. const Eigen::PlainObjectBase<DerivedV>& V,
  312. const Eigen::PlainObjectBase<DerivedF>& F,
  313. const Eigen::PlainObjectBase<DerivedP>& P,
  314. Eigen::PlainObjectBase<DerivedB>& inside) {
  315. Eigen::VectorXi I(F.rows());
  316. I.setLinSpaced(F.rows(), 0, F.rows()-1);
  317. igl::cgal::points_inside_component(V, F, I, P, inside);
  318. }
  319. #ifdef IGL_STATIC_LIBRARY
  320. // Explicit template specialization
  321. template void igl::cgal::points_inside_component< Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix< int, -1, -1, 0, -1, -1>, Eigen::Matrix< int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix< int, -1, -1, 0, -1, -1> > ( Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix< int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix< int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix< int, -1, -1, 0, -1, -1> >&);
  322. template void igl::cgal::points_inside_component< Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix< int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix< int, -1, -1, 0, -1, -1> > ( Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix< int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix< int, -1, -1, 0, -1, -1> >&);
  323. template void igl::cgal::points_inside_component<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  324. template void igl::cgal::points_inside_component<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  325. template void igl::cgal::points_inside_component<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  326. #endif