split_nonmanifold.cpp 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2022 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 "split_nonmanifold.h"
  9. #include "unique_edge_map.h"
  10. #include "connected_components.h"
  11. #include "unique.h"
  12. #include "sort.h"
  13. #include "triangle_triangle_adjacency.h"
  14. #include "is_edge_manifold.h"
  15. #include <unordered_map>
  16. #include <cassert>
  17. #include <type_traits>
  18. #include "is_vertex_manifold.h"
  19. #include "matlab_format.h"
  20. #include <iostream>
  21. #include <unordered_set>
  22. #include <utility>
  23. template <
  24. typename DerivedF,
  25. typename DerivedSF,
  26. typename DerivedSVI
  27. >
  28. IGL_INLINE void igl::split_nonmanifold(
  29. const Eigen::MatrixBase<DerivedF> & F,
  30. Eigen::PlainObjectBase <DerivedSF> & SF,
  31. Eigen::PlainObjectBase <DerivedSVI> & SVI)
  32. {
  33. using Scalar = typename DerivedSF::Scalar;
  34. // Scalar must allow negative values
  35. static_assert(std::is_signed<Scalar>::value,"Scalar must be signed");
  36. using MatrixX2I = Eigen::Matrix<Scalar,Eigen::Dynamic,2>;
  37. using MatrixX3I = Eigen::Matrix<Scalar,Eigen::Dynamic,3>;
  38. using VectorXI = Eigen::Matrix< Scalar,Eigen::Dynamic,1>;
  39. MatrixX2I E,uE;
  40. VectorXI EMAP,uEC,uEE;
  41. igl::unique_edge_map(F,E,uE,EMAP,uEC,uEE);
  42. // Let's assume the most convenient connectivity data structure and worry
  43. // about performance later
  44. // There are always 3#F "corners".
  45. //
  46. // V[c] = v means that corner c is mapped to new-vertex v
  47. // Start with all corners mapped to singleton new-vertices
  48. Eigen::VectorXi V = Eigen::VectorXi::LinSpaced(F.size(),0,F.size()-1);
  49. // Convenience map so that CF(f,i) = V[c] = v where c is the ith corner of
  50. // face f.
  51. Eigen::Map<Eigen::MatrixXi> CF = Eigen::Map<Eigen::MatrixXi>(V.data(),F.rows(),F.cols());
  52. // C[v][j] = c means that c is the jth corner in the group of corners at
  53. // new-vertex v. As we merge these, we will clear "dead" new-vertices.
  54. std::vector<std::vector<int> > C(F.size());
  55. for(int i = 0;i<F.size();i++) { C[i] = {i}; }
  56. const int m = F.rows();
  57. // O(S) where S = |star(v)|
  58. // @param[in] v new-vertex index
  59. // @return list of face indices incident on new-vertex v
  60. const auto star = [&](const int v)->std::vector<int>
  61. {
  62. std::vector<int> faces(C[v].size());
  63. for(int i = 0;i<C[v].size();i++)
  64. {
  65. faces[i] = C[v][i]%m;
  66. }
  67. return faces;
  68. };
  69. // O(S) where S = |star(v)|
  70. // @param[in] v new-vertex index
  71. // @return list of half-edge indices incident on new-vertex v
  72. const auto nonmanifold_edge_star = [&](const int v)->std::vector<int>
  73. {
  74. std::vector<int> edges;
  75. // loop over edges opposite corners of v
  76. for(int e : C[v])
  77. {
  78. const int f = e%m;
  79. for(int j = 1;j<3;j++)
  80. {
  81. // next edge
  82. const int e1 = (e+j*m)%(3*m);
  83. const int u1 = EMAP(e1);
  84. if(uEC(u1+1)-uEC(u1) > 2)
  85. {
  86. edges.push_back(e1);
  87. }
  88. }
  89. }
  90. return edges;
  91. };
  92. // O(S) where S = |star(v)|
  93. const std::function<void(
  94. Eigen::VectorXi &,
  95. std::vector<std::vector<int> > &,
  96. const int, const int)> merge_vertex =
  97. [&merge_vertex](Eigen::VectorXi & V,
  98. std::vector<std::vector<int> > & C,
  99. const int u, const int v)
  100. {
  101. if(u == v) { return; }
  102. if(u > v) { merge_vertex(V,C,v,u); return; }
  103. assert(u < v);
  104. // Consider each corner in v
  105. for(const int c : C[v])
  106. {
  107. V[c] = u;
  108. }
  109. // Merge C[v] into C[u]
  110. C[u].insert(C[u].end(),C[v].begin(),C[v].end());
  111. C[v].clear();
  112. };
  113. // O(S) where S is the size of the star of e's first vertex.
  114. // This could probably be O(1) with careful bookkeeping
  115. const auto is_boundary = [&](const int e)->bool
  116. {
  117. // e----d
  118. // \ |
  119. // \f₁↑
  120. // \ |
  121. // s
  122. const int s = (e+1*m)%(3*m);
  123. const int d = (e+2*m)%(3*m);
  124. const int f = e%m;
  125. const int vs = V[s];
  126. const int vd = V[d];
  127. // Consider every face in the star of s
  128. for(const int g : star(vs))
  129. {
  130. if(g == f) { continue; }
  131. // Consider each edge in g
  132. for(int i = 0;i<3;i++)
  133. {
  134. const int a = (g+(i+1)*m)%(3*m);
  135. const int b = (g+(i+2)*m)%(3*m);
  136. // Is that edge the same as e?
  137. if(V[a] == vd && V[b] == vs) { return false; }
  138. if(V[a] == vs && V[b] == vd) { return false; }
  139. }
  140. }
  141. return true;
  142. };
  143. // Ω(m) and probably O(m log m) or worse.
  144. // This should take in the candidate merge edge pair, extract the submesh and
  145. // just check if that's manifold. Then it would be O(S) where S is the size of
  146. // biggest star of the edges' vertices.
  147. //
  148. // My guess is that is_edge_manifold is O(m) but is_vertex_manifold is
  149. // O(max(F))
  150. const auto is_manifold = [](Eigen::MatrixXi F)->bool
  151. {
  152. Eigen::Array<bool,Eigen::Dynamic,1> referenced =
  153. Eigen::Array<bool,Eigen::Dynamic,1>::Zero(F.maxCoeff()+1,1);
  154. for(int i = 0;i<F.size();i++)
  155. {
  156. referenced(F(i)) = true;
  157. }
  158. Eigen::Array<bool,Eigen::Dynamic,1> VM;
  159. igl::is_vertex_manifold(F,VM);
  160. for(int i = 0;i<VM.size();i++)
  161. {
  162. if(referenced(i) && !VM(i))
  163. {
  164. return false;
  165. }
  166. }
  167. return igl::is_edge_manifold(F);
  168. };
  169. // Ω(S) where S is the largest star of (vs1,vd2) or (vd1,vs2)
  170. // I think that is_vertex/edge_manifold(L) is O(|L| log |L|) so I think that
  171. // should make this O(|S| log |S|) with some gross constants because of all
  172. // the copying and sorting things into different data structures.
  173. //
  174. // merging edges (vs1,vd2) and (vd1,vs2) requires merging vertices (vs1→vd1) and
  175. // (vd2→vd2).
  176. //
  177. // Merging vertices (a→b) will change and only change the stars of a and b.
  178. // That is, some vertex c ≠ a,b will have the sam star before and after.
  179. //
  180. // Whether a vertex is singular depends entirely on its star.
  181. //
  182. // Therefore, the only vertices we need to check for non-manifoldness are
  183. // vs=(vs1,vd2) and vd=(vd1,vs2).
  184. const auto simulated_merge_is_manifold =
  185. [&](
  186. const int vs1, const int vd2,
  187. const int vd1, const int vs2)->bool
  188. {
  189. // all_faces[i] = f means that f is the ith face in the list of stars.
  190. std::vector<int> all_faces;
  191. for(int v : {vs1,vd2,vd1,vs2})
  192. {
  193. std::vector<int> star_v = star(v);
  194. all_faces.insert(all_faces.end(),star_v.begin(),star_v.end());
  195. }
  196. // unique_faces[l] = f means that f is the lth unique face in the list of
  197. // stars.
  198. std::vector<int> unique_faces;
  199. std::vector<size_t> _, local;
  200. igl::unique(all_faces,unique_faces,_,local);
  201. Eigen::MatrixXi L(unique_faces.size(),3);
  202. // collect local faces
  203. for(int l = 0;l<unique_faces.size();l++)
  204. {
  205. L.row(l) = CF.row(unique_faces[l]);
  206. }
  207. {
  208. int f = 0;
  209. const auto merge_local = [&](const int v1, const int v2)
  210. {
  211. const int u = std::min(v1,v2);
  212. for(const int v : {v1,v2})
  213. {
  214. for(const int c : C[v])
  215. {
  216. const int i = c/m;
  217. L(local[f++],i) = u;
  218. }
  219. }
  220. };
  221. // must match order {vs1,vd2,vd1,vs2} above
  222. merge_local(vs1,vd2);
  223. merge_local(vd1,vs2);
  224. }
  225. // remove unreferenced vertices by mapping each index in L to a unique
  226. // index between 0 and size(unique(L))
  227. std::unordered_map<int,int> M;
  228. for(int & i : L.reshaped())
  229. {
  230. if(M.find(i) == M.end())
  231. {
  232. M[i] = M.size();
  233. }
  234. i = M[i];
  235. }
  236. // Only need to check if the two vertices being merged are manifold
  237. Eigen::Array<bool,Eigen::Dynamic,1> VM;
  238. const int vs = std::min(vs1,vd2);
  239. const int vd = std::min(vd1,vs2);
  240. igl::is_vertex_manifold(L,VM);
  241. if(!VM(M[vs])) {
  242. return false;
  243. }
  244. if(!VM(M[vd])) {
  245. return false;
  246. }
  247. // Probably only need to check incident edges in star, but this also
  248. // checks link
  249. return igl::is_edge_manifold(L);
  250. };
  251. const auto merge_edge = [&](const int e1, const int e2)
  252. {
  253. // Ideally we would track whether an edge is a boundary so we can just
  254. // assert these. But because of "implied stitches" it's not necessarily just
  255. // e1 and e2 which become non-boundary when e1 and e2 are merged.
  256. //assert(is_boundary(e1));
  257. //assert(is_boundary(e2));
  258. if(!is_boundary(e1) || !is_boundary(e2)) { return false; }
  259. assert(e1 != e2);
  260. if(EMAP(e1) != EMAP(e2)) { return false; }
  261. assert(EMAP(e1) == EMAP(e2));
  262. const int u = EMAP(e1);
  263. const bool consistent = E(e1,0) == E(e2,1);
  264. // skip if inconsistently oriented
  265. if(!consistent) { return false; }
  266. // The code below is assuming merging consistently oriented edges
  267. if(E(e1,1) != E(e2,0))
  268. {
  269. }
  270. assert(E(e1,1) == E(e2,0));
  271. //
  272. // e1--d1 s2--e2
  273. // \ | | /
  274. // \f₁↑ ↓f₂/
  275. // \ | | /
  276. // s1 d2
  277. //
  278. //
  279. // "Cutting and Stitching: Converting Sets of Polygons to Manifold
  280. // Surfaces" [Guéziec et al. 2001]
  281. const int s1 = (e1+1*m)%(3*m);
  282. const int d1 = (e1+2*m)%(3*m);
  283. #ifndef NDEBUG
  284. {
  285. const int f1 = e1 % m;
  286. const int i1 = e1 / m;
  287. const int s1_test = f1 + ((i1+1)%3)*m;
  288. const int d1_test = f1 + ((i1+2)%3)*m;
  289. assert(s1 == s1_test);
  290. assert(d1 == d1_test);
  291. }
  292. #endif
  293. int s2 = (e2+1*m)%(3*m);
  294. int d2 = (e2+2*m)%(3*m);
  295. const int vs1 = V[s1];
  296. const int vd2 = V[d2];
  297. const int vd1 = V[d1];
  298. const int vs2 = V[s2];
  299. #ifdef IGL_SPLIT_NONMANIFOLD_DEBUG
  300. const auto simulated_merge_is_manifold_old = [&]()->bool
  301. {
  302. Eigen::VectorXi V_copy = V;
  303. std::vector<std::vector<int> > C_copy = C;
  304. merge_vertex(V_copy,C_copy,vs1,vd2);
  305. merge_vertex(V_copy,C_copy,vd1,vs2);
  306. Eigen::Map<Eigen::MatrixXi> CF_copy =
  307. Eigen::Map<Eigen::MatrixXi>(V_copy.data(),CF.rows(),CF.cols());
  308. if(!is_manifold(CF_copy)) { return false; }
  309. return true;
  310. };
  311. const bool ret_old = simulated_merge_is_manifold_old();
  312. const bool ret = simulated_merge_is_manifold(vs1,vd2,vd1,vs2);
  313. if(ret != ret_old)
  314. {
  315. assert(false);
  316. }
  317. #endif
  318. // I claim this is completely unnecessary if the unique edge was originally
  319. // manifold.
  320. //
  321. // I also hypothesize that this is unnecessary when conducting depth-first
  322. // traversals starting at a successful merge.
  323. //
  324. // That is, we never need to call this in the current algorithm.
  325. const int edge_valence = uEC(u+1)-uEC(u);
  326. assert(edge_valence >= 2);
  327. if(edge_valence>2 && !simulated_merge_is_manifold(vs1,vd2,vd1,vs2))
  328. {
  329. return false;
  330. }
  331. // Now we can merge
  332. merge_vertex(V,C,vs1,vd2);
  333. merge_vertex(V,C,vd1,vs2);
  334. return true;
  335. };
  336. // Consider each unique edge in the original mesh
  337. // number of faces incident on each unique edge
  338. Eigen::VectorXi D = uEC.tail(uEC.rows()-1)-uEC.head(uEC.rows()-1);
  339. Eigen::VectorXi uI;
  340. {
  341. Eigen::VectorXi sD;
  342. igl::sort(D,1,true,sD,uI);
  343. }
  344. const std::function<void(const int)> dfs = [&](const int e)
  345. {
  346. // we just successfully merged e, find all other non-manifold edges sharing
  347. // a current vertex with e and try to merge it too.
  348. const int s = (e+1*m)%(3*m);
  349. const int d = (e+2*m)%(3*m);
  350. for(const int c : {s,d})
  351. {
  352. const int v = V[c];
  353. std::vector<int> nme = nonmanifold_edge_star(v);
  354. // My thinking is that this must be size 0 or 2.
  355. //
  356. // But this seems very not true...
  357. for(int i = 0;i<nme.size();i++)
  358. {
  359. const int e1 = nme[i];
  360. for(int j = i+1;j<nme.size();j++)
  361. {
  362. const int e2 = nme[j];
  363. if(merge_edge(e1,e2))
  364. {
  365. dfs(e2);
  366. }
  367. }
  368. }
  369. }
  370. };
  371. // Every edge starts as a boundary
  372. for(auto u : uI)
  373. {
  374. // if boundary skip
  375. if(uEC(u+1)-uEC(u) == 1) { continue; }
  376. for(int j = uEC(u);j<uEC(u+1);j++)
  377. {
  378. const int e1 = uEE(j);
  379. for(int k = j+1;k<uEC(u+1);k++)
  380. {
  381. const int e2 = uEE(k);
  382. if(merge_edge(e1,e2))
  383. {
  384. // for non-manifold edges, launch search from e1 and e2
  385. if(uEC(u+1)-uEC(u) > 2)
  386. {
  387. dfs(e1);
  388. }
  389. break;
  390. }
  391. }
  392. }
  393. }
  394. // Ideally we'd do this so that all duplicated vertices end up at the end
  395. // rather than scrambling the whole mesh.
  396. {
  397. SVI.resize(F.size());
  398. std::vector<bool> marked(F.size());
  399. VectorXI J = VectorXI::Constant(F.size(),-1);
  400. SF.resize(F.rows(),F.cols());
  401. {
  402. int nv = 0;
  403. for(int f = 0;f<m;f++)
  404. {
  405. for(int i = 0;i<3;i++)
  406. {
  407. const int c = CF(f,i);
  408. if(J(c) == -1)
  409. {
  410. J(c) = nv;
  411. SVI(nv) = F(f,i);
  412. nv++;
  413. }
  414. SF(f,i) = J(c);
  415. }
  416. }
  417. SVI.conservativeResize(nv);
  418. }
  419. }
  420. }
  421. template <
  422. typename DerivedV,
  423. typename DerivedF,
  424. typename DerivedSV,
  425. typename DerivedSF,
  426. typename DerivedSVI
  427. >
  428. IGL_INLINE void igl::split_nonmanifold(
  429. const Eigen::MatrixBase<DerivedV> & V,
  430. const Eigen::MatrixBase<DerivedF> & F,
  431. Eigen::PlainObjectBase <DerivedSV> & SV,
  432. Eigen::PlainObjectBase <DerivedSF> & SF,
  433. Eigen::PlainObjectBase <DerivedSVI> & SVI)
  434. {
  435. igl::split_nonmanifold(F,SF,SVI);
  436. SV = V(SVI.derived(),Eigen::placeholders::all);
  437. }
  438. #ifdef IGL_STATIC_LIBRARY
  439. template void igl::split_nonmanifold<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::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<double, -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> >&);
  440. #endif