trim_with_solid.cpp 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 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 "trim_with_solid.h"
  9. #include "assign.h"
  10. #include "intersect_other.h"
  11. #include "remesh_self_intersections.h"
  12. #include "point_solid_signed_squared_distance.h"
  13. #include "../../extract_manifold_patches.h"
  14. #include "../../connected_components.h"
  15. #include "../../facet_adjacency_matrix.h"
  16. #include "../../list_to_matrix.h"
  17. #include "../../find.h"
  18. #include "../../get_seconds.h"
  19. #include "../../barycenter.h"
  20. #include "../../remove_unreferenced.h"
  21. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  22. #include <vector>
  23. template <
  24. typename DerivedVA,
  25. typename DerivedFA,
  26. typename DerivedVB,
  27. typename DerivedFB,
  28. typename DerivedV,
  29. typename DerivedF,
  30. typename DerivedD,
  31. typename DerivedJ>
  32. IGL_INLINE void igl::copyleft::cgal::trim_with_solid(
  33. const Eigen::PlainObjectBase<DerivedVA> & VA,
  34. const Eigen::PlainObjectBase<DerivedFA> & FA,
  35. const Eigen::PlainObjectBase<DerivedVB> & VB,
  36. const Eigen::PlainObjectBase<DerivedFB> & FB,
  37. Eigen::PlainObjectBase<DerivedV> & Vd,
  38. Eigen::PlainObjectBase<DerivedF> & F,
  39. Eigen::PlainObjectBase<DerivedD> & D,
  40. Eigen::PlainObjectBase<DerivedJ> & J)
  41. {
  42. return trim_with_solid(
  43. VA,FA,VB,FB,TrimWithSolidMethod::CHECK_EACH_PATCH,
  44. Vd,F,D,J);
  45. }
  46. template <
  47. typename DerivedVA,
  48. typename DerivedFA,
  49. typename DerivedVB,
  50. typename DerivedFB,
  51. typename DerivedV,
  52. typename DerivedF,
  53. typename DerivedD,
  54. typename DerivedJ>
  55. IGL_INLINE void igl::copyleft::cgal::trim_with_solid(
  56. const Eigen::PlainObjectBase<DerivedVA> & VA,
  57. const Eigen::PlainObjectBase<DerivedFA> & FA,
  58. const Eigen::PlainObjectBase<DerivedVB> & VB,
  59. const Eigen::PlainObjectBase<DerivedFB> & FB,
  60. const TrimWithSolidMethod method,
  61. Eigen::PlainObjectBase<DerivedV> & Vd,
  62. Eigen::PlainObjectBase<DerivedF> & F,
  63. Eigen::PlainObjectBase<DerivedD> & D,
  64. Eigen::PlainObjectBase<DerivedJ> & J)
  65. {
  66. // resolve intersections using exact representation
  67. typedef Eigen::Matrix<CGAL::Epeck::FT,Eigen::Dynamic,3> MatrixX3E;
  68. typedef Eigen::Matrix<CGAL::Epeck::FT,Eigen::Dynamic,1> VectorXE;
  69. // Previously this used intersect_other to resolve intersections between A and
  70. // B hoping that it'd merge the two into a mesh where patches are separated by
  71. // non-manifold edges. This happens most of the time but sometimes the
  72. // new triangulations on faces of A don't match those on faces of B.
  73. // Specifically it seems you can get T junctions:
  74. //
  75. // /|\
  76. // / | \
  77. // / B| \
  78. // ---| A ⋅
  79. // \ B| /
  80. // \ | /
  81. // \|/
  82. // Probably intersect_other should not be attempting to output a single mesh
  83. // (i.e., when detect_only=false).
  84. //
  85. // # Alternative 1)
  86. //
  87. // Just call point_solid_signed_squared_distance for each output face.
  88. // Obviously O(#output-faces) calls to point_solid_signed_squared_distance.
  89. // But we get to use intersect_other to avoid finding and remeshing
  90. // self-intersections in A
  91. //
  92. // # Alternative 2)
  93. //
  94. // Use SelfIntersectMesh to really get a
  95. // single mesh with non-manifold edges. _But_ this would resolve any existing
  96. // self-intersections in (VA,FA) which is not requested. An idea to "undo"
  97. // this resolution is to undo any intersections _involving_ faces between A,B.
  98. //
  99. // This results in O(#patches) calls to point_solid_signed_squared_distance
  100. // but calls remeshes on the order of O(#self-intersections-in-A)
  101. //
  102. // # Alterative 3)
  103. //
  104. // Use intersect_other but then create an adjacency matrix based on facets
  105. // that share an edge but have a dissimilar J value. This will likely result
  106. // in lots of tiny patches along the intersection belt. So let's say it has
  107. // O(#A-B-intersection) calls to point_solid_signed_squared_distance
  108. //
  109. // If point_solid_signed_squared_distance turns out to me costly then 1) is
  110. // out and we should do 2) or 3).
  111. //
  112. // Indeed. point_solid_signed_squared_distance is a major bottleneck for some
  113. // examples.
  114. //
  115. MatrixX3E V;
  116. const auto set_D_via_patches =
  117. [&V,&F,&D,&VB,&FB](const int num_patches, const Eigen::VectorXi & P)
  118. {
  119. // Aggregate representative query points for each patch
  120. std::vector<bool> flag(num_patches);
  121. std::vector<std::vector<CGAL::Epeck::FT> > vQ;
  122. Eigen::VectorXi P2Q(num_patches);
  123. for(int f = 0;f<P.rows();f++)
  124. {
  125. const auto p = P(f);
  126. // if not yet processed this patch
  127. if(!flag[p])
  128. {
  129. P2Q(p) = vQ.size();
  130. std::vector<CGAL::Epeck::FT> q = {
  131. (V(F(f,0),0)+ V(F(f,1),0)+ V(F(f,2),0))/3.,
  132. (V(F(f,0),1)+ V(F(f,1),1)+ V(F(f,2),1))/3.,
  133. (V(F(f,0),2)+ V(F(f,1),2)+ V(F(f,2),2))/3.};
  134. vQ.emplace_back(q);
  135. flag[p] = true;
  136. }
  137. }
  138. MatrixX3E Q;
  139. igl::list_to_matrix(vQ,Q);
  140. VectorXE SP;
  141. point_solid_signed_squared_distance(Q,VB,FB,SP);
  142. Eigen::Matrix<bool,Eigen::Dynamic,1> DP = SP.array()>0;
  143. // distribute flag to all faces
  144. D.resize(F.rows());
  145. for(int f = 0;f<F.rows();f++)
  146. {
  147. D(f) = DP(P2Q(P(f)));
  148. }
  149. };
  150. switch(method)
  151. {
  152. case CHECK_EACH_PATCH:
  153. case CHECK_EACH_FACE:
  154. {
  155. Eigen::MatrixXi _1;
  156. Eigen::VectorXi _2;
  157. // Intersect A and B meshes and stitch together new faces
  158. igl::copyleft::cgal::intersect_other(
  159. VA,FA,VB,FB,{false,false,true},_1,V,F,J,_2);
  160. const auto keep = igl::find( (J.array()<FA.rows()).eval() );
  161. F = F(keep,Eigen::all).eval();
  162. J = J(keep).eval();
  163. {
  164. Eigen::VectorXi _;
  165. igl::remove_unreferenced(decltype(V)(V),decltype(F)(F),V,F,_);
  166. }
  167. switch(method)
  168. {
  169. default: /*unreachable*/ break;
  170. case CHECK_EACH_PATCH:
  171. {
  172. Eigen::SparseMatrix<bool> A;
  173. igl::facet_adjacency_matrix(F,A);
  174. for(int i = 0; i < A.outerSize(); i++)
  175. {
  176. for(decltype(A)::InnerIterator it(A,i); it; ++it)
  177. {
  178. const int a = it.row();
  179. const int b = it.col();
  180. if(J(a) == J(b))
  181. {
  182. A.coeffRef(a,b) = false;
  183. }
  184. }
  185. }
  186. A.prune(false);
  187. Eigen::VectorXi P,K;
  188. const int num_patches = igl::connected_components(A,P,K);
  189. set_D_via_patches(num_patches,P);
  190. break;
  191. }
  192. case CHECK_EACH_FACE:
  193. {
  194. MatrixX3E Q;
  195. igl::barycenter(V,F,Q);
  196. VectorXE SP;
  197. point_solid_signed_squared_distance(Q,VB,FB,SP);
  198. D = (SP.array()>0).template cast<typename DerivedD::Scalar>();
  199. break;
  200. }
  201. }
  202. break;
  203. }
  204. case RESOLVE_BOTH_AND_RESTORE_THEN_CHECK_EACH_PATCH:
  205. {
  206. RemeshSelfIntersectionsParam params;
  207. // This is somewhat dubious but appears to work. The stitch_all flag is
  208. // poorly documented.
  209. params.stitch_all = false;
  210. {
  211. Eigen::MatrixXi IF;
  212. Eigen::Matrix<typename DerivedVA::Scalar,Eigen::Dynamic,3> VAB(VA.rows() + VB.rows(),3);
  213. VAB << VA,VB;
  214. Eigen::Matrix<typename DerivedFA::Scalar,Eigen::Dynamic,3> FAB(FA.rows() + FB.rows(),3);
  215. FAB << FA,FB.array()+VA.rows();
  216. /// Sigh. Can't use this because of how it calls remove_unreferenced
  217. // remesh_self_intersections(VAB,FAB,params,V,F,IF,J);
  218. {
  219. Eigen::VectorXi I;
  220. igl::copyleft::cgal::remesh_self_intersections(
  221. VAB, FAB, params, V, F, IF, J, I);
  222. // Undo self-intersection remeshing in FA
  223. {
  224. Eigen::Array<bool,Eigen::Dynamic,1> avoids_B =
  225. Eigen::Array<bool,Eigen::Dynamic,1>::Constant(FA.rows(),1,true);
  226. for(int p = 0;p<IF.rows();p++)
  227. {
  228. if(IF(p,0) >= FA.rows() || IF(p,1) >= FA.rows())
  229. {
  230. if(IF(p,0) < FA.rows()){ avoids_B[IF(p,0)] = false; }
  231. if(IF(p,1) < FA.rows()){ avoids_B[IF(p,1)] = false; }
  232. }
  233. }
  234. // Find first entry for each in J
  235. Eigen::VectorXi first = Eigen::VectorXi::Constant(FA.rows(),1,-1);
  236. for(int j = 0;j<J.rows();j++)
  237. {
  238. if(J(j) < FA.rows() && first[J(j)] == -1)
  239. {
  240. first[J(j)] = j;
  241. // restore original face at this first entry
  242. if(avoids_B[J(j)])
  243. {
  244. F.row(j) = FA.row(J(j));
  245. }
  246. }
  247. }
  248. // Maybe this cannot happen for co-planar?
  249. assert(first.minCoeff() >= 0 && "Every face should be found");
  250. std::vector<int> keep;
  251. for(int f = 0;f<F.rows();f++)
  252. {
  253. if(J(f)>=FA.rows() || !avoids_B[J(f)] || first[J(f)] == f)
  254. {
  255. keep.push_back(f);
  256. }
  257. }
  258. F = F(keep,Eigen::all).eval();
  259. J = J(keep).eval();
  260. }
  261. // Don't do this until the very end:
  262. assert(I.size() == V.rows());
  263. // Merge coinciding vertices into non-manifold vertices.
  264. std::for_each(F.data(),F.data()+F.size(),[&I](typename DerivedF::Scalar & a){a=I[a];});
  265. }
  266. }
  267. // Partition result into manifold patches
  268. Eigen::VectorXi P;
  269. const int num_patches = igl::extract_manifold_patches(F,P);
  270. // only keep faces from A
  271. Eigen::Array<bool,Eigen::Dynamic,1> A = J.array()< FA.rows();
  272. const auto AI = igl::find(A);
  273. F = F(AI,Eigen::all).eval();
  274. J = J(AI).eval();
  275. P = P(AI).eval();
  276. set_D_via_patches(num_patches,P);
  277. break;
  278. }
  279. }
  280. {
  281. Eigen::VectorXi _;
  282. igl::remove_unreferenced(MatrixX3E(V),DerivedF(F),V,F,_);
  283. }
  284. assign(V,Vd);
  285. }
  286. #ifdef IGL_STATIC_LIBRARY
  287. // Explicit template instantiation
  288. template void igl::copyleft::cgal::trim_with_solid<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<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::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> > 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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  289. template void igl::copyleft::cgal::trim_with_solid<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<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Array<bool, -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<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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Array<bool, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  290. template void igl::copyleft::cgal::trim_with_solid<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<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Array<bool, -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<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::copyleft::cgal::TrimWithSolidMethod, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Array<bool, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  291. #endif