order_facets_around_edge.cpp 16 KB

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