order_facets_around_edge.cpp 16 KB

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