closest_facet.cpp 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360
  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. //
  9. #include "closest_facet.h"
  10. #include <vector>
  11. #include <stdexcept>
  12. #include <CGAL/AABB_tree.h>
  13. #include <CGAL/AABB_traits.h>
  14. #include <CGAL/AABB_triangle_primitive.h>
  15. #include <CGAL/intersections.h>
  16. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  17. #include "order_facets_around_edge.h"
  18. template<
  19. typename DerivedV,
  20. typename DerivedF,
  21. typename DerivedI,
  22. typename DerivedP,
  23. typename DerivedR,
  24. typename DerivedS >
  25. IGL_INLINE void igl::copyleft::cgal::closest_facet(
  26. const Eigen::PlainObjectBase<DerivedV>& V,
  27. const Eigen::PlainObjectBase<DerivedF>& F,
  28. const Eigen::PlainObjectBase<DerivedI>& I,
  29. const Eigen::PlainObjectBase<DerivedP>& P,
  30. Eigen::PlainObjectBase<DerivedR>& R,
  31. Eigen::PlainObjectBase<DerivedS>& S) {
  32. typedef CGAL::Exact_predicates_exact_constructions_kernel Kernel;
  33. typedef Kernel::Point_3 Point_3;
  34. typedef Kernel::Plane_3 Plane_3;
  35. typedef Kernel::Segment_3 Segment_3;
  36. typedef Kernel::Triangle_3 Triangle;
  37. typedef std::vector<Triangle>::iterator Iterator;
  38. typedef CGAL::AABB_triangle_primitive<Kernel, Iterator> Primitive;
  39. typedef CGAL::AABB_traits<Kernel, Primitive> AABB_triangle_traits;
  40. typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
  41. if (F.rows() <= 0 || I.rows() <= 0) {
  42. throw std::runtime_error(
  43. "Closest facet cannot be computed on empty mesh.");
  44. }
  45. const size_t num_faces = I.rows();
  46. std::vector<Triangle> triangles;
  47. for (size_t i=0; i<num_faces; i++) {
  48. const Eigen::Vector3i f = F.row(I(i, 0));
  49. triangles.emplace_back(
  50. Point_3(V(f[0], 0), V(f[0], 1), V(f[0], 2)),
  51. Point_3(V(f[1], 0), V(f[1], 1), V(f[1], 2)),
  52. Point_3(V(f[2], 0), V(f[2], 1), V(f[2], 2)));
  53. if (triangles.back().is_degenerate()) {
  54. throw std::runtime_error(
  55. "Input facet components contains degenerated triangles");
  56. }
  57. }
  58. Tree tree(triangles.begin(), triangles.end());
  59. tree.accelerate_distance_queries();
  60. auto on_the_positive_side = [&](size_t fid, const Point_3& p) {
  61. const auto& f = F.row(fid).eval();
  62. Point_3 v0(V(f[0], 0), V(f[0], 1), V(f[0], 2));
  63. Point_3 v1(V(f[1], 0), V(f[1], 1), V(f[1], 2));
  64. Point_3 v2(V(f[2], 0), V(f[2], 1), V(f[2], 2));
  65. auto ori = CGAL::orientation(v0, v1, v2, p);
  66. switch (ori) {
  67. case CGAL::POSITIVE:
  68. return true;
  69. case CGAL::NEGATIVE:
  70. return false;
  71. case CGAL::COPLANAR:
  72. throw std::runtime_error(
  73. "It seems input mesh contains self intersection");
  74. default:
  75. throw std::runtime_error("Unknown CGAL state.");
  76. }
  77. return false;
  78. };
  79. auto get_orientation = [&](size_t fid, size_t s, size_t d) -> bool {
  80. const auto& f = F.row(fid);
  81. if ((size_t)f[0] == s && (size_t)f[1] == d) return false;
  82. else if ((size_t)f[1] == s && (size_t)f[2] == d) return false;
  83. else if ((size_t)f[2] == s && (size_t)f[0] == d) return false;
  84. else if ((size_t)f[0] == d && (size_t)f[1] == s) return true;
  85. else if ((size_t)f[1] == d && (size_t)f[2] == s) return true;
  86. else if ((size_t)f[2] == d && (size_t)f[0] == s) return true;
  87. else {
  88. throw std::runtime_error(
  89. "Cannot compute orientation due to incorrect connectivity");
  90. return false;
  91. }
  92. };
  93. auto index_to_signed_index = [&](size_t index, bool ori) -> int{
  94. return (index+1) * (ori? 1:-1);
  95. };
  96. //auto signed_index_to_index = [&](int signed_index) -> size_t {
  97. // return abs(signed_index) - 1;
  98. //};
  99. enum ElementType { VERTEX, EDGE, FACE };
  100. auto determine_element_type = [&](const Point_3& p, const size_t fid,
  101. size_t& element_index) {
  102. const auto& tri = triangles[fid];
  103. const Point_3 p0 = tri[0];
  104. const Point_3 p1 = tri[1];
  105. const Point_3 p2 = tri[2];
  106. if (p == p0) { element_index = 0; return VERTEX; }
  107. if (p == p1) { element_index = 1; return VERTEX; }
  108. if (p == p2) { element_index = 2; return VERTEX; }
  109. if (CGAL::collinear(p0, p1, p)) { element_index = 2; return EDGE; }
  110. if (CGAL::collinear(p1, p2, p)) { element_index = 0; return EDGE; }
  111. if (CGAL::collinear(p2, p0, p)) { element_index = 1; return EDGE; }
  112. element_index = 0;
  113. return FACE;
  114. };
  115. auto process_edge_case = [&](
  116. size_t query_idx,
  117. const size_t s, const size_t d,
  118. size_t preferred_facet,
  119. bool& orientation) {
  120. Point_3 mid_edge_point(
  121. (V(s,0) + V(d,0)) * 0.5,
  122. (V(s,1) + V(d,1)) * 0.5,
  123. (V(s,2) + V(d,2)) * 0.5);
  124. Point_3 query_point(
  125. P(query_idx, 0),
  126. P(query_idx, 1),
  127. P(query_idx, 2));
  128. std::vector<Tree::Primitive_id> intersected_faces;
  129. tree.all_intersected_primitives(Segment_3(mid_edge_point, query_point),
  130. std::back_inserter(intersected_faces));
  131. const size_t num_intersected_faces = intersected_faces.size();
  132. std::vector<size_t> intersected_face_indices(num_intersected_faces);
  133. std::vector<int> intersected_face_signed_indices(num_intersected_faces);
  134. std::transform(intersected_faces.begin(),
  135. intersected_faces.end(),
  136. intersected_face_indices.begin(),
  137. [&](const Tree::Primitive_id& itr) -> size_t
  138. { return I(itr-triangles.begin(), 0); });
  139. std::transform(
  140. intersected_face_indices.begin(),
  141. intersected_face_indices.end(),
  142. intersected_face_signed_indices.begin(),
  143. [&](size_t index) {
  144. return index_to_signed_index(
  145. index, get_orientation(index, s,d));
  146. });
  147. assert(num_intersected_faces >= 1);
  148. if (num_intersected_faces == 1) {
  149. // The edge must be a boundary edge. Thus, the orientation can be
  150. // simply determined by checking if the query point is on the
  151. // positive side of the facet.
  152. const size_t fid = intersected_face_indices[0];
  153. orientation = on_the_positive_side(fid, query_point);
  154. return fid;
  155. }
  156. Eigen::VectorXi order;
  157. DerivedP pivot = P.row(query_idx).eval();
  158. igl::copyleft::cgal::order_facets_around_edge(V, F, s, d,
  159. intersected_face_signed_indices,
  160. pivot, order);
  161. // Although first and last are equivalent, make the choice based on
  162. // preferred_facet.
  163. const size_t first = order[0];
  164. const size_t last = order[num_intersected_faces-1];
  165. if (intersected_face_indices[first] == preferred_facet) {
  166. orientation = intersected_face_signed_indices[first] < 0;
  167. return intersected_face_indices[first];
  168. } else if (intersected_face_indices[last] == preferred_facet) {
  169. orientation = intersected_face_signed_indices[last] > 0;
  170. return intersected_face_indices[last];
  171. } else {
  172. orientation = intersected_face_signed_indices[order[0]] < 0;
  173. return intersected_face_indices[order[0]];
  174. }
  175. };
  176. auto process_face_case = [&](
  177. const size_t query_idx, const size_t fid, bool& orientation) {
  178. const auto& f = F.row(I(fid, 0));
  179. return process_edge_case(query_idx, f[0], f[1], I(fid, 0), orientation);
  180. };
  181. auto process_vertex_case = [&](const size_t query_idx, size_t s,
  182. size_t preferred_facet, bool& orientation) {
  183. Point_3 closest_point(V(s, 0), V(s, 1), V(s, 2));
  184. Point_3 query_point(P(query_idx, 0), P(query_idx, 1), P(query_idx, 2));
  185. std::vector<Tree::Primitive_id> intersected_faces;
  186. tree.all_intersected_primitives(Segment_3(closest_point, query_point),
  187. std::back_inserter(intersected_faces));
  188. const size_t num_intersected_faces = intersected_faces.size();
  189. std::vector<size_t> intersected_face_indices(num_intersected_faces);
  190. std::transform(intersected_faces.begin(),
  191. intersected_faces.end(),
  192. intersected_face_indices.begin(),
  193. [&](const Tree::Primitive_id& itr) -> size_t
  194. { return I(itr-triangles.begin(), 0); });
  195. std::set<size_t> adj_vertices_set;
  196. for (auto fid : intersected_face_indices) {
  197. const auto& f = F.row(fid);
  198. if ((size_t)f[0] != s) adj_vertices_set.insert(f[0]);
  199. if ((size_t)f[1] != s) adj_vertices_set.insert(f[1]);
  200. if ((size_t)f[2] != s) adj_vertices_set.insert(f[2]);
  201. }
  202. const size_t num_adj_vertices = adj_vertices_set.size();
  203. std::vector<size_t> adj_vertices(num_adj_vertices);
  204. std::copy(adj_vertices_set.begin(), adj_vertices_set.end(),
  205. adj_vertices.begin());
  206. std::vector<Point_3> adj_points;
  207. for (size_t vid : adj_vertices) {
  208. adj_points.emplace_back(V(vid,0), V(vid,1), V(vid,2));
  209. }
  210. // A plane is on the exterior if all adj_points lies on or to
  211. // one side of the plane.
  212. auto is_on_exterior = [&](const Plane_3& separator) {
  213. size_t positive=0;
  214. size_t negative=0;
  215. size_t coplanar=0;
  216. for (const auto& point : adj_points) {
  217. switch(separator.oriented_side(point)) {
  218. case CGAL::ON_POSITIVE_SIDE:
  219. positive++;
  220. break;
  221. case CGAL::ON_NEGATIVE_SIDE:
  222. negative++;
  223. break;
  224. case CGAL::ON_ORIENTED_BOUNDARY:
  225. coplanar++;
  226. break;
  227. default:
  228. throw "Unknown plane-point orientation";
  229. }
  230. }
  231. auto query_orientation = separator.oriented_side(query_point);
  232. bool r = (positive == 0 && query_orientation == CGAL::POSITIVE)
  233. || (negative == 0 && query_orientation == CGAL::NEGATIVE);
  234. return r;
  235. };
  236. size_t d = std::numeric_limits<size_t>::max();
  237. for (size_t i=0; i<num_adj_vertices; i++) {
  238. const size_t vi = adj_vertices[i];
  239. for (size_t j=i+1; j<num_adj_vertices; j++) {
  240. Plane_3 separator(closest_point, adj_points[i], adj_points[j]);
  241. if (separator.is_degenerate()) {
  242. throw std::runtime_error(
  243. "Input mesh contains degenerated faces");
  244. }
  245. if (is_on_exterior(separator)) {
  246. d = vi;
  247. assert(!CGAL::collinear(
  248. query_point, adj_points[i], closest_point));
  249. break;
  250. }
  251. }
  252. }
  253. assert(d != std::numeric_limits<size_t>::max());
  254. return process_edge_case(query_idx, s, d, preferred_facet, orientation);
  255. };
  256. const size_t num_queries = P.rows();
  257. R.resize(num_queries, 1);
  258. S.resize(num_queries, 1);
  259. for (size_t i=0; i<num_queries; i++) {
  260. const Point_3 query(P(i,0), P(i,1), P(i,2));
  261. auto projection = tree.closest_point_and_primitive(query);
  262. const Point_3 closest_point = projection.first;
  263. size_t fid = projection.second - triangles.begin();
  264. bool fid_ori = false;
  265. // Gether all facets sharing the closest point.
  266. std::vector<Tree::Primitive_id> intersected_faces;
  267. tree.all_intersected_primitives(Segment_3(closest_point, query),
  268. std::back_inserter(intersected_faces));
  269. const size_t num_intersected_faces = intersected_faces.size();
  270. std::vector<size_t> intersected_face_indices(num_intersected_faces);
  271. std::transform(intersected_faces.begin(),
  272. intersected_faces.end(),
  273. intersected_face_indices.begin(),
  274. [&](const Tree::Primitive_id& itr) -> size_t
  275. { return I(itr-triangles.begin(), 0); });
  276. size_t element_index;
  277. auto element_type = determine_element_type(closest_point, fid,
  278. element_index);
  279. switch(element_type) {
  280. case VERTEX:
  281. {
  282. const auto& f = F.row(I(fid, 0));
  283. const size_t s = f[element_index];
  284. fid = process_vertex_case(i, s, I(fid, 0), fid_ori);
  285. }
  286. break;
  287. case EDGE:
  288. {
  289. const auto& f = F.row(I(fid, 0));
  290. const size_t s = f[(element_index+1)%3];
  291. const size_t d = f[(element_index+2)%3];
  292. fid = process_edge_case(i, s, d, I(fid, 0), fid_ori);
  293. }
  294. break;
  295. case FACE:
  296. {
  297. fid = process_face_case(i, fid, fid_ori);
  298. }
  299. break;
  300. default:
  301. throw std::runtime_error("Unknown element type.");
  302. }
  303. R(i,0) = fid;
  304. S(i,0) = fid_ori;
  305. }
  306. }
  307. template<
  308. typename DerivedV,
  309. typename DerivedF,
  310. typename DerivedP,
  311. typename DerivedR,
  312. typename DerivedS >
  313. IGL_INLINE void igl::copyleft::cgal::closest_facet(
  314. const Eigen::PlainObjectBase<DerivedV>& V,
  315. const Eigen::PlainObjectBase<DerivedF>& F,
  316. const Eigen::PlainObjectBase<DerivedP>& P,
  317. Eigen::PlainObjectBase<DerivedR>& R,
  318. Eigen::PlainObjectBase<DerivedS>& S) {
  319. const size_t num_faces = F.rows();
  320. Eigen::VectorXi I(num_faces);
  321. I.setLinSpaced(num_faces, 0, num_faces-1);
  322. igl::copyleft::cgal::closest_facet(V, F, I, P, R, S);
  323. }
  324. #ifdef IGL_STATIC_LIBRARY
  325. template void igl::copyleft::cgal::closest_facet<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -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::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<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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  326. #endif