order_facets_around_edge.cpp 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426
  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, bool debug)
  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. if (debug) {
  188. std::cout << "tie positive: " << std::endl;
  189. for (auto& f : tie_positive_oriented) {
  190. std::cout << get_face_index(f) << " ";
  191. }
  192. std::cout << std::endl;
  193. std::cout << "positive side: " << std::endl;
  194. for (auto& f : positive_side) {
  195. std::cout << get_face_index(f) << " ";
  196. }
  197. std::cout << std::endl;
  198. std::cout << "tie negative: " << std::endl;
  199. for (auto& f : tie_negative_oriented) {
  200. std::cout << get_face_index(f) << " ";
  201. }
  202. std::cout << std::endl;
  203. std::cout << "negative side: " << std::endl;
  204. for (auto& f : negative_side) {
  205. std::cout << get_face_index(f) << " ";
  206. }
  207. std::cout << std::endl;
  208. }
  209. Eigen::PlainObjectBase<DerivedI> positive_order, negative_order;
  210. order_facets_around_edge(V, F, s, d, positive_side, positive_order, debug);
  211. order_facets_around_edge(V, F, s, d, negative_side, negative_order, debug);
  212. std::vector<size_t> tie_positive_order = index_sort(tie_positive_oriented);
  213. std::vector<size_t> tie_negative_order = index_sort(tie_negative_oriented);
  214. // Copy results into order vector.
  215. const size_t tie_positive_size = tie_positive_oriented.size();
  216. const size_t tie_negative_size = tie_negative_oriented.size();
  217. const size_t positive_size = positive_order.size();
  218. const size_t negative_size = negative_order.size();
  219. order.resize(
  220. tie_positive_size + positive_size + tie_negative_size + negative_size,1);
  221. size_t count=0;
  222. for (size_t i=0; i<tie_positive_size; i++)
  223. {
  224. order(count+i, 0) = tie_positive_oriented_index[tie_positive_order[i]];
  225. }
  226. count += tie_positive_size;
  227. for (size_t i=0; i<negative_size; i++)
  228. {
  229. order(count+i, 0) = negative_side_index[negative_order(i, 0)];
  230. }
  231. count += negative_size;
  232. for (size_t i=0; i<tie_negative_size; i++)
  233. {
  234. order(count+i, 0) = tie_negative_oriented_index[tie_negative_order[i]];
  235. }
  236. count += tie_negative_size;
  237. for (size_t i=0; i<positive_size; i++)
  238. {
  239. order(count+i, 0) = positive_side_index[positive_order(i, 0)];
  240. }
  241. count += positive_size;
  242. assert(count == num_adj_faces);
  243. // Find the correct start point.
  244. size_t start_idx = 0;
  245. for (size_t i=0; i<num_adj_faces; i++)
  246. {
  247. const Point_3& p_a = opposite_vertices[order(i, 0)];
  248. const Point_3& p_b =
  249. opposite_vertices[order((i+1)%num_adj_faces, 0)];
  250. auto orientation = CGAL::orientation(p_s, p_d, p_a, p_b);
  251. if (orientation == CGAL::POSITIVE)
  252. {
  253. // Angle between triangle (p_s, p_d, p_a) and (p_s, p_d, p_b) is
  254. // more than 180 degrees.
  255. start_idx = (i+1)%num_adj_faces;
  256. break;
  257. } else if (orientation == CGAL::COPLANAR &&
  258. Plane_3(p_s, p_d, p_a).orthogonal_direction() !=
  259. Plane_3(p_s, p_d, p_b).orthogonal_direction())
  260. {
  261. // All 4 points are coplanar, but p_a and p_b are on each side of
  262. // the edge (p_s, p_d). This means the angle between triangle
  263. // (p_s, p_d, p_a) and (p_s, p_d, p_b) is exactly 180 degrees.
  264. start_idx = (i+1)%num_adj_faces;
  265. break;
  266. }
  267. }
  268. DerivedI circular_order = order;
  269. for (size_t i=0; i<num_adj_faces; i++)
  270. {
  271. order(i, 0) = circular_order((start_idx + i)%num_adj_faces, 0);
  272. }
  273. }
  274. template<
  275. typename DerivedV,
  276. typename DerivedF,
  277. typename DerivedI>
  278. IGL_INLINE
  279. void igl::cgal::order_facets_around_edge(
  280. const Eigen::PlainObjectBase<DerivedV>& V,
  281. const Eigen::PlainObjectBase<DerivedF>& F,
  282. size_t s,
  283. size_t d,
  284. const std::vector<int>& adj_faces,
  285. const Eigen::PlainObjectBase<DerivedV>& pivot_point,
  286. Eigen::PlainObjectBase<DerivedI>& order)
  287. {
  288. assert(V.cols() == 3);
  289. assert(F.cols() == 3);
  290. assert(pivot_point.cols() == 3);
  291. auto signed_index_to_index = [&](int signed_idx)
  292. {
  293. return abs(signed_idx) -1;
  294. };
  295. auto get_opposite_vertex_index = [&](size_t fid)
  296. {
  297. typedef typename DerivedF::Scalar Index;
  298. if (F(fid, 0) != (Index)s && F(fid, 0) != (Index)d) return F(fid, 0);
  299. if (F(fid, 1) != (Index)s && F(fid, 1) != (Index)d) return F(fid, 1);
  300. if (F(fid, 2) != (Index)s && F(fid, 2) != (Index)d) return F(fid, 2);
  301. assert(false);
  302. // avoid warning
  303. return -1;
  304. };
  305. {
  306. // Check if s, d and pivot are collinear.
  307. typedef CGAL::Exact_predicates_exact_constructions_kernel K;
  308. K::Point_3 ps(V(s,0), V(s,1), V(s,2));
  309. K::Point_3 pd(V(d,0), V(d,1), V(d,2));
  310. K::Point_3 pp(pivot_point(0,0), pivot_point(0,1), pivot_point(0,2));
  311. if (CGAL::collinear(ps, pd, pp)) {
  312. throw std::runtime_error(
  313. "Pivot point is collinear with the outer edge!");
  314. }
  315. }
  316. const size_t N = adj_faces.size();
  317. const size_t num_faces = N + 1; // N adj faces + 1 pivot face
  318. // Because face indices are used for tie breaking, the original face indices
  319. // in the new faces array must be ascending.
  320. auto comp = [&](int i, int j) {
  321. return signed_index_to_index(adj_faces[i]) <
  322. signed_index_to_index(adj_faces[j]);
  323. };
  324. std::vector<size_t> adj_order(N);
  325. for (size_t i=0; i<N; i++) adj_order[i] = i;
  326. std::sort(adj_order.begin(), adj_order.end(), comp);
  327. DerivedV vertices(num_faces + 2, 3);
  328. for (size_t i=0; i<N; i++) {
  329. const size_t fid = signed_index_to_index(adj_faces[adj_order[i]]);
  330. vertices.row(i) = V.row(get_opposite_vertex_index(fid));
  331. }
  332. vertices.row(N ) = pivot_point;
  333. vertices.row(N+1) = V.row(s);
  334. vertices.row(N+2) = V.row(d);
  335. DerivedF faces(num_faces, 3);
  336. for (size_t i=0; i<N; i++)
  337. {
  338. if (adj_faces[adj_order[i]] < 0) {
  339. faces(i,0) = N+1; // s
  340. faces(i,1) = N+2; // d
  341. faces(i,2) = i ;
  342. } else {
  343. faces(i,0) = N+2; // d
  344. faces(i,1) = N+1; // s
  345. faces(i,2) = i ;
  346. }
  347. }
  348. // Last face is the pivot face.
  349. faces(N, 0) = N+1;
  350. faces(N, 1) = N+2;
  351. faces(N, 2) = N;
  352. std::vector<int> adj_faces_with_pivot(num_faces);
  353. for (size_t i=0; i<num_faces; i++)
  354. {
  355. if ((size_t)faces(i,0) == N+1 && (size_t)faces(i,1) == N+2)
  356. {
  357. adj_faces_with_pivot[i] = int(i+1) * -1;
  358. } else
  359. {
  360. adj_faces_with_pivot[i] = int(i+1);
  361. }
  362. }
  363. DerivedI order_with_pivot;
  364. igl::cgal::order_facets_around_edge(
  365. vertices, faces, N+1, N+2,
  366. adj_faces_with_pivot, order_with_pivot);
  367. assert((size_t)order_with_pivot.size() == num_faces);
  368. order.resize(N);
  369. size_t pivot_index = num_faces + 1;
  370. for (size_t i=0; i<num_faces; i++)
  371. {
  372. if ((size_t)order_with_pivot[i] == N)
  373. {
  374. pivot_index = i;
  375. break;
  376. }
  377. }
  378. assert(pivot_index < num_faces);
  379. for (size_t i=0; i<N; i++)
  380. {
  381. order[i] = adj_order[order_with_pivot[(pivot_index+i+1)%num_faces]];
  382. }
  383. }
  384. #ifdef IGL_STATIC_LIBRARY
  385. // Explicit template specialization
  386. template void
  387. 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::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  388. 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::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  389. 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::vector<int, std::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> >&);
  390. 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::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  391. 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::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  392. template void igl::cgal::order_facets_around_edge<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&, unsigned long, unsigned long, std::__1::vector<int, std::__1::allocator<int> > 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> >&);
  393. #endif