// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2023 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 "isolines_intrinsic.h" #include "edge_crossings.h" #include "cat.h" #include "unique_edge_map.h" #ifndef NDEBUG # include "is_edge_manifold.h" # include "is_vertex_manifold.h" #endif #include #include template < typename DerivedF, typename DerivedS, typename Derivedvals, typename DerivediB, typename DerivediFI, typename DerivediE, typename DerivedI> void igl::isolines_intrinsic( const Eigen::MatrixBase & F, const Eigen::MatrixBase & S, const Eigen::MatrixBase & vals, Eigen::PlainObjectBase & iB, Eigen::PlainObjectBase & iFI, Eigen::PlainObjectBase & iE, Eigen::PlainObjectBase & I) { Eigen::MatrixXi uE; Eigen::VectorXi EMAP,uEC,uEE; { Eigen::MatrixXi E; igl::unique_edge_map(F,E,uE,EMAP,uEC,uEE); } { std::vector viB(vals.size()); std::vector viFI(vals.size()); std::vector viE(vals.size()); std::vector vI(vals.size()); int num_vertices = 0; for(int j = 0;j void igl::isolines_intrinsic( const Eigen::MatrixBase & F, const Eigen::MatrixBase & S, const Eigen::MatrixBase & uE, const Eigen::MatrixBase & EMAP, const Eigen::MatrixBase & uEC, const Eigen::MatrixBase & uEE, const typename DerivedS::Scalar val, Eigen::PlainObjectBase & iB, Eigen::PlainObjectBase & iFI, Eigen::PlainObjectBase & iE) { using Scalar = typename DerivedS::Scalar; std::unordered_map uE2I; Eigen::Matrix T; igl::edge_crossings(uE,S,val,uE2I,T); iB.resize(uE2I.size(),F.cols()); iFI.resize(uE2I.size()); Eigen::VectorXi U(uE2I.size()); for(auto & pair : uE2I) { const int u = pair.first; const int w = pair.second; // first face incident on uE(u,:) const int e = uEE(uEC(u)); const int f = e % F.rows(); const int k = e / F.rows(); const bool flip = uE(u,0) != F(f,(k+1)%3); const double t = T(w); iB(w,k) = 0; iB(w,(k+1)%3) = flip? t:1-t; iB(w,(k+2)%3) = flip?1-t:t; iFI(w) = f; U(w) = u; } // Vertex crossings std::unordered_map V2I; { const auto add_vertex_crossing = [&iB,&iFI](const int k, const int i, const int j) { if(k >= iB.rows()) { iB.conservativeResize(2*iB.rows()+1,Eigen::NoChange); iFI.conservativeResize(2*iFI.rows()+1,Eigen::NoChange); } iFI(k) = i; iB.row(k) << 0,0,0; iB(k,j) = 1; }; int k = iB.rows(); for(int i = 0;i= iE.rows()) { iE.conservativeResize(2*iE.rows()+1,Eigen::NoChange); } iE.row(k) << i,j; }; { int r = 0; for(int f = 0;f < F.rows();f++) { // find first crossing edge int i; for(i = 0;i<3;i++) { if(uE2I.find(EMAP(f+F.rows()*i)) != uE2I.end()) { break; } } int j; for(j = i+1;j<3;j++) { if(uE2I.find(EMAP(f+F.rows()*j)) != uE2I.end()) { break; } } if(j<3) { // Connect two edge crossings. // other vertex const int k = 3-i-j; const int wi = uE2I[EMAP(f+F.rows()*i)]; const int wj = uE2I[EMAP(f+F.rows()*j)]; // flip orientation based on triangle gradient bool flip = S(F(f,k)) < val; flip = k%2? !flip:flip; if(flip) { set_row(r++,wi,wj); }else { set_row(r++,wj,wi); } }else if(i<3) { // The only valid vertex crossing is the opposite vertex const int v = F(f,i); // Is it a crossing? assert(V2I.find(v) != V2I.end()); //if(V2I.find(v) != V2I.end()) { const int wv = V2I[v]; const int wi = uE2I[EMAP(f+F.rows()*i)]; const bool flip = S(F(f,(i+1)%3)) > val; if(flip) { set_row(r++,wi,wv); }else { set_row(r++,wv,wi); } } }else { // Could have 2 vertex crossings. We're only interested if there're exactly two and if the other vertex is "above". int i = 0; for(i = 0;i<3;i++) { if(S(F(f,i)) == val) { break; } } int j; for(j = i+1;j<3;j++) { if(S(F(f,j)) == val) { break; } } if(j<3) { // check if the third is a crossing. const int k = 3-i-j; // Triangle is constant on the val. Skip. if(S(F(f,k)) == val){ continue; } // Is this a boundary edge? const int u = EMAP(f+F.rows()*k); const int count = uEC(u+1)-uEC(u); if( count == 1 || S(F(f,k)) > val) { const int wi = V2I[F(f,i)]; const int wj = V2I[F(f,j)]; bool flip = S(F(f,k)) < val; flip = k%2 ? !flip:flip; if(flip) { set_row(r++,wj,wi); }else { set_row(r++,wi,wj); } } } } } iE.conservativeResize(r,Eigen::NoChange); } #ifndef NDEBUG if(igl::is_vertex_manifold(F) && igl::is_edge_manifold(F)) { // Check that every vertex has one in one out Eigen::VectorXi in_count = Eigen::VectorXi::Zero(iB.rows()); Eigen::VectorXi out_count = Eigen::VectorXi::Zero(iB.rows()); for(int e = 0;e, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&, Eigen::PlainObjectBase >&); #endif