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