order_facets_around_edge.cpp 16 KB

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