outer_hull_legacy.cpp 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <[email protected]>
  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 "outer_hull_legacy.h"
  9. #include "extract_cells.h"
  10. #include "remesh_self_intersections.h"
  11. #include "assign.h"
  12. #include "../../remove_unreferenced.h"
  13. #include <CGAL/AABB_tree.h>
  14. #include <CGAL/AABB_traits.h>
  15. #include <CGAL/AABB_triangle_primitive.h>
  16. #include <CGAL/intersections.h>
  17. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  18. #include "points_inside_component.h"
  19. #include "order_facets_around_edges.h"
  20. #include "outer_facet.h"
  21. #include "../../sortrows.h"
  22. #include "../../facet_components.h"
  23. #include "../../winding_number.h"
  24. #include "../../triangle_triangle_adjacency.h"
  25. #include "../../unique_edge_map.h"
  26. #include "../../barycenter.h"
  27. #include "../../per_face_normals.h"
  28. #include "../../PlainMatrix.h"
  29. #include "../../sort_angles.h"
  30. #include <Eigen/Geometry>
  31. #include <vector>
  32. #include <map>
  33. #include <queue>
  34. #include <iostream>
  35. #include <type_traits>
  36. #include <CGAL/number_utils.h>
  37. //#define IGL_OUTER_HULL_DEBUG
  38. template <
  39. typename DerivedV,
  40. typename DerivedF,
  41. typename DerivedG,
  42. typename DerivedJ,
  43. typename Derivedflip>
  44. IGL_INLINE void igl::copyleft::cgal::outer_hull_legacy(
  45. const Eigen::MatrixBase<DerivedV> & V,
  46. const Eigen::MatrixBase<DerivedF> & F,
  47. Eigen::PlainObjectBase<DerivedG> & G,
  48. Eigen::PlainObjectBase<DerivedJ> & J,
  49. Eigen::PlainObjectBase<Derivedflip> & flip)
  50. {
  51. #ifdef IGL_OUTER_HULL_DEBUG
  52. std::cerr << "Extracting outer hull" << std::endl;
  53. #endif
  54. using namespace Eigen;
  55. using namespace std;
  56. typedef typename DerivedF::Index Index;
  57. Matrix<Index,DerivedF::RowsAtCompileTime,1> C;
  58. typedef Matrix<typename DerivedV::Scalar,Dynamic,DerivedV::ColsAtCompileTime> MatrixXV;
  59. //typedef Matrix<typename DerivedF::Scalar,Dynamic,DerivedF::ColsAtCompileTime> MatrixXF;
  60. typedef Matrix<typename DerivedG::Scalar,Dynamic,DerivedG::ColsAtCompileTime> MatrixXG;
  61. typedef Matrix<typename DerivedJ::Scalar,Dynamic,DerivedJ::ColsAtCompileTime> MatrixXJ;
  62. const Index m = F.rows();
  63. // UNUSED:
  64. //const auto & duplicate_simplex = [&F](const int f, const int g)->bool
  65. //{
  66. // return
  67. // (F(f,0) == F(g,0) && F(f,1) == F(g,1) && F(f,2) == F(g,2)) ||
  68. // (F(f,1) == F(g,0) && F(f,2) == F(g,1) && F(f,0) == F(g,2)) ||
  69. // (F(f,2) == F(g,0) && F(f,0) == F(g,1) && F(f,1) == F(g,2)) ||
  70. // (F(f,0) == F(g,2) && F(f,1) == F(g,1) && F(f,2) == F(g,0)) ||
  71. // (F(f,1) == F(g,2) && F(f,2) == F(g,1) && F(f,0) == F(g,0)) ||
  72. // (F(f,2) == F(g,2) && F(f,0) == F(g,1) && F(f,1) == F(g,0));
  73. //};
  74. #ifdef IGL_OUTER_HULL_DEBUG
  75. cout<<"outer hull..."<<endl;
  76. #endif
  77. #ifdef IGL_OUTER_HULL_DEBUG
  78. cout<<"edge map..."<<endl;
  79. #endif
  80. typedef Matrix<typename DerivedF::Scalar,Dynamic,2> MatrixX2I;
  81. typedef Matrix<typename DerivedF::Index,Dynamic,1> VectorXI;
  82. //typedef Matrix<typename DerivedV::Scalar, 3, 1> Vector3F;
  83. MatrixX2I E,uE;
  84. VectorXI EMAP;
  85. vector<vector<typename DerivedF::Index> > uE2E;
  86. unique_edge_map(F,E,uE,EMAP,uE2E);
  87. #ifdef IGL_OUTER_HULL_DEBUG
  88. for (size_t ui=0; ui<uE.rows(); ui++) {
  89. std::cout << ui << ": " << uE2E[ui].size() << " -- (";
  90. for (size_t i=0; i<uE2E[ui].size(); i++) {
  91. std::cout << uE2E[ui][i] << ", ";
  92. }
  93. std::cout << ")" << std::endl;
  94. }
  95. #endif
  96. std::vector<std::vector<typename DerivedF::Index> > uE2oE;
  97. std::vector<std::vector<bool> > uE2C;
  98. order_facets_around_edges(V, F, uE, uE2E, uE2oE, uE2C);
  99. uE2E = uE2oE;
  100. VectorXI diIM(3*m);
  101. for (auto ue : uE2E) {
  102. for (size_t i=0; i<ue.size(); i++) {
  103. auto fe = ue[i];
  104. diIM[fe] = i;
  105. }
  106. }
  107. vector<vector<vector<Index > > > TT,_1;
  108. triangle_triangle_adjacency(E,EMAP,uE2E,false,TT,_1);
  109. VectorXI counts;
  110. #ifdef IGL_OUTER_HULL_DEBUG
  111. cout<<"facet components..."<<endl;
  112. #endif
  113. facet_components(TT,C,counts);
  114. assert(C.maxCoeff()+1 == counts.rows());
  115. const size_t ncc = counts.rows();
  116. G.resize(0,F.cols());
  117. J.resize(0,1);
  118. flip.setConstant(m,1,false);
  119. #ifdef IGL_OUTER_HULL_DEBUG
  120. cout<<"reindex..."<<endl;
  121. #endif
  122. // H contains list of faces on outer hull;
  123. vector<bool> FH(m,false);
  124. vector<bool> EH(3*m,false);
  125. vector<MatrixXG> vG(ncc);
  126. vector<MatrixXJ> vJ(ncc);
  127. vector<MatrixXJ> vIM(ncc);
  128. //size_t face_count = 0;
  129. for(size_t id = 0;id<ncc;id++)
  130. {
  131. vIM[id].resize(counts[id],1);
  132. }
  133. // current index into each IM
  134. vector<size_t> g(ncc,0);
  135. // place order of each face in its respective component
  136. for(Index f = 0;f<m;f++)
  137. {
  138. vIM[C(f)](g[C(f)]++) = f;
  139. }
  140. #ifdef IGL_OUTER_HULL_DEBUG
  141. cout<<"barycenters..."<<endl;
  142. #endif
  143. // assumes that "resolve" has handled any coplanar cases correctly and nearly
  144. // coplanar cases can be sorted based on barycenter.
  145. MatrixXV BC;
  146. barycenter(V,F,BC);
  147. #ifdef IGL_OUTER_HULL_DEBUG
  148. cout<<"loop over CCs (="<<ncc<<")..."<<endl;
  149. #endif
  150. for(Index id = 0;id<(Index)ncc;id++)
  151. {
  152. auto & IM = vIM[id];
  153. // starting face that's guaranteed to be on the outer hull and in this
  154. // component
  155. int f;
  156. bool f_flip;
  157. #ifdef IGL_OUTER_HULL_DEBUG
  158. cout<<"outer facet..."<<endl;
  159. #endif
  160. igl::copyleft::cgal::outer_facet(V,F,IM,f,f_flip);
  161. #ifdef IGL_OUTER_HULL_DEBUG
  162. cout<<"outer facet: "<<f<<endl;
  163. //cout << V.row(F(f, 0)) << std::endl;
  164. //cout << V.row(F(f, 1)) << std::endl;
  165. //cout << V.row(F(f, 2)) << std::endl;
  166. #endif
  167. int FHcount = 1;
  168. FH[f] = true;
  169. // Q contains list of face edges to continue traversing upong
  170. queue<int> Q;
  171. Q.push(f+0*m);
  172. Q.push(f+1*m);
  173. Q.push(f+2*m);
  174. flip(f) = f_flip;
  175. //std::cout << "face " << face_count++ << ": " << f << std::endl;
  176. //std::cout << "f " << F.row(f).array()+1 << std::endl;
  177. //cout<<"flip("<<f<<") = "<<(flip(f)?"true":"false")<<endl;
  178. #ifdef IGL_OUTER_HULL_DEBUG
  179. cout<<"BFS..."<<endl;
  180. #endif
  181. while(!Q.empty())
  182. {
  183. // face-edge
  184. const int e = Q.front();
  185. Q.pop();
  186. // face
  187. const int f = e%m;
  188. // corner
  189. const int c = e/m;
  190. #ifdef IGL_OUTER_HULL_DEBUG
  191. std::cout << "edge: " << e << ", ue: " << EMAP(e) << std::endl;
  192. std::cout << "face: " << f << std::endl;
  193. std::cout << "corner: " << c << std::endl;
  194. std::cout << "consistent: " << uE2C[EMAP(e)][diIM[e]] << std::endl;
  195. #endif
  196. // Should never see edge again...
  197. if(EH[e] == true)
  198. {
  199. continue;
  200. }
  201. EH[e] = true;
  202. // source of edge according to f
  203. const int fs = flip(f)?F(f,(c+2)%3):F(f,(c+1)%3);
  204. // destination of edge according to f
  205. const int fd = flip(f)?F(f,(c+1)%3):F(f,(c+2)%3);
  206. // edge valence
  207. const size_t val = uE2E[EMAP(e)].size();
  208. #ifdef IGL_OUTER_HULL_DEBUG
  209. //std::cout << "vd: " << V.row(fd) << std::endl;
  210. //std::cout << "vs: " << V.row(fs) << std::endl;
  211. //std::cout << "edge: " << V.row(fd) - V.row(fs) << std::endl;
  212. for (size_t i=0; i<val; i++) {
  213. if (i == diIM(e)) {
  214. std::cout << "* ";
  215. } else {
  216. std::cout << " ";
  217. }
  218. std::cout << i << ": "
  219. << " (e: " << uE2E[EMAP(e)][i] << ", f: "
  220. << uE2E[EMAP(e)][i] % m * (uE2C[EMAP(e)][i] ? 1:-1) << ")" << std::endl;
  221. }
  222. #endif
  223. // is edge consistent with edge of face used for sorting
  224. const int e_cons = (uE2C[EMAP(e)][diIM(e)] ? 1: -1);
  225. int nfei = -1;
  226. // Loop once around trying to find suitable next face
  227. for(size_t step = 1; step<val+2;step++)
  228. {
  229. const int nfei_new = (diIM(e) + 2*val + e_cons*step*(flip(f)?-1:1))%val;
  230. const int nf = uE2E[EMAP(e)][nfei_new] % m;
  231. {
  232. #ifdef IGL_OUTER_HULL_DEBUG
  233. //cout<<"Next facet: "<<(f+1)<<" --> "<<(nf+1)<<", |"<<
  234. // di[EMAP(e)][diIM(e)]<<" - "<<di[EMAP(e)][nfei_new]<<"| = "<<
  235. // abs(di[EMAP(e)][diIM(e)] - di[EMAP(e)][nfei_new])
  236. // <<endl;
  237. #endif
  238. // Only use this face if not already seen
  239. if(!FH[nf])
  240. {
  241. nfei = nfei_new;
  242. //} else {
  243. // std::cout << "skipping face " << nfei_new << " because it is seen before"
  244. // << std::endl;
  245. }
  246. break;
  247. //} else {
  248. // std::cout << di[EMAP(e)][diIM(e)].transpose() << std::endl;
  249. // std::cout << di[EMAP(e)][diIM(nfei_new)].transpose() << std::endl;
  250. // std::cout << "skipping face " << nfei_new << " with identical dihedral angle"
  251. // << std::endl;
  252. }
  253. //#ifdef IGL_OUTER_HULL_DEBUG
  254. // cout<<"Skipping co-planar facet: "<<(f+1)<<" --> "<<(nf+1)<<endl;
  255. //#endif
  256. }
  257. int max_ne = -1;
  258. if(nfei >= 0)
  259. {
  260. max_ne = uE2E[EMAP(e)][nfei];
  261. }
  262. if(max_ne>=0)
  263. {
  264. // face of neighbor
  265. const int nf = max_ne%m;
  266. #ifdef IGL_OUTER_HULL_DEBUG
  267. if(!FH[nf])
  268. {
  269. // first time seeing face
  270. cout<<(f+1)<<" --> "<<(nf+1)<<endl;
  271. }
  272. #endif
  273. FH[nf] = true;
  274. //std::cout << "face " << face_count++ << ": " << nf << std::endl;
  275. //std::cout << "f " << F.row(nf).array()+1 << std::endl;
  276. FHcount++;
  277. // corner of neighbor
  278. const int nc = max_ne/m;
  279. const int nd = F(nf,(nc+2)%3);
  280. const bool cons = (flip(f)?fd:fs) == nd;
  281. flip(nf) = (cons ? flip(f) : !flip(f));
  282. //cout<<"flip("<<nf<<") = "<<(flip(nf)?"true":"false")<<endl;
  283. const int ne1 = nf+((nc+1)%3)*m;
  284. const int ne2 = nf+((nc+2)%3)*m;
  285. if(!EH[ne1])
  286. {
  287. Q.push(ne1);
  288. }
  289. if(!EH[ne2])
  290. {
  291. Q.push(ne2);
  292. }
  293. }
  294. }
  295. {
  296. vG[id].resize(FHcount,3);
  297. vJ[id].resize(FHcount,1);
  298. //nG += FHcount;
  299. size_t h = 0;
  300. assert(counts(id) == IM.rows());
  301. for(int i = 0;i<counts(id);i++)
  302. {
  303. const size_t f = IM(i);
  304. //if(f_flip)
  305. //{
  306. // flip(f) = !flip(f);
  307. //}
  308. if(FH[f])
  309. {
  310. vG[id].row(h) = (flip(f)?F.row(f).reverse().eval():F.row(f));
  311. vJ[id](h,0) = f;
  312. h++;
  313. }
  314. }
  315. assert((int)h == FHcount);
  316. }
  317. }
  318. // Is A inside B? Assuming A and B are consistently oriented but closed and
  319. // non-intersecting.
  320. const auto & has_overlapping_bbox = [](
  321. const Eigen::MatrixBase<DerivedV> & V,
  322. const MatrixXG & A,
  323. const MatrixXG & B)->bool
  324. {
  325. const auto & bounding_box = [](
  326. const Eigen::MatrixBase<DerivedV> & V,
  327. const MatrixXG & F)->
  328. PlainMatrix<DerivedV,2,3>
  329. {
  330. PlainMatrix<DerivedV,2,3> BB(2,3);
  331. BB<<
  332. 1e26,1e26,1e26,
  333. -1e26,-1e26,-1e26;
  334. const size_t m = F.rows();
  335. for(size_t f = 0;f<m;f++)
  336. {
  337. for(size_t c = 0;c<3;c++)
  338. {
  339. const auto & vfc = V.row(F(f,c)).eval();
  340. BB(0,0) = std::min(BB(0,0), vfc(0,0));
  341. BB(0,1) = std::min(BB(0,1), vfc(0,1));
  342. BB(0,2) = std::min(BB(0,2), vfc(0,2));
  343. BB(1,0) = std::max(BB(1,0), vfc(0,0));
  344. BB(1,1) = std::max(BB(1,1), vfc(0,1));
  345. BB(1,2) = std::max(BB(1,2), vfc(0,2));
  346. }
  347. }
  348. return BB;
  349. };
  350. // A lot of the time we're dealing with unrelated, distant components: cull
  351. // them.
  352. PlainMatrix<DerivedV,2,3> ABB = bounding_box(V,A);
  353. PlainMatrix<DerivedV,2,3> BBB = bounding_box(V,B);
  354. if( (BBB.row(0)-ABB.row(1)).maxCoeff()>0 ||
  355. (ABB.row(0)-BBB.row(1)).maxCoeff()>0 )
  356. {
  357. // bounding boxes do not overlap
  358. return false;
  359. } else {
  360. return true;
  361. }
  362. };
  363. // Reject components which are completely inside other components
  364. vector<bool> keep(ncc,true);
  365. size_t nG = 0;
  366. // This is O( ncc * ncc * m)
  367. for(size_t id = 0;id<ncc;id++)
  368. {
  369. if (!keep[id]) continue;
  370. std::vector<size_t> unresolved;
  371. for(size_t oid = 0;oid<ncc;oid++)
  372. {
  373. if(id == oid || !keep[oid])
  374. {
  375. continue;
  376. }
  377. if (has_overlapping_bbox(V, vG[id], vG[oid])) {
  378. unresolved.push_back(oid);
  379. }
  380. }
  381. const size_t num_unresolved_components = unresolved.size();
  382. PlainMatrix<DerivedV,Eigen::Dynamic,3> query_points(num_unresolved_components, 3);
  383. for (size_t i=0; i<num_unresolved_components; i++) {
  384. const size_t oid = unresolved[i];
  385. PlainMatrix<DerivedF,1> f = vG[oid].row(0);
  386. query_points(i,0) = (V(f(0,0), 0) + V(f(0,1), 0) + V(f(0,2), 0))/3.0;
  387. query_points(i,1) = (V(f(0,0), 1) + V(f(0,1), 1) + V(f(0,2), 1))/3.0;
  388. query_points(i,2) = (V(f(0,0), 2) + V(f(0,1), 2) + V(f(0,2), 2))/3.0;
  389. }
  390. Eigen::VectorXi inside;
  391. igl::copyleft::cgal::points_inside_component(V, vG[id], query_points, inside);
  392. assert((size_t)inside.size() == num_unresolved_components);
  393. for (size_t i=0; i<num_unresolved_components; i++) {
  394. if (inside(i, 0)) {
  395. const size_t oid = unresolved[i];
  396. keep[oid] = false;
  397. }
  398. }
  399. }
  400. for (size_t id = 0; id<ncc; id++) {
  401. if (keep[id]) {
  402. nG += vJ[id].rows();
  403. }
  404. }
  405. // collect G and J across components
  406. G.resize(nG,3);
  407. J.resize(nG,1);
  408. {
  409. size_t off = 0;
  410. for(Index id = 0;id<(Index)ncc;id++)
  411. {
  412. if(keep[id])
  413. {
  414. assert(vG[id].rows() == vJ[id].rows());
  415. G.block(off,0,vG[id].rows(),vG[id].cols()) = vG[id];
  416. J.block(off,0,vJ[id].rows(),vJ[id].cols()) = vJ[id];
  417. off += vG[id].rows();
  418. }
  419. }
  420. }
  421. }
  422. #ifdef IGL_STATIC_LIBRARY
  423. // Explicit template instantiation
  424. // generated by autoexplicit.sh
  425. template void igl::copyleft::cgal::outer_hull_legacy<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  426. template void igl::copyleft::cgal::outer_hull_legacy< 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::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > &);
  427. template void igl::copyleft::cgal::outer_hull_legacy<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::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  428. #ifdef WIN32
  429. #endif
  430. #endif