trim_with_solid.cpp 13 KB

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