// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2016 Alec Jacobson // // This Source Code Form is subject to the terms of the Mozilla Public License // v. 2.0. If a copy of the MPL was not distributed with this file, You can // obtain one at http://mozilla.org/MPL/2.0/. #include "trim_with_solid.h" #include "assign.h" #include "intersect_other.h" #include "remesh_self_intersections.h" #include "point_solid_signed_squared_distance.h" #include "../../extract_manifold_patches.h" #include "../../connected_components.h" #include "../../facet_adjacency_matrix.h" #include "../../list_to_matrix.h" #include "../../find.h" #include "../../get_seconds.h" #include "../../barycenter.h" #include "../../remove_unreferenced.h" #include #include template < typename DerivedVA, typename DerivedFA, typename DerivedVB, typename DerivedFB, typename DerivedV, typename DerivedF, typename DerivedD, typename DerivedJ> IGL_INLINE void igl::copyleft::cgal::trim_with_solid( const Eigen::PlainObjectBase & VA, const Eigen::PlainObjectBase & FA, const Eigen::PlainObjectBase & VB, const Eigen::PlainObjectBase & FB, Eigen::PlainObjectBase & Vd, Eigen::PlainObjectBase & F, Eigen::PlainObjectBase & D, Eigen::PlainObjectBase & J) { return trim_with_solid( VA,FA,VB,FB,TrimWithSolidMethod::CHECK_EACH_PATCH, Vd,F,D,J); } template < typename DerivedVA, typename DerivedFA, typename DerivedVB, typename DerivedFB, typename DerivedV, typename DerivedF, typename DerivedD, typename DerivedJ> IGL_INLINE void igl::copyleft::cgal::trim_with_solid( const Eigen::PlainObjectBase & VA, const Eigen::PlainObjectBase & FA, const Eigen::PlainObjectBase & VB, const Eigen::PlainObjectBase & FB, const TrimWithSolidMethod method, Eigen::PlainObjectBase & Vd, Eigen::PlainObjectBase & F, Eigen::PlainObjectBase & D, Eigen::PlainObjectBase & J) { // resolve intersections using exact representation typedef Eigen::Matrix MatrixX3E; typedef Eigen::Matrix VectorXE; // Previously this used intersect_other to resolve intersections between A and // B hoping that it'd merge the two into a mesh where patches are separated by // non-manifold edges. This happens most of the time but sometimes the // new triangulations on faces of A don't match those on faces of B. // Specifically it seems you can get T junctions: // // /|\ // / | \ // / B| \ // ---| A ⋅ // \ B| / // \ | / // \|/ // Probably intersect_other should not be attempting to output a single mesh // (i.e., when detect_only=false). // // # Alternative 1) // // Just call point_solid_signed_squared_distance for each output face. // Obviously O(#output-faces) calls to point_solid_signed_squared_distance. // But we get to use intersect_other to avoid finding and remeshing // self-intersections in A // // # Alternative 2) // // Use SelfIntersectMesh to really get a // single mesh with non-manifold edges. _But_ this would resolve any existing // self-intersections in (VA,FA) which is not requested. An idea to "undo" // this resolution is to undo any intersections _involving_ faces between A,B. // // This results in O(#patches) calls to point_solid_signed_squared_distance // but calls remeshes on the order of O(#self-intersections-in-A) // // # Alterative 3) // // Use intersect_other but then create an adjacency matrix based on facets // that share an edge but have a dissimilar J value. This will likely result // in lots of tiny patches along the intersection belt. So let's say it has // O(#A-B-intersection) calls to point_solid_signed_squared_distance // // If point_solid_signed_squared_distance turns out to me costly then 1) is // out and we should do 2) or 3). // // Indeed. point_solid_signed_squared_distance is a major bottleneck for some // examples. // MatrixX3E V; const auto set_D_via_patches = [&V,&F,&D,&VB,&FB](const int num_patches, const Eigen::VectorXi & P) { // Aggregate representative query points for each patch std::vector flag(num_patches); std::vector > vQ; Eigen::VectorXi P2Q(num_patches); for(int f = 0;f q = { (V(F(f,0),0)+ V(F(f,1),0)+ V(F(f,2),0))/3., (V(F(f,0),1)+ V(F(f,1),1)+ V(F(f,2),1))/3., (V(F(f,0),2)+ V(F(f,1),2)+ V(F(f,2),2))/3.}; vQ.emplace_back(q); flag[p] = true; } } MatrixX3E Q; igl::list_to_matrix(vQ,Q); VectorXE SP; point_solid_signed_squared_distance(Q,VB,FB,SP); Eigen::Matrix DP = SP.array()>0; // distribute flag to all faces D.resize(F.rows()); for(int f = 0;f A; igl::facet_adjacency_matrix(F,A); for(int i = 0; i < A.outerSize(); i++) { for(decltype(A)::InnerIterator it(A,i); it; ++it) { const int a = it.row(); const int b = it.col(); if(J(a) == J(b)) { A.coeffRef(a,b) = false; } } } A.prune(false); Eigen::VectorXi P,K; const int num_patches = igl::connected_components(A,P,K); set_D_via_patches(num_patches,P); break; } case CHECK_EACH_FACE: { MatrixX3E Q; igl::barycenter(V,F,Q); VectorXE SP; point_solid_signed_squared_distance(Q,VB,FB,SP); D = (SP.array()>0).template cast(); break; } } break; } case RESOLVE_BOTH_AND_RESTORE_THEN_CHECK_EACH_PATCH: { RemeshSelfIntersectionsParam params; // This is somewhat dubious but appears to work. The stitch_all flag is // poorly documented. params.stitch_all = false; { Eigen::MatrixXi IF; Eigen::Matrix VAB(VA.rows() + VB.rows(),3); VAB << VA,VB; Eigen::Matrix FAB(FA.rows() + FB.rows(),3); FAB << FA,FB.array()+VA.rows(); /// Sigh. Can't use this because of how it calls remove_unreferenced // remesh_self_intersections(VAB,FAB,params,V,F,IF,J); { Eigen::VectorXi I; igl::copyleft::cgal::remesh_self_intersections( VAB, FAB, params, V, F, IF, J, I); // Undo self-intersection remeshing in FA { Eigen::Array avoids_B = Eigen::Array::Constant(FA.rows(),1,true); for(int p = 0;p= FA.rows() || IF(p,1) >= FA.rows()) { if(IF(p,0) < FA.rows()){ avoids_B[IF(p,0)] = false; } if(IF(p,1) < FA.rows()){ avoids_B[IF(p,1)] = false; } } } // Find first entry for each in J Eigen::VectorXi first = Eigen::VectorXi::Constant(FA.rows(),1,-1); for(int j = 0;j= 0 && "Every face should be found"); std::vector keep; for(int f = 0;f=FA.rows() || !avoids_B[J(f)] || first[J(f)] == f) { keep.push_back(f); } } F = F(keep,Eigen::all).eval(); J = J(keep).eval(); } // Don't do this until the very end: assert(I.size() == V.rows()); // Merge coinciding vertices into non-manifold vertices. std::for_each(F.data(),F.data()+F.size(),[&I](typename DerivedF::Scalar & a){a=I[a];}); } } // Partition result into manifold patches Eigen::VectorXi P; const int num_patches = igl::extract_manifold_patches(F,P); // only keep faces from A Eigen::Array A = J.array()< FA.rows(); const auto AI = igl::find(A); F = F(AI,Eigen::all).eval(); J = J(AI).eval(); P = P(AI).eval(); set_D_via_patches(num_patches,P); break; } } { Eigen::VectorXi _; igl::remove_unreferenced(MatrixX3E(V),DerivedF(F),V,F,_); } assign(V,Vd); } #ifdef IGL_STATIC_LIBRARY // Explicit template instantiation template void igl::copyleft::cgal::trim_with_solid, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&); template void igl::copyleft::cgal::trim_with_solid, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Array, Eigen::Matrix >(Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&); template void igl::copyleft::cgal::trim_with_solid, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Array, Eigen::Matrix >(Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, igl::copyleft::cgal::TrimWithSolidMethod, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&); #endif