points_inside_component.cpp 20 KB

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