order_facets_around_edges.cpp 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "order_facets_around_edges.h"
  9. #include "order_facets_around_edge.h"
  10. #include "../sort_angles.h"
  11. #include <type_traits>
  12. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  13. template<
  14. typename DerivedV,
  15. typename DerivedF,
  16. typename DerivedN,
  17. typename DeriveduE,
  18. typename uE2EType,
  19. typename uE2oEType,
  20. typename uE2CType >
  21. IGL_INLINE
  22. typename std::enable_if<!std::is_same<typename DerivedV::Scalar,
  23. typename CGAL::Exact_predicates_exact_constructions_kernel::FT>::value, void>::type
  24. igl::cgal::order_facets_around_edges(
  25. const Eigen::PlainObjectBase<DerivedV>& V,
  26. const Eigen::PlainObjectBase<DerivedF>& F,
  27. const Eigen::PlainObjectBase<DerivedN>& N,
  28. const Eigen::PlainObjectBase<DeriveduE>& uE,
  29. const std::vector<std::vector<uE2EType> >& uE2E,
  30. std::vector<std::vector<uE2oEType> >& uE2oE,
  31. std::vector<std::vector<uE2CType > >& uE2C ) {
  32. typedef Eigen::Matrix<typename DerivedN::Scalar, 3, 1> Vector3F;
  33. const typename DerivedV::Scalar EPS = 1e-12;
  34. const size_t num_faces = F.rows();
  35. const size_t num_undirected_edges = uE.rows();
  36. auto edge_index_to_face_index = [&](size_t ei) { return ei % num_faces; };
  37. auto edge_index_to_corner_index = [&](size_t ei) { return ei / num_faces; };
  38. uE2oE.resize(num_undirected_edges);
  39. uE2C.resize(num_undirected_edges);
  40. for(size_t ui = 0;ui<num_undirected_edges;ui++)
  41. {
  42. const auto& adj_edges = uE2E[ui];
  43. const size_t edge_valance = adj_edges.size();
  44. assert(edge_valance > 0);
  45. const auto ref_edge = adj_edges[0];
  46. const auto ref_face = edge_index_to_face_index(ref_edge);
  47. Vector3F ref_normal = N.row(ref_face);
  48. const auto ref_corner_o = edge_index_to_corner_index(ref_edge);
  49. const auto ref_corner_s = (ref_corner_o+1)%3;
  50. const auto ref_corner_d = (ref_corner_o+2)%3;
  51. const typename DerivedF::Scalar o = F(ref_face, ref_corner_o);
  52. const typename DerivedF::Scalar s = F(ref_face, ref_corner_s);
  53. const typename DerivedF::Scalar d = F(ref_face, ref_corner_d);
  54. Vector3F edge = V.row(d) - V.row(s);
  55. auto edge_len = edge.norm();
  56. bool degenerated = edge_len < EPS;
  57. if (degenerated) {
  58. if (edge_valance <= 2) {
  59. // There is only one way to order 2 or less faces.
  60. edge.setZero();
  61. } else {
  62. edge.setZero();
  63. Eigen::Matrix<typename DerivedN::Scalar, Eigen::Dynamic, 3>
  64. normals(edge_valance, 3);
  65. for (size_t fei=0; fei<edge_valance; fei++) {
  66. const auto fe = adj_edges[fei];
  67. const auto f = edge_index_to_face_index(fe);
  68. normals.row(fei) = N.row(f);
  69. }
  70. for (size_t i=0; i<edge_valance; i++) {
  71. size_t j = (i+1) % edge_valance;
  72. Vector3F ni = normals.row(i);
  73. Vector3F nj = normals.row(j);
  74. edge = ni.cross(nj);
  75. edge_len = edge.norm();
  76. if (edge_len >= EPS) {
  77. edge.normalize();
  78. break;
  79. }
  80. }
  81. // Ensure edge direction are consistent with reference face.
  82. Vector3F in_face_vec = V.row(o) - V.row(s);
  83. if (edge.cross(in_face_vec).dot(ref_normal) < 0) {
  84. edge *= -1;
  85. }
  86. if (edge.norm() < EPS) {
  87. std::cerr << "=====================================" << std::endl;
  88. std::cerr << " ui: " << ui << std::endl;
  89. std::cerr << "edge: " << ref_edge << std::endl;
  90. std::cerr << "face: " << ref_face << std::endl;
  91. std::cerr << " vs: " << V.row(s) << std::endl;
  92. std::cerr << " vd: " << V.row(d) << std::endl;
  93. std::cerr << "adj face normals: " << std::endl;
  94. std::cerr << normals << std::endl;
  95. std::cerr << "Very degenerated case detected:" << std::endl;
  96. std::cerr << "Near zero edge surrounded by "
  97. << edge_valance << " neearly colinear faces" <<
  98. std::endl;
  99. std::cerr << "=====================================" << std::endl;
  100. }
  101. }
  102. } else {
  103. edge.normalize();
  104. }
  105. Eigen::MatrixXd angle_data(edge_valance, 3);
  106. std::vector<bool> cons(edge_valance);
  107. for (size_t fei=0; fei<edge_valance; fei++) {
  108. const auto fe = adj_edges[fei];
  109. const auto f = edge_index_to_face_index(fe);
  110. const auto c = edge_index_to_corner_index(fe);
  111. cons[fei] = (d == F(f, (c+1)%3));
  112. assert( cons[fei] || (d == F(f,(c+2)%3)));
  113. assert(!cons[fei] || (s == F(f,(c+2)%3)));
  114. assert(!cons[fei] || (d == F(f,(c+1)%3)));
  115. Vector3F n = N.row(f);
  116. angle_data(fei, 0) = ref_normal.cross(n).dot(edge);
  117. angle_data(fei, 1) = ref_normal.dot(n);
  118. if (cons[fei]) {
  119. angle_data(fei, 0) *= -1;
  120. angle_data(fei, 1) *= -1;
  121. }
  122. angle_data(fei, 0) *= -1; // Sort clockwise.
  123. angle_data(fei, 2) = (cons[fei]?1.:-1.)*(f+1);
  124. }
  125. Eigen::VectorXi order;
  126. igl::sort_angles(angle_data, order);
  127. auto& ordered_edges = uE2oE[ui];
  128. auto& consistency = uE2C[ui];
  129. ordered_edges.resize(edge_valance);
  130. consistency.resize(edge_valance);
  131. for (size_t fei=0; fei<edge_valance; fei++) {
  132. ordered_edges[fei] = adj_edges[order[fei]];
  133. consistency[fei] = cons[order[fei]];
  134. }
  135. }
  136. }
  137. template<
  138. typename DerivedV,
  139. typename DerivedF,
  140. typename DerivedN,
  141. typename DeriveduE,
  142. typename uE2EType,
  143. typename uE2oEType,
  144. typename uE2CType >
  145. IGL_INLINE
  146. typename std::enable_if<std::is_same<typename DerivedV::Scalar,
  147. typename CGAL::Exact_predicates_exact_constructions_kernel::FT>::value, void>::type
  148. igl::cgal::order_facets_around_edges(
  149. const Eigen::PlainObjectBase<DerivedV>& V,
  150. const Eigen::PlainObjectBase<DerivedF>& F,
  151. const Eigen::PlainObjectBase<DerivedN>& N,
  152. const Eigen::PlainObjectBase<DeriveduE>& uE,
  153. const std::vector<std::vector<uE2EType> >& uE2E,
  154. std::vector<std::vector<uE2oEType> >& uE2oE,
  155. std::vector<std::vector<uE2CType > >& uE2C ) {
  156. typedef Eigen::Matrix<typename DerivedN::Scalar, 3, 1> Vector3F;
  157. typedef Eigen::Matrix<typename DerivedV::Scalar, 3, 1> Vector3E;
  158. const typename DerivedV::Scalar EPS = 1e-12;
  159. const size_t num_faces = F.rows();
  160. const size_t num_undirected_edges = uE.rows();
  161. auto edge_index_to_face_index = [&](size_t ei) { return ei % num_faces; };
  162. auto edge_index_to_corner_index = [&](size_t ei) { return ei / num_faces; };
  163. uE2oE.resize(num_undirected_edges);
  164. uE2C.resize(num_undirected_edges);
  165. for(size_t ui = 0;ui<num_undirected_edges;ui++)
  166. {
  167. const auto& adj_edges = uE2E[ui];
  168. const size_t edge_valance = adj_edges.size();
  169. assert(edge_valance > 0);
  170. const auto ref_edge = adj_edges[0];
  171. const auto ref_face = edge_index_to_face_index(ref_edge);
  172. Vector3F ref_normal = N.row(ref_face);
  173. const auto ref_corner_o = edge_index_to_corner_index(ref_edge);
  174. const auto ref_corner_s = (ref_corner_o+1)%3;
  175. const auto ref_corner_d = (ref_corner_o+2)%3;
  176. const typename DerivedF::Scalar o = F(ref_face, ref_corner_o);
  177. const typename DerivedF::Scalar s = F(ref_face, ref_corner_s);
  178. const typename DerivedF::Scalar d = F(ref_face, ref_corner_d);
  179. Vector3E exact_edge = V.row(d) - V.row(s);
  180. exact_edge.array() /= exact_edge.squaredNorm();
  181. Vector3F edge(
  182. CGAL::to_double(exact_edge[0]),
  183. CGAL::to_double(exact_edge[1]),
  184. CGAL::to_double(exact_edge[2]));
  185. edge.normalize();
  186. Eigen::MatrixXd angle_data(edge_valance, 3);
  187. std::vector<bool> cons(edge_valance);
  188. for (size_t fei=0; fei<edge_valance; fei++) {
  189. const auto fe = adj_edges[fei];
  190. const auto f = edge_index_to_face_index(fe);
  191. const auto c = edge_index_to_corner_index(fe);
  192. cons[fei] = (d == F(f, (c+1)%3));
  193. assert( cons[fei] || (d == F(f,(c+2)%3)));
  194. assert(!cons[fei] || (s == F(f,(c+2)%3)));
  195. assert(!cons[fei] || (d == F(f,(c+1)%3)));
  196. Vector3F n = N.row(f);
  197. angle_data(fei, 0) = ref_normal.cross(n).dot(edge);
  198. angle_data(fei, 1) = ref_normal.dot(n);
  199. if (cons[fei]) {
  200. angle_data(fei, 0) *= -1;
  201. angle_data(fei, 1) *= -1;
  202. }
  203. angle_data(fei, 0) *= -1; // Sort clockwise.
  204. angle_data(fei, 2) = (cons[fei]?1.:-1.)*(f+1);
  205. }
  206. Eigen::VectorXi order;
  207. igl::sort_angles(angle_data, order);
  208. auto& ordered_edges = uE2oE[ui];
  209. auto& consistency = uE2C[ui];
  210. ordered_edges.resize(edge_valance);
  211. consistency.resize(edge_valance);
  212. for (size_t fei=0; fei<edge_valance; fei++) {
  213. ordered_edges[fei] = adj_edges[order[fei]];
  214. consistency[fei] = cons[order[fei]];
  215. }
  216. }
  217. }
  218. template<
  219. typename DerivedV,
  220. typename DerivedF,
  221. typename DeriveduE,
  222. typename uE2EType,
  223. typename uE2oEType,
  224. typename uE2CType >
  225. IGL_INLINE void igl::cgal::order_facets_around_edges(
  226. const Eigen::PlainObjectBase<DerivedV>& V,
  227. const Eigen::PlainObjectBase<DerivedF>& F,
  228. const Eigen::PlainObjectBase<DeriveduE>& uE,
  229. const std::vector<std::vector<uE2EType> >& uE2E,
  230. std::vector<std::vector<uE2oEType> >& uE2oE,
  231. std::vector<std::vector<uE2CType > >& uE2C ) {
  232. //typedef Eigen::Matrix<typename DerivedV::Scalar, 3, 1> Vector3E;
  233. const size_t num_faces = F.rows();
  234. const size_t num_undirected_edges = uE.rows();
  235. auto edge_index_to_face_index = [&](size_t ei) { return ei % num_faces; };
  236. auto edge_index_to_corner_index = [&](size_t ei) { return ei / num_faces; };
  237. uE2oE.resize(num_undirected_edges);
  238. uE2C.resize(num_undirected_edges);
  239. for(size_t ui = 0;ui<num_undirected_edges;ui++)
  240. {
  241. const auto& adj_edges = uE2E[ui];
  242. const size_t edge_valance = adj_edges.size();
  243. assert(edge_valance > 0);
  244. const auto ref_edge = adj_edges[0];
  245. const auto ref_face = edge_index_to_face_index(ref_edge);
  246. const auto ref_corner_o = edge_index_to_corner_index(ref_edge);
  247. const auto ref_corner_s = (ref_corner_o+1)%3;
  248. const auto ref_corner_d = (ref_corner_o+2)%3;
  249. //const typename DerivedF::Scalar o = F(ref_face, ref_corner_o);
  250. const typename DerivedF::Scalar s = F(ref_face, ref_corner_s);
  251. const typename DerivedF::Scalar d = F(ref_face, ref_corner_d);
  252. std::vector<bool> cons(edge_valance);
  253. std::vector<int> adj_faces(edge_valance);
  254. for (size_t fei=0; fei<edge_valance; fei++) {
  255. const auto fe = adj_edges[fei];
  256. const auto f = edge_index_to_face_index(fe);
  257. const auto c = edge_index_to_corner_index(fe);
  258. cons[fei] = (d == F(f, (c+1)%3));
  259. adj_faces[fei] = (f+1) * (cons[fei] ? 1:-1);
  260. assert( cons[fei] || (d == F(f,(c+2)%3)));
  261. assert(!cons[fei] || (s == F(f,(c+2)%3)));
  262. assert(!cons[fei] || (d == F(f,(c+1)%3)));
  263. }
  264. Eigen::VectorXi order;
  265. order_facets_around_edge(V, F, s, d, adj_faces, order);
  266. assert(order.size() == edge_valance);
  267. auto& ordered_edges = uE2oE[ui];
  268. auto& consistency = uE2C[ui];
  269. ordered_edges.resize(edge_valance);
  270. consistency.resize(edge_valance);
  271. for (size_t fei=0; fei<edge_valance; fei++) {
  272. ordered_edges[fei] = adj_edges[order[fei]];
  273. consistency[fei] = cons[order[fei]];
  274. }
  275. }
  276. }
  277. #ifdef IGL_STATIC_LIBRARY
  278. // Explicit template specialization
  279. template void igl::cgal::order_facets_around_edges<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, 2, 0, -1, 2>, long, long, bool>(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&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 2, 0, -1, 2> > const&, std::__1::vector<std::__1::vector<long, std::__1::allocator<long> >, std::__1::allocator<std::__1::vector<long, std::__1::allocator<long> > > > const&, std::__1::vector<std::__1::vector<long, std::__1::allocator<long> >, std::__1::allocator<std::__1::vector<long, std::__1::allocator<long> > > >&, std::__1::vector<std::__1::vector<bool, std::__1::allocator<bool> >, std::__1::allocator<std::__1::vector<bool, std::__1::allocator<bool> > > >&);
  280. template void igl::cgal::order_facets_around_edges<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 2, 0, -1, 2>, long, long, bool>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 2, 0, -1, 2> > const&, std::__1::vector<std::__1::vector<long, std::__1::allocator<long> >, std::__1::allocator<std::__1::vector<long, std::__1::allocator<long> > > > const&, std::__1::vector<std::__1::vector<long, std::__1::allocator<long> >, std::__1::allocator<std::__1::vector<long, std::__1::allocator<long> > > >&, std::__1::vector<std::__1::vector<bool, std::__1::allocator<bool> >, std::__1::allocator<std::__1::vector<bool, std::__1::allocator<bool> > > >&);
  281. template void igl::cgal::order_facets_around_edges<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 2, 0, -1, 2>, long, long, bool>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 2, 0, -1, 2> > const&, std::__1::vector<std::__1::vector<long, std::__1::allocator<long> >, std::__1::allocator<std::__1::vector<long, std::__1::allocator<long> > > > const&, std::__1::vector<std::__1::vector<long, std::__1::allocator<long> >, std::__1::allocator<std::__1::vector<long, std::__1::allocator<long> > > >&, std::__1::vector<std::__1::vector<bool, std::__1::allocator<bool> >, std::__1::allocator<std::__1::vector<bool, std::__1::allocator<bool> > > >&);
  282. #endif