outer_hull.cpp 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643
  1. #include "outer_hull.h"
  2. #include "outer_facet.h"
  3. #include "sortrows.h"
  4. #include "facet_components.h"
  5. #include "winding_number.h"
  6. #include "triangle_triangle_adjacency.h"
  7. #include "unique_edge_map.h"
  8. #include "barycenter.h"
  9. #include "per_face_normals.h"
  10. #include "writePLY.h"
  11. #include "sort_angles.h"
  12. #include <Eigen/Geometry>
  13. #include <vector>
  14. #include <map>
  15. #include <queue>
  16. #include <iostream>
  17. #define IGL_OUTER_HULL_DEBUG
  18. template <
  19. typename DerivedV,
  20. typename DerivedF,
  21. typename DerivedN,
  22. typename DerivedG,
  23. typename DerivedJ,
  24. typename Derivedflip>
  25. IGL_INLINE void igl::outer_hull(
  26. const Eigen::PlainObjectBase<DerivedV> & V,
  27. const Eigen::PlainObjectBase<DerivedF> & F,
  28. const Eigen::PlainObjectBase<DerivedN> & N,
  29. Eigen::PlainObjectBase<DerivedG> & G,
  30. Eigen::PlainObjectBase<DerivedJ> & J,
  31. Eigen::PlainObjectBase<Derivedflip> & flip)
  32. {
  33. #ifdef IGL_OUTER_HULL_DEBUG
  34. std::cerr << "Extracting outer hull" << std::endl;
  35. writePLY("outer_hull_input.ply", V, F);
  36. #endif
  37. using namespace Eigen;
  38. using namespace std;
  39. typedef typename DerivedF::Index Index;
  40. Matrix<Index,DerivedF::RowsAtCompileTime,1> C;
  41. typedef Matrix<typename DerivedV::Scalar,Dynamic,DerivedV::ColsAtCompileTime> MatrixXV;
  42. typedef Matrix<typename DerivedF::Scalar,Dynamic,DerivedF::ColsAtCompileTime> MatrixXF;
  43. typedef Matrix<typename DerivedG::Scalar,Dynamic,DerivedG::ColsAtCompileTime> MatrixXG;
  44. typedef Matrix<typename DerivedJ::Scalar,Dynamic,DerivedJ::ColsAtCompileTime> MatrixXJ;
  45. typedef Matrix<typename DerivedN::Scalar,1,3> RowVector3N;
  46. const Index m = F.rows();
  47. const auto & duplicate_simplex = [&F](const int f, const int g)->bool
  48. {
  49. return
  50. (F(f,0) == F(g,0) && F(f,1) == F(g,1) && F(f,2) == F(g,2)) ||
  51. (F(f,1) == F(g,0) && F(f,2) == F(g,1) && F(f,0) == F(g,2)) ||
  52. (F(f,2) == F(g,0) && F(f,0) == F(g,1) && F(f,1) == F(g,2)) ||
  53. (F(f,0) == F(g,2) && F(f,1) == F(g,1) && F(f,2) == F(g,0)) ||
  54. (F(f,1) == F(g,2) && F(f,2) == F(g,1) && F(f,0) == F(g,0)) ||
  55. (F(f,2) == F(g,2) && F(f,0) == F(g,1) && F(f,1) == F(g,0));
  56. };
  57. #ifdef IGL_OUTER_HULL_DEBUG
  58. cout<<"outer hull..."<<endl;
  59. #endif
  60. #ifdef IGL_OUTER_HULL_DEBUG
  61. cout<<"edge map..."<<endl;
  62. #endif
  63. typedef Matrix<typename DerivedF::Scalar,Dynamic,2> MatrixX2I;
  64. typedef Matrix<typename DerivedF::Index,Dynamic,1> VectorXI;
  65. typedef Matrix<typename DerivedV::Scalar, 3, 1> Vector3F;
  66. MatrixX2I E,uE;
  67. VectorXI EMAP;
  68. vector<vector<typename DerivedF::Index> > uE2E;
  69. unique_edge_map(F,E,uE,EMAP,uE2E);
  70. #ifdef IGL_OUTER_HULL_DEBUG
  71. for (size_t ui=0; ui<uE.rows(); ui++) {
  72. std::cout << ui << ": " << uE2E[ui].size() << " -- (";
  73. for (size_t i=0; i<uE2E[ui].size(); i++) {
  74. std::cout << uE2E[ui][i] << ", ";
  75. }
  76. std::cout << ")" << std::endl;
  77. }
  78. #endif
  79. // TODO:
  80. // uE --> face-edge index, sorted CCW around edge according to normal
  81. // uE --> sorted order index
  82. // uE --> bool, whether needed to flip face to make "consistent" with unique
  83. // edge
  84. // Place order of each half-edge in its corresponding sorted list around edge
  85. VectorXI diIM(3*m);
  86. // Whether face's edge used for sorting is consistent with unique edge
  87. VectorXI dicons(3*m);
  88. // dihedral angles of faces around edge with face of edge in dicons
  89. vector<vector<typename DerivedV::Scalar> > di(uE2E.size());
  90. // For each list of face-edges incide on a unique edge
  91. for(size_t ui = 0;ui<(size_t)uE.rows();ui++)
  92. {
  93. // Base normal vector to orient against
  94. const auto fe0 = uE2E[ui][0];
  95. const RowVector3N & eVp = N.row(fe0%m);
  96. MatrixXd di_I(uE2E[ui].size(),3);
  97. const typename DerivedF::Scalar o = F(fe0%m, fe0/m);
  98. const typename DerivedF::Scalar d = F(fe0%m,((fe0/m)+2)%3);
  99. const typename DerivedF::Scalar s = F(fe0%m,((fe0/m)+1)%3);
  100. // Edge vector
  101. auto eV = (V.row(d)-V.row(s)).normalized();
  102. auto edge_len = (V.row(d) - V.row(s)).norm();
  103. auto edge_valance = uE2E[ui].size();
  104. assert(edge_valance % 2 == 0);
  105. bool degenerated = !eV.allFinite() || edge_len < 1e-12;
  106. #ifdef IGL_OUTER_HULL_DEBUG
  107. if (degenerated && edge_valance > 2) {
  108. cerr.precision(30);
  109. std::cerr << ui << ": " << (V.row(d) - V.row(s)).norm() << std::endl;
  110. std::cerr << "Edge valance: " << edge_valance << std::endl;
  111. std::cerr << V.row(d) << std::endl;
  112. std::cerr << V.row(s) << std::endl;
  113. }
  114. #endif
  115. if (degenerated) {
  116. const size_t num_adj_faces = uE2E[ui].size();
  117. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 3>
  118. normals(num_adj_faces, 3);
  119. for (size_t fei=0; fei<num_adj_faces; fei++) {
  120. const auto & fe = uE2E[ui][fei];
  121. const auto f = fe % m;
  122. const RowVector3N & n = N.row(f);
  123. normals.row(fei) = n;
  124. }
  125. for (size_t i=0; i<num_adj_faces; i++) {
  126. size_t j = (i+1) % num_adj_faces;
  127. eV = normals.row(i).cross(normals.row(j));
  128. auto length = eV.norm();
  129. if (length > 1e-12) {
  130. eV /= length;
  131. break;
  132. }
  133. }
  134. }
  135. if (!eV.allFinite() || eV.norm() < 1e-12) {
  136. //cerr << "This is bad... all adj face normals are colinear" << std::endl;
  137. eV.setZero();
  138. }
  139. if (degenerated){
  140. // Adjust edge direction.
  141. Vector3F in_face_vec = V.row(o) - V.row(s);
  142. Vector3F edge = eV;
  143. if (edge.cross(in_face_vec).dot(eVp) < 0) {
  144. #ifdef IGL_OUTER_HULL_DEBUG
  145. cerr << "Flipping edge..." << std::endl;
  146. #endif
  147. eV *= -1;
  148. }
  149. //cerr << "Resolved: " << eV << std::endl;
  150. }
  151. vector<bool> cons(uE2E[ui].size());
  152. // Loop over incident face edges
  153. for(size_t fei = 0;fei<uE2E[ui].size();fei++)
  154. {
  155. const auto & fe = uE2E[ui][fei];
  156. const auto f = fe % m;
  157. const auto c = fe / m;
  158. // source should match destination to be consistent
  159. cons[fei] = (d == F(f,(c+1)%3));
  160. assert( cons[fei] || (d == F(f,(c+2)%3)));
  161. assert(!cons[fei] || (s == F(f,(c+2)%3)));
  162. assert(!cons[fei] || (d == F(f,(c+1)%3)));
  163. // Angle between n and f
  164. const RowVector3N & n = N.row(f);
  165. di_I(fei, 0) = eVp.cross(n).dot(eV);
  166. di_I(fei, 1) = eVp.dot(n);
  167. assert(di_I(fei,0) == di_I(fei,0) && "NaN Alert!");
  168. assert(di_I(fei,1) == di_I(fei,1) && "NaN Alert!");
  169. if (cons[fei]) {
  170. di_I(fei, 0) *= -1;
  171. di_I(fei, 1) *= -1;
  172. }
  173. di_I(fei, 0) *= -1; // Sort clockwise.
  174. // This signing is very important to make sure different edges sort
  175. // duplicate faces the same way, regardless of their orientations
  176. di_I(fei,2) = (cons[fei]?1.:-1.)*(f+1);
  177. }
  178. #if 0
  179. // Despite the effort to get stable normals the atan2 up doesn't
  180. // compute (exactly) -θ for -n if it computes θ for n. So just
  181. // explicitly check if there's a duplicate face
  182. // Shitty O(val^2) implementation
  183. for(size_t fei = 0;fei<uE2E[ui].size();fei++)
  184. {
  185. const auto & fe = uE2E[ui][fei];
  186. const auto f = fe % m;
  187. for(size_t gei = fei+1;gei<uE2E[ui].size();gei++)
  188. {
  189. const auto & ge = uE2E[ui][gei];
  190. const auto g = ge % m;
  191. if(duplicate_simplex(f,g))
  192. {
  193. #ifdef IGL_OUTER_HULL_DEBUG
  194. cout<<"Forcing duplicate: "<<(f+1)<<","<<(g+1)<<endl;
  195. #endif
  196. di_I(gei,0) = di_I(fei,0);
  197. }
  198. }
  199. }
  200. #endif
  201. VectorXi IM;
  202. //igl::sort(di[ui],true,di[ui],IM);
  203. // Sort, but break ties using "signed index" to ensure that duplicates
  204. // always show up in same order.
  205. igl::sort_angles(di_I, IM);
  206. vector<typename DerivedF::Index> temp = uE2E[ui];
  207. std::cout.precision(20);
  208. std::cout << "sorted" << std::endl;
  209. for(size_t fei = 0;fei<uE2E[ui].size();fei++)
  210. {
  211. std::cout << di_I.row(IM(fei)) << std::endl;
  212. uE2E[ui][fei] = temp[IM(fei)];
  213. const auto & fe = uE2E[ui][fei];
  214. diIM(fe) = fei;
  215. dicons(fe) = cons[IM(fei)];
  216. }
  217. //MatrixXd s_di_I;
  218. //igl::sortrows(di_I,true,s_di_I,IM);
  219. //di[ui].resize(uE2E[ui].size());
  220. //for(size_t i = 0;i<di[ui].size();i++)
  221. //{
  222. // di[ui][i] = s_di_I(i,0);
  223. //}
  224. //// copy old list
  225. //vector<typename DerivedF::Index> temp = uE2E[ui];
  226. //for(size_t fei = 0;fei<uE2E[ui].size();fei++)
  227. //{
  228. // uE2E[ui][fei] = temp[IM(fei)];
  229. // const auto & fe = uE2E[ui][fei];
  230. // diIM(fe) = fei;
  231. // dicons(fe) = cons[IM(fei)];
  232. //}
  233. }
  234. vector<vector<vector<Index > > > TT,_1;
  235. triangle_triangle_adjacency(E,EMAP,uE2E,false,TT,_1);
  236. VectorXI counts;
  237. #ifdef IGL_OUTER_HULL_DEBUG
  238. cout<<"facet components..."<<endl;
  239. #endif
  240. facet_components(TT,C,counts);
  241. assert(C.maxCoeff()+1 == counts.rows());
  242. const size_t ncc = counts.rows();
  243. G.resize(0,F.cols());
  244. J.resize(0,1);
  245. flip.setConstant(m,1,false);
  246. #ifdef IGL_OUTER_HULL_DEBUG
  247. cout<<"reindex..."<<endl;
  248. #endif
  249. // H contains list of faces on outer hull;
  250. vector<bool> FH(m,false);
  251. vector<bool> EH(3*m,false);
  252. vector<MatrixXG> vG(ncc);
  253. vector<MatrixXJ> vJ(ncc);
  254. vector<MatrixXJ> vIM(ncc);
  255. for(size_t id = 0;id<ncc;id++)
  256. {
  257. vIM[id].resize(counts[id],1);
  258. }
  259. // current index into each IM
  260. vector<size_t> g(ncc,0);
  261. // place order of each face in its respective component
  262. for(Index f = 0;f<m;f++)
  263. {
  264. vIM[C(f)](g[C(f)]++) = f;
  265. }
  266. #ifdef IGL_OUTER_HULL_DEBUG
  267. cout<<"barycenters..."<<endl;
  268. #endif
  269. // assumes that "resolve" has handled any coplanar cases correctly and nearly
  270. // coplanar cases can be sorted based on barycenter.
  271. MatrixXV BC;
  272. barycenter(V,F,BC);
  273. #ifdef IGL_OUTER_HULL_DEBUG
  274. cout<<"loop over CCs (="<<ncc<<")..."<<endl;
  275. #endif
  276. for(Index id = 0;id<(Index)ncc;id++)
  277. {
  278. auto & IM = vIM[id];
  279. // starting face that's guaranteed to be on the outer hull and in this
  280. // component
  281. int f;
  282. bool f_flip;
  283. #ifdef IGL_OUTER_HULL_DEBUG
  284. cout<<"outer facet..."<<endl;
  285. #endif
  286. outer_facet(V,F,N,IM,f,f_flip);
  287. #ifdef IGL_OUTER_HULL_DEBUG
  288. cout<<"outer facet: "<<f<<endl;
  289. #endif
  290. int FHcount = 1;
  291. FH[f] = true;
  292. // Q contains list of face edges to continue traversing upong
  293. queue<int> Q;
  294. Q.push(f+0*m);
  295. Q.push(f+1*m);
  296. Q.push(f+2*m);
  297. flip(f) = f_flip;
  298. //cout<<"flip("<<f<<") = "<<(flip(f)?"true":"false")<<endl;
  299. #ifdef IGL_OUTER_HULL_DEBUG
  300. cout<<"BFS..."<<endl;
  301. #endif
  302. while(!Q.empty())
  303. {
  304. // face-edge
  305. const int e = Q.front();
  306. Q.pop();
  307. // face
  308. const int f = e%m;
  309. // corner
  310. const int c = e/m;
  311. // Should never see edge again...
  312. if(EH[e] == true)
  313. {
  314. continue;
  315. }
  316. EH[e] = true;
  317. // source of edge according to f
  318. const int fs = flip(f)?F(f,(c+2)%3):F(f,(c+1)%3);
  319. // destination of edge according to f
  320. const int fd = flip(f)?F(f,(c+1)%3):F(f,(c+2)%3);
  321. // edge valence
  322. const size_t val = uE2E[EMAP(e)].size();
  323. //// find overlapping face-edges
  324. //const auto & neighbors = uE2E[EMAP(e)];
  325. //// normal after possible flipping
  326. //const auto & fN = (flip(f)?-1.:1.)*N.row(f);
  327. //// Edge vector according to f's (flipped) orientation.
  328. ////const auto & eV = (V.row(fd)-V.row(fs)).normalized();
  329. //#warning "EXPERIMENTAL, DO NOT USE"
  330. //// THIS IS WRONG! The first face is---after sorting---no longer the face
  331. //// used for orienting the sort.
  332. //const auto ui = EMAP(e);
  333. //const auto fe0 = uE2E[ui][0];
  334. //const auto es = F(fe0%m,((fe0/m)+1)%3);
  335. // is edge consistent with edge of face used for sorting
  336. const int e_cons = (dicons(e) ? 1: -1);
  337. int nfei = -1;
  338. // Loop once around trying to find suitable next face
  339. for(size_t step = 1; step<val+2;step++)
  340. {
  341. const int nfei_new = (diIM(e) + 2*val + e_cons*step*(flip(f)?-1:1))%val;
  342. const int nf = uE2E[EMAP(e)][nfei_new] % m;
  343. // Don't consider faces with identical dihedral angles
  344. //if((di[EMAP(e)][diIM(e)] != di[EMAP(e)][nfei_new]))
  345. //#warning "THIS IS HACK, FIX ME"
  346. // if( abs(di[EMAP(e)][diIM(e)] - di[EMAP(e)][nfei_new]) < 1e-16 )
  347. {
  348. #ifdef IGL_OUTER_HULL_DEBUG
  349. //cout<<"Next facet: "<<(f+1)<<" --> "<<(nf+1)<<", |"<<
  350. // di[EMAP(e)][diIM(e)]<<" - "<<di[EMAP(e)][nfei_new]<<"| = "<<
  351. // abs(di[EMAP(e)][diIM(e)] - di[EMAP(e)][nfei_new])
  352. // <<endl;
  353. #endif
  354. // Only use this face if not already seen
  355. if(!FH[nf])
  356. {
  357. nfei = nfei_new;
  358. }
  359. break;
  360. }
  361. //#ifdef IGL_OUTER_HULL_DEBUG
  362. // cout<<"Skipping co-planar facet: "<<(f+1)<<" --> "<<(nf+1)<<endl;
  363. //#endif
  364. }
  365. int max_ne = -1;
  366. //// Loop over and find max dihedral angle
  367. //typename DerivedV::Scalar max_di = -1;
  368. //for(const auto & ne : neighbors)
  369. //{
  370. // const int nf = ne%m;
  371. // if(nf == f)
  372. // {
  373. // continue;
  374. // }
  375. // // Corner of neighbor
  376. // const int nc = ne/m;
  377. // // Is neighbor oriented consistently with (flipped) f?
  378. // //const int ns = F(nf,(nc+1)%3);
  379. // const int nd = F(nf,(nc+2)%3);
  380. // const bool cons = (flip(f)?fd:fs) == nd;
  381. // // Normal after possibly flipping to match flip or orientation of f
  382. // const auto & nN = (cons? (flip(f)?-1:1.) : (flip(f)?1.:-1.) )*N.row(nf);
  383. // // Angle between n and f
  384. // const auto & ndi = M_PI - atan2( fN.cross(nN).dot(eV), fN.dot(nN));
  385. // if(ndi>=max_di)
  386. // {
  387. // max_ne = ne;
  388. // max_di = ndi;
  389. // }
  390. //}
  391. ////cout<<(max_ne != max_ne_2)<<" =?= "<<e_cons<<endl;
  392. //if(max_ne != max_ne_2)
  393. //{
  394. // cout<<(f+1)<<" ---> "<<(max_ne%m)+1<<" != "<<(max_ne_2%m)+1<<" ... "<<e_cons<<" "<<flip(f)<<endl;
  395. // typename DerivedV::Scalar max_di = -1;
  396. // for(size_t nei = 0;nei<neighbors.size();nei++)
  397. // {
  398. // const auto & ne = neighbors[nei];
  399. // const int nf = ne%m;
  400. // if(nf == f)
  401. // {
  402. // cout<<" "<<(ne%m)+1<<":\t"<<0<<"\t"<<di[EMAP[e]][nei]<<" "<<diIM(ne)<<endl;
  403. // continue;
  404. // }
  405. // // Corner of neighbor
  406. // const int nc = ne/m;
  407. // // Is neighbor oriented consistently with (flipped) f?
  408. // //const int ns = F(nf,(nc+1)%3);
  409. // const int nd = F(nf,(nc+2)%3);
  410. // const bool cons = (flip(f)?fd:fs) == nd;
  411. // // Normal after possibly flipping to match flip or orientation of f
  412. // const auto & nN = (cons? (flip(f)?-1:1.) : (flip(f)?1.:-1.) )*N.row(nf);
  413. // // Angle between n and f
  414. // const auto & ndi = M_PI - atan2( fN.cross(nN).dot(eV), fN.dot(nN));
  415. // cout<<" "<<(ne%m)+1<<":\t"<<ndi<<"\t"<<di[EMAP[e]][nei]<<" "<<diIM(ne)<<endl;
  416. // if(ndi>=max_di)
  417. // {
  418. // max_ne = ne;
  419. // max_di = ndi;
  420. // }
  421. // }
  422. //}
  423. if(nfei >= 0)
  424. {
  425. max_ne = uE2E[EMAP(e)][nfei];
  426. }
  427. if(max_ne>=0)
  428. {
  429. // face of neighbor
  430. const int nf = max_ne%m;
  431. #ifdef IGL_OUTER_HULL_DEBUG
  432. if(!FH[nf])
  433. {
  434. // first time seeing face
  435. cout<<(f+1)<<" --> "<<(nf+1)<<endl;
  436. }
  437. #endif
  438. FH[nf] = true;
  439. FHcount++;
  440. // corner of neighbor
  441. const int nc = max_ne/m;
  442. const int nd = F(nf,(nc+2)%3);
  443. const bool cons = (flip(f)?fd:fs) == nd;
  444. flip(nf) = (cons ? flip(f) : !flip(f));
  445. //cout<<"flip("<<nf<<") = "<<(flip(nf)?"true":"false")<<endl;
  446. const int ne1 = nf+((nc+1)%3)*m;
  447. const int ne2 = nf+((nc+2)%3)*m;
  448. if(!EH[ne1])
  449. {
  450. Q.push(ne1);
  451. }
  452. if(!EH[ne2])
  453. {
  454. Q.push(ne2);
  455. }
  456. }
  457. }
  458. {
  459. vG[id].resize(FHcount,3);
  460. vJ[id].resize(FHcount,1);
  461. //nG += FHcount;
  462. size_t h = 0;
  463. assert(counts(id) == IM.rows());
  464. for(int i = 0;i<counts(id);i++)
  465. {
  466. const size_t f = IM(i);
  467. //if(f_flip)
  468. //{
  469. // flip(f) = !flip(f);
  470. //}
  471. if(FH[f])
  472. {
  473. vG[id].row(h) = (flip(f)?F.row(f).reverse().eval():F.row(f));
  474. vJ[id](h,0) = f;
  475. h++;
  476. }
  477. }
  478. assert((int)h == FHcount);
  479. }
  480. }
  481. // Is A inside B? Assuming A and B are consistently oriented but closed and
  482. // non-intersecting.
  483. const auto & is_component_inside_other = [](
  484. const Eigen::PlainObjectBase<DerivedV> & V,
  485. const MatrixXV & BC,
  486. const MatrixXG & A,
  487. const MatrixXJ & AJ,
  488. const MatrixXG & B)->bool
  489. {
  490. const auto & bounding_box = [](
  491. const Eigen::PlainObjectBase<DerivedV> & V,
  492. const MatrixXG & F)->
  493. MatrixXV
  494. {
  495. MatrixXV BB(2,3);
  496. BB<<
  497. 1e26,1e26,1e26,
  498. -1e26,-1e26,-1e26;
  499. const size_t m = F.rows();
  500. for(size_t f = 0;f<m;f++)
  501. {
  502. for(size_t c = 0;c<3;c++)
  503. {
  504. const auto & vfc = V.row(F(f,c));
  505. BB.row(0) = BB.row(0).array().min(vfc.array()).eval();
  506. BB.row(1) = BB.row(1).array().max(vfc.array()).eval();
  507. }
  508. }
  509. return BB;
  510. };
  511. // A lot of the time we're dealing with unrelated, distant components: cull
  512. // them.
  513. MatrixXV ABB = bounding_box(V,A);
  514. MatrixXV BBB = bounding_box(V,B);
  515. if( (BBB.row(0)-ABB.row(1)).maxCoeff()>0 ||
  516. (ABB.row(0)-BBB.row(1)).maxCoeff()>0 )
  517. {
  518. // bounding boxes do not overlap
  519. return false;
  520. }
  521. ////////////////////////////////////////////////////////////////////////
  522. // POTENTIAL ROBUSTNESS WEAK AREA
  523. ////////////////////////////////////////////////////////////////////////
  524. //
  525. // q could be so close (<~1e-16) to B that the winding number is not a robust way to
  526. // determine inside/outsideness. We could try to find a _better_ q which is
  527. // farther away, but couldn't they all be bad?
  528. MatrixXV q = BC.row(AJ(0));
  529. // In a perfect world, it's enough to test a single point.
  530. double w;
  531. // winding_number_3 expects colmajor
  532. const typename DerivedV::Scalar * Vdata;
  533. Vdata = V.data();
  534. Matrix<
  535. typename DerivedV::Scalar,
  536. DerivedV::RowsAtCompileTime,
  537. DerivedV::ColsAtCompileTime,
  538. ColMajor> Vcol;
  539. if(DerivedV::IsRowMajor)
  540. {
  541. // copy to convert to colmajor
  542. Vcol = V;
  543. Vdata = Vcol.data();
  544. }
  545. winding_number_3(
  546. Vdata,V.rows(),
  547. B.data(),B.rows(),
  548. q.data(),1,&w);
  549. return fabs(w)>0.5;
  550. };
  551. // Reject components which are completely inside other components
  552. vector<bool> keep(ncc,true);
  553. size_t nG = 0;
  554. // This is O( ncc * ncc * m)
  555. for(size_t id = 0;id<ncc;id++)
  556. {
  557. for(size_t oid = 0;oid<ncc;oid++)
  558. {
  559. if(id == oid)
  560. {
  561. continue;
  562. }
  563. const bool inside = is_component_inside_other(V,BC,vG[id],vJ[id],vG[oid]);
  564. #ifdef IGL_OUTER_HULL_DEBUG
  565. cout<<id<<" is inside "<<oid<<" ? "<<inside<<endl;
  566. #endif
  567. keep[id] = keep[id] && !inside;
  568. }
  569. if(keep[id])
  570. {
  571. nG += vJ[id].rows();
  572. }
  573. }
  574. // collect G and J across components
  575. G.resize(nG,3);
  576. J.resize(nG,1);
  577. {
  578. size_t off = 0;
  579. for(Index id = 0;id<(Index)ncc;id++)
  580. {
  581. if(keep[id])
  582. {
  583. assert(vG[id].rows() == vJ[id].rows());
  584. G.block(off,0,vG[id].rows(),vG[id].cols()) = vG[id];
  585. J.block(off,0,vJ[id].rows(),vJ[id].cols()) = vJ[id];
  586. off += vG[id].rows();
  587. }
  588. }
  589. }
  590. }
  591. template <
  592. typename DerivedV,
  593. typename DerivedF,
  594. typename DerivedG,
  595. typename DerivedJ,
  596. typename Derivedflip>
  597. IGL_INLINE void igl::outer_hull(
  598. const Eigen::PlainObjectBase<DerivedV> & V,
  599. const Eigen::PlainObjectBase<DerivedF> & F,
  600. Eigen::PlainObjectBase<DerivedG> & G,
  601. Eigen::PlainObjectBase<DerivedJ> & J,
  602. Eigen::PlainObjectBase<Derivedflip> & flip)
  603. {
  604. Eigen::Matrix<typename DerivedV::Scalar,DerivedF::RowsAtCompileTime,3> N;
  605. per_face_normals_stable(V,F,N);
  606. return outer_hull(V,F,N,G,J,flip);
  607. }
  608. #ifdef IGL_STATIC_LIBRARY
  609. // Explicit template specialization
  610. template void igl::outer_hull<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<bool, -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&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<bool, -1, 1, 0, -1, 1> >&);
  611. template void igl::outer_hull<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::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, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  612. template void igl::outer_hull<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, 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<bool, -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&, Eigen::PlainObjectBase<Eigen::Matrix<double, -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<bool, -1, 1, 0, -1, 1> >&);
  613. #endif