order_facets_around_edge.cpp 17 KB

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