// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2017 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 "straighten_seams.h" #include "LinSpaced.h" #include "on_boundary.h" #include "sparse.h" #include "max.h" #include "count.h" #include "any.h" #include "slice_mask.h" #include "slice_into.h" #include "unique_simplices.h" #include "adjacency_matrix.h" #include "setxor.h" #include "edges_to_path.h" #include "ramer_douglas_peucker.h" #include "vertex_components.h" #include "list_to_matrix.h" #include "ears.h" #include "slice.h" #include "sum.h" #include "find.h" #include template < typename DerivedV, typename DerivedF, typename DerivedVT, typename DerivedFT, typename Scalar, typename DerivedUE, typename DerivedUT, typename DerivedOT> IGL_INLINE void igl::straighten_seams( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, const Eigen::MatrixBase & VT, const Eigen::MatrixBase & FT, const Scalar tol, Eigen::PlainObjectBase & UE, Eigen::PlainObjectBase & UT, Eigen::PlainObjectBase & OT) { using namespace Eigen; // number of faces assert(FT.rows() == F.rows() && "#FT must == #F"); assert(F.cols() == 3 && "F should contain triangles"); assert(FT.cols() == 3 && "FT should contain triangles"); const int m = F.rows(); // Boundary edges of the texture map and 3d meshes Array _; Array BT,BF; on_boundary(FT,_,BT); on_boundary(F,_,BF); assert((!((BF && (BT!=true)).any())) && "Not dealing with boundaries of mesh that get 'stitched' in texture mesh"); typedef Matrix MatrixX2I; const MatrixX2I ET = (MatrixX2I(FT.rows()*3,2) <vBT = Map >(BT.data(),BT.size(),1); ArrayvBF = Map >(BF.data(),BF.size(),1); MatrixX2I OF; slice_mask(ET,vBT,1,OT); slice_mask(EF,vBT,1,OF); VectorXi OFMAP; slice_mask(EFMAP,vBT,1,OFMAP); // Two boundary edges on the texture-mapping are "equivalent" to each other on // the 3D-mesh if their 3D-mesh vertex indices match SparseMatrix OEQ; { SparseMatrix OEQR; sparse( igl::LinSpaced(OT.rows(),0,OT.rows()-1), OFMAP, Array::Ones(OT.rows(),1), OT.rows(), m*3, OEQR); OEQ = OEQR * OEQR.transpose(); // Remove diagonal OEQ.prune([](const int r, const int c, const bool)->bool{return r!=c;}); } // For each edge in OT, for each endpoint, how many _other_ texture-vertices // are images of all the 3d-mesh vertices in F who map from "corners" in F/FT // mapping to this endpoint. // // Adjacency matrix between 3d-vertices and texture-vertices SparseMatrix V2VT; sparse( F, FT, Array::Ones(F.rows(),F.cols()), V.rows(), VT.rows(), V2VT); // For each 3d-vertex count how many different texture-coordinates its getting // from different incident corners VectorXi DV; count(V2VT,2,DV); VectorXi M,I; max(V2VT,1,M,I); assert( (M.array() == 1).all() ); VectorXi DT; // Map counts onto texture-vertices slice(DV,I,1,DT); // Boundary in 3D && UV Array BTF; slice_mask(vBF, vBT, 1, BTF); // Texture-vertex is "sharp" if incident on "half-"edge that is not a // boundary in the 3D mesh but is a boundary in the texture-mesh AND is not // "cut cleanly" (the vertex is mapped to exactly 2 locations) Array SV = Array::Zero(VT.rows(),1); //std::cout<<"#SV: "< CL = DT.array()==2; SparseMatrix VTOT; { Eigen::MatrixXi I = igl::LinSpaced(OT.rows(),0,OT.rows()-1).replicate(1,2); sparse( OT, I, Array::Ones(OT.rows(),OT.cols()), VT.rows(), OT.rows(), VTOT); Array cuts; count( (VTOT*OEQ).eval(), 2, cuts); CL = (CL && (cuts.array() == 2)).eval(); } //std::cout<<"#CL: "< earT = Array::Zero(VT.rows(),1); for(int e = 0;e A; adjacency_matrix(FT,A); earT = (earT || (A*earT.matrix()).array()).eval(); //std::cout<<"#earT: "< V2VTSV,V2VTC; slice_mask(V2VT,SV,2,V2VTSV); Array Cb; any(V2VTSV,2,Cb); slice_mask(V2VT,Cb,1,V2VTC); any(V2VTC,1,SV); } //std::cout<<"#SV: "< OTVT = VTOT.transpose(); int nc; ArrayXi C; { // Doesn't Compile on older Eigen: //SparseMatrix A = OTVT * (!SV).matrix().asDiagonal() * VTOT; SparseMatrix A = OTVT * (SV!=true).matrix().asDiagonal() * VTOT; vertex_components(A,C); nc = C.maxCoeff()+1; } //std::cout<<"nc: "< > vUE; // loop over each component std::vector done(nc,false); for(int c = 0;c OEQIc; slice(OEQ,Ic,1,OEQIc); Eigen::VectorXi N; sum(OEQIc,2,N); const int ncopies = N(0)+1; assert((N.array() == ncopies-1).all()); assert((ncopies == 1 || ncopies == 2) && "Not dealing with non-manifold meshes"); Eigen::VectorXi vpath,epath,eend; typedef Eigen::Matrix MatrixX2S; switch(ncopies) { case 1: { MatrixX2I OTIc; slice(OT,Ic,1,OTIc); edges_to_path(OTIc,vpath,epath,eend); Array SVvpath; slice(SV,vpath,1,SVvpath); assert( (vpath(0) != vpath(vpath.size()-1) || !SVvpath.any()) && "Not dealing with 1-loops touching 'sharp' corners"); // simple open boundary MatrixX2S PI; slice(VT,vpath,1,PI); const Scalar bbd = (PI.colwise().maxCoeff() - PI.colwise().minCoeff()).norm(); // Do not collapse boundaries to fewer than 3 vertices const bool allow_boundary_collapse = false; assert(PI.size() >= 2); const bool is_closed = PI(0) == PI(PI.size()-1); assert(!is_closed || vpath.size() >= 4); Scalar eff_tol = std::min(tol,2.); VectorXi UIc; while(true) { MatrixX2S UPI,UTvpath; ramer_douglas_peucker(PI,eff_tol*bbd,UPI,UIc,UTvpath); slice_into(UTvpath,vpath,1,UT); if(!is_closed || allow_boundary_collapse) { break; } if(UPI.rows()>=4) { break; } eff_tol = eff_tol*0.5; } for(int i = 0;i IV; SparseMatrix OEQIcT = OEQIc.transpose().eval(); find(OEQIcT,Icc,II,IV); assert(II.size() == Ic.size() && (II.array() == igl::LinSpaced(Ic.size(),0,Ic.size()-1).array()).all()); assert(Icc.size() == Ic.size()); const int cc = C(Icc(0)); Eigen::VectorXi CIcc; slice(C,Icc,1,CIcc); assert((CIcc.array() == cc).all()); assert(!done[cc]); done[cc] = true; } Array flipped; { MatrixX2I OFIc,OFIcc; slice(OF,Ic,1,OFIc); slice(OF,Icc,1,OFIcc); Eigen::VectorXi XOR,IA,IB; setxor(OFIc,OFIcc,XOR,IA,IB); assert(XOR.size() == 0); flipped = OFIc.array().col(0) != OFIcc.array().col(0); } if(Ic.size() == 1) { // No change to UT vUE.push_back({OT(Ic(0),0),OT(Ic(0),1)}); assert(Icc.size() == 1); vUE.push_back({OT(Icc(0),flipped(0)?1:0),OT(Icc(0),flipped(0)?0:1)}); }else { MatrixX2I OTIc; slice(OT,Ic,1,OTIc); edges_to_path(OTIc,vpath,epath,eend); // Flip endpoints if needed for(int e = 0;e PI(vpath.size(),VT.cols()*2); for(int p = 0;p UPI,SI; VectorXi UIc; ramer_douglas_peucker(PI,tol*bbd,UPI,UIc,SI); slice_into(SI.leftCols (VT.cols()), vpath,1,UT); slice_into(SI.rightCols(VT.cols()),vpathc,1,UT); for(int i = 0;i, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, double, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, double, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&); #endif