order_facets_around_edge.cpp 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403
  1. #include "order_facets_around_edge.h"
  2. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  3. #include <stdexcept>
  4. namespace igl
  5. {
  6. namespace cgal
  7. {
  8. namespace order_facets_around_edges_helper
  9. {
  10. template<typename T>
  11. std::vector<size_t> index_sort(const std::vector<T>& data)
  12. {
  13. const size_t len = data.size();
  14. std::vector<size_t> order(len);
  15. for (size_t i=0; i<len; i++)
  16. {
  17. order[i] = i;
  18. }
  19. auto comp = [&](size_t i, size_t j)
  20. {
  21. return data[i] < data[j];
  22. };
  23. std::sort(order.begin(), order.end(), comp);
  24. return order;
  25. }
  26. }
  27. }
  28. }
  29. // adj_faces contains signed index starting from +- 1.
  30. template<
  31. typename DerivedV,
  32. typename DerivedF,
  33. typename DerivedI >
  34. void igl::cgal::order_facets_around_edge(
  35. const Eigen::PlainObjectBase<DerivedV>& V,
  36. const Eigen::PlainObjectBase<DerivedF>& F,
  37. size_t s,
  38. size_t d,
  39. const std::vector<int>& adj_faces,
  40. Eigen::PlainObjectBase<DerivedI>& order)
  41. {
  42. using namespace igl::cgal::order_facets_around_edges_helper;
  43. // Although we only need exact predicates in the algorithm,
  44. // exact constructions are needed to avoid degeneracies due to
  45. // casting to double.
  46. typedef CGAL::Exact_predicates_exact_constructions_kernel K;
  47. typedef K::Point_3 Point_3;
  48. typedef K::Plane_3 Plane_3;
  49. auto get_face_index = [&](int adj_f)->size_t
  50. {
  51. return abs(adj_f) - 1;
  52. };
  53. auto get_opposite_vertex = [&](size_t fid)->size_t
  54. {
  55. typedef typename DerivedF::Scalar Index;
  56. if (F(fid, 0) != (Index)s && F(fid, 0) != (Index)d) return F(fid, 0);
  57. if (F(fid, 1) != (Index)s && F(fid, 1) != (Index)d) return F(fid, 1);
  58. if (F(fid, 2) != (Index)s && F(fid, 2) != (Index)d) return F(fid, 2);
  59. assert(false);
  60. return -1;
  61. };
  62. // Handle base cases
  63. if (adj_faces.size() == 0)
  64. {
  65. order.resize(0, 1);
  66. return;
  67. } else if (adj_faces.size() == 1)
  68. {
  69. order.resize(1, 1);
  70. order(0, 0) = 0;
  71. return;
  72. } else if (adj_faces.size() == 2)
  73. {
  74. const size_t o1 = get_opposite_vertex(get_face_index(adj_faces[0]));
  75. const size_t o2 = get_opposite_vertex(get_face_index(adj_faces[1]));
  76. const Point_3 ps(V(s, 0), V(s, 1), V(s, 2));
  77. const Point_3 pd(V(d, 0), V(d, 1), V(d, 2));
  78. const Point_3 p1(V(o1, 0), V(o1, 1), V(o1, 2));
  79. const Point_3 p2(V(o2, 0), V(o2, 1), V(o2, 2));
  80. order.resize(2, 1);
  81. switch (CGAL::orientation(ps, pd, p1, p2))
  82. {
  83. case CGAL::POSITIVE:
  84. order(0, 0) = 1;
  85. order(1, 0) = 0;
  86. break;
  87. case CGAL::NEGATIVE:
  88. order(0, 0) = 0;
  89. order(1, 0) = 1;
  90. break;
  91. case CGAL::COPLANAR:
  92. {
  93. switch (CGAL::coplanar_orientation(ps, pd, p1, p2)) {
  94. case CGAL::POSITIVE:
  95. // Duplicated face, use index to break tie.
  96. order(0, 0) = adj_faces[0] < adj_faces[1] ? 0:1;
  97. order(1, 0) = adj_faces[0] < adj_faces[1] ? 1:0;
  98. break;
  99. case CGAL::NEGATIVE:
  100. // Coplanar faces, one on each side of the edge.
  101. // It is equally valid to order them (0, 1) or (1, 0).
  102. // I cannot think of any reason to prefer one to the
  103. // other. So just use (0, 1) ordering by default.
  104. order(0, 0) = 0;
  105. order(1, 0) = 1;
  106. break;
  107. case CGAL::COLLINEAR:
  108. std::cerr << "Degenerated triangle detected." <<
  109. std::endl;
  110. assert(false);
  111. break;
  112. default:
  113. assert(false);
  114. }
  115. }
  116. break;
  117. default:
  118. assert(false);
  119. }
  120. return;
  121. }
  122. const size_t num_adj_faces = adj_faces.size();
  123. const size_t o = get_opposite_vertex( get_face_index(adj_faces[0]));
  124. const Point_3 p_s(V(s, 0), V(s, 1), V(s, 2));
  125. const Point_3 p_d(V(d, 0), V(d, 1), V(d, 2));
  126. const Point_3 p_o(V(o, 0), V(o, 1), V(o, 2));
  127. const Plane_3 separator(p_s, p_d, p_o);
  128. if (separator.is_degenerate()) {
  129. throw std::runtime_error(
  130. "Cannot order facets around edge due to degenerated facets");
  131. }
  132. std::vector<Point_3> opposite_vertices;
  133. for (size_t i=0; i<num_adj_faces; i++)
  134. {
  135. const size_t o = get_opposite_vertex( get_face_index(adj_faces[i]));
  136. opposite_vertices.emplace_back(
  137. V(o, 0), V(o, 1), V(o, 2));
  138. }
  139. std::vector<int> positive_side;
  140. std::vector<int> negative_side;
  141. std::vector<int> tie_positive_oriented;
  142. std::vector<int> tie_negative_oriented;
  143. std::vector<size_t> positive_side_index;
  144. std::vector<size_t> negative_side_index;
  145. std::vector<size_t> tie_positive_oriented_index;
  146. std::vector<size_t> tie_negative_oriented_index;
  147. for (size_t i=0; i<num_adj_faces; i++)
  148. {
  149. const int f = adj_faces[i];
  150. const Point_3& p_a = opposite_vertices[i];
  151. auto orientation = separator.oriented_side(p_a);
  152. switch (orientation) {
  153. case CGAL::ON_POSITIVE_SIDE:
  154. positive_side.push_back(f);
  155. positive_side_index.push_back(i);
  156. break;
  157. case CGAL::ON_NEGATIVE_SIDE:
  158. negative_side.push_back(f);
  159. negative_side_index.push_back(i);
  160. break;
  161. case CGAL::ON_ORIENTED_BOUNDARY:
  162. {
  163. auto inplane_orientation = CGAL::coplanar_orientation(
  164. p_s, p_d, p_o, p_a);
  165. switch (inplane_orientation) {
  166. case CGAL::POSITIVE:
  167. tie_positive_oriented.push_back(f);
  168. tie_positive_oriented_index.push_back(i);
  169. break;
  170. case CGAL::NEGATIVE:
  171. tie_negative_oriented.push_back(f);
  172. tie_negative_oriented_index.push_back(i);
  173. break;
  174. case CGAL::COLLINEAR:
  175. default:
  176. throw std::runtime_error(
  177. "Degenerated facet detected.");
  178. break;
  179. }
  180. }
  181. break;
  182. default:
  183. // Should not be here.
  184. throw std::runtime_error("Unknown CGAL state detected.");
  185. }
  186. }
  187. Eigen::PlainObjectBase<DerivedI> positive_order, negative_order;
  188. order_facets_around_edge(V, F, s, d, positive_side, positive_order);
  189. order_facets_around_edge(V, F, s, d, negative_side, negative_order);
  190. std::vector<size_t> tie_positive_order = index_sort(tie_positive_oriented);
  191. std::vector<size_t> tie_negative_order = index_sort(tie_negative_oriented);
  192. // Copy results into order vector.
  193. const size_t tie_positive_size = tie_positive_oriented.size();
  194. const size_t tie_negative_size = tie_negative_oriented.size();
  195. const size_t positive_size = positive_order.size();
  196. const size_t negative_size = negative_order.size();
  197. order.resize(
  198. tie_positive_size + positive_size + tie_negative_size + negative_size,1);
  199. size_t count=0;
  200. for (size_t i=0; i<tie_positive_size; i++)
  201. {
  202. order(count+i, 0) = tie_positive_oriented_index[tie_positive_order[i]];
  203. }
  204. count += tie_positive_size;
  205. for (size_t i=0; i<negative_size; i++)
  206. {
  207. order(count+i, 0) = negative_side_index[negative_order(i, 0)];
  208. }
  209. count += negative_size;
  210. for (size_t i=0; i<tie_negative_size; i++)
  211. {
  212. order(count+i, 0) = tie_negative_oriented_index[tie_negative_order[i]];
  213. }
  214. count += tie_negative_size;
  215. for (size_t i=0; i<positive_size; i++)
  216. {
  217. order(count+i, 0) = positive_side_index[positive_order(i, 0)];
  218. }
  219. count += positive_size;
  220. assert(count == num_adj_faces);
  221. // Find the correct start point.
  222. size_t start_idx = 0;
  223. for (size_t i=0; i<num_adj_faces; i++)
  224. {
  225. const Point_3& p_a = opposite_vertices[order(i, 0)];
  226. const Point_3& p_b =
  227. opposite_vertices[order((i+1)%num_adj_faces, 0)];
  228. auto orientation = CGAL::orientation(p_s, p_d, p_a, p_b);
  229. if (orientation == CGAL::POSITIVE)
  230. {
  231. // Angle between triangle (p_s, p_d, p_a) and (p_s, p_d, p_b) is
  232. // more than 180 degrees.
  233. start_idx = (i+1)%num_adj_faces;
  234. break;
  235. } else if (orientation == CGAL::COPLANAR &&
  236. Plane_3(p_s, p_d, p_a).orthogonal_direction() !=
  237. Plane_3(p_s, p_d, p_b).orthogonal_direction())
  238. {
  239. // All 4 points are coplanar, but p_a and p_b are on each side of
  240. // the edge (p_s, p_d). This means the angle between triangle
  241. // (p_s, p_d, p_a) and (p_s, p_d, p_b) is exactly 180 degrees.
  242. start_idx = (i+1)%num_adj_faces;
  243. break;
  244. }
  245. }
  246. DerivedI circular_order = order;
  247. for (size_t i=0; i<num_adj_faces; i++)
  248. {
  249. order(i, 0) = circular_order((start_idx + i)%num_adj_faces, 0);
  250. }
  251. }
  252. template<
  253. typename DerivedV,
  254. typename DerivedF,
  255. typename DerivedI>
  256. IGL_INLINE
  257. void igl::cgal::order_facets_around_edge(
  258. const Eigen::PlainObjectBase<DerivedV>& V,
  259. const Eigen::PlainObjectBase<DerivedF>& F,
  260. size_t s,
  261. size_t d,
  262. const std::vector<int>& adj_faces,
  263. const Eigen::PlainObjectBase<DerivedV>& pivot_point,
  264. Eigen::PlainObjectBase<DerivedI>& order)
  265. {
  266. assert(V.cols() == 3);
  267. assert(F.cols() == 3);
  268. assert(pivot_point.cols() == 3);
  269. auto signed_index_to_index = [&](int signed_idx)
  270. {
  271. return abs(signed_idx) -1;
  272. };
  273. auto get_opposite_vertex_index = [&](size_t fid)
  274. {
  275. typedef typename DerivedF::Scalar Index;
  276. if (F(fid, 0) != (Index)s && F(fid, 0) != (Index)d) return F(fid, 0);
  277. if (F(fid, 1) != (Index)s && F(fid, 1) != (Index)d) return F(fid, 1);
  278. if (F(fid, 2) != (Index)s && F(fid, 2) != (Index)d) return F(fid, 2);
  279. assert(false);
  280. // avoid warning
  281. return -1;
  282. };
  283. {
  284. // Check if s, d and pivot are collinear.
  285. typedef CGAL::Exact_predicates_exact_constructions_kernel K;
  286. K::Point_3 ps(V(s,0), V(s,1), V(s,2));
  287. K::Point_3 pd(V(d,0), V(d,1), V(d,2));
  288. K::Point_3 pp(pivot_point(0,0), pivot_point(0,1), pivot_point(0,2));
  289. if (CGAL::collinear(ps, pd, pp)) {
  290. throw std::runtime_error(
  291. "Pivot point is collinear with the outer edge!");
  292. }
  293. }
  294. const size_t N = adj_faces.size();
  295. const size_t num_faces = N + 1; // N adj faces + 1 pivot face
  296. // Because face indices are used for tie breaking, the original face indices
  297. // in the new faces array must be ascending.
  298. auto comp = [&](int i, int j) {
  299. return signed_index_to_index(adj_faces[i]) <
  300. signed_index_to_index(adj_faces[j]);
  301. };
  302. std::vector<size_t> adj_order(N);
  303. for (size_t i=0; i<N; i++) adj_order[i] = i;
  304. std::sort(adj_order.begin(), adj_order.end(), comp);
  305. DerivedV vertices(num_faces + 2, 3);
  306. for (size_t i=0; i<N; i++) {
  307. const size_t fid = signed_index_to_index(adj_faces[adj_order[i]]);
  308. vertices.row(i) = V.row(get_opposite_vertex_index(fid));
  309. }
  310. vertices.row(N ) = pivot_point;
  311. vertices.row(N+1) = V.row(s);
  312. vertices.row(N+2) = V.row(d);
  313. DerivedF faces(num_faces, 3);
  314. for (size_t i=0; i<N; i++)
  315. {
  316. if (adj_faces[adj_order[i]] < 0) {
  317. faces(i,0) = N+1; // s
  318. faces(i,1) = N+2; // d
  319. faces(i,2) = i ;
  320. } else {
  321. faces(i,0) = N+2; // d
  322. faces(i,1) = N+1; // s
  323. faces(i,2) = i ;
  324. }
  325. }
  326. // Last face is the pivot face.
  327. faces(N, 0) = N+1;
  328. faces(N, 1) = N+2;
  329. faces(N, 2) = N;
  330. std::vector<int> adj_faces_with_pivot(num_faces);
  331. for (size_t i=0; i<num_faces; i++)
  332. {
  333. if (faces(i,0) == N+1 && faces(i,1) == N+2)
  334. {
  335. adj_faces_with_pivot[i] = int(i+1) * -1;
  336. } else
  337. {
  338. adj_faces_with_pivot[i] = int(i+1);
  339. }
  340. }
  341. DerivedI order_with_pivot;
  342. igl::cgal::order_facets_around_edge(
  343. vertices, faces, N+1, N+2,
  344. adj_faces_with_pivot, order_with_pivot);
  345. assert(order_with_pivot.size() == num_faces);
  346. order.resize(N);
  347. size_t pivot_index = num_faces + 1;
  348. for (size_t i=0; i<num_faces; i++)
  349. {
  350. if (order_with_pivot[i] == N)
  351. {
  352. pivot_index = i;
  353. break;
  354. }
  355. }
  356. assert(pivot_index < num_faces);
  357. for (size_t i=0; i<N; i++)
  358. {
  359. order[i] = adj_order[order_with_pivot[(pivot_index+i+1)%num_faces]];
  360. }
  361. }
  362. #ifdef IGL_STATIC_LIBRARY
  363. // Explicit template specialization
  364. template void igl::cgal::order_facets_around_edge<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  365. template void igl::cgal::order_facets_around_edge<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  366. template void igl::cgal::order_facets_around_edge<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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  367. template void igl::cgal::order_facets_around_edge<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > 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> >&);
  368. template void igl::cgal::order_facets_around_edge<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  369. template void igl::cgal::order_facets_around_edge<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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  370. #endif