split_nonmanifold.cpp 13 KB

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