outer_hull.cpp 24 KB

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