// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2013 Olga Diamanti, 2015 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/. #ifndef IGL_CONJUGATE_FF_SOLVER_DATA_H #define IGL_CONJUGATE_FF_SOLVER_DATA_H #include "igl_inline.h" #include #include namespace igl { // Data class for the Conjugate Frame Field Solver template class ConjugateFFSolverData { public: const Eigen::PlainObjectBase &V; int numV; const Eigen::PlainObjectBase &F; int numF; Eigen::MatrixXi EV; int numE; Eigen::MatrixXi F2E; Eigen::MatrixXi E2F; Eigen::VectorXd K; Eigen::VectorXi isBorderEdge; int numInteriorEdges; Eigen::Matrix E2F_int; Eigen::VectorXi indInteriorToFull; Eigen::VectorXi indFullToInterior; Eigen::PlainObjectBase B1, B2, FN; Eigen::Matrix kmin, kmax; Eigen::Matrix dmin, dmax; Eigen::Matrix dmin3, dmax3; Eigen::VectorXd nonPlanarityMeasure; Eigen::SparseMatrix > planarityWeight; //conjugacy matrix std::vector > H; //conjugacy matrix eigenvectors and (scaled) eigenvalues std::vector > UH; std::vector > s; //laplacians Eigen::SparseMatrix> DDA, DDB; private: IGL_INLINE void computeCurvatureAndPrincipals(); IGL_INLINE void precomputeConjugacyStuff(); IGL_INLINE void computeLaplacians(); IGL_INLINE void computek(); IGL_INLINE void computeCoefficientLaplacian(int n, Eigen::SparseMatrix > &D); IGL_INLINE void precomputeInteriorEdges(); public: IGL_INLINE ConjugateFFSolverData(const Eigen::PlainObjectBase &_V, const Eigen::PlainObjectBase &_F); IGL_INLINE void evaluateConjugacy(const Eigen::Matrix &pvU, const Eigen::Matrix &pvV, Eigen::Matrix &conjValues) const ; }; } #include #include #include #include #include #include template IGL_INLINE igl::ConjugateFFSolverData:: ConjugateFFSolverData(const Eigen::PlainObjectBase &_V, const Eigen::PlainObjectBase &_F): V(_V), numV(_V.rows()), F(_F), numF(_F.rows()) { igl::edge_topology(V,F,EV,F2E,E2F); numE = EV.rows(); precomputeInteriorEdges(); igl::local_basis(V,F,B1,B2,FN); computek(); computeLaplacians(); computeCurvatureAndPrincipals(); precomputeConjugacyStuff(); }; template IGL_INLINE void igl::ConjugateFFSolverData::computeCurvatureAndPrincipals() { Eigen::MatrixXd VCBary; Eigen::MatrixXi FCBary; VCBary.setZero(numV+numF,3); FCBary.setZero(3*numF,3); igl::false_barycentric_subdivision(V, F, VCBary, FCBary); Eigen::MatrixXd dmax3_,dmin3_; igl::principal_curvature(VCBary, FCBary, dmax3_, dmin3_, kmax, kmin, 5,true); dmax3 = dmax3_.bottomRows(numF); dmin3 = dmin3_.bottomRows(numF); kmax = kmax.bottomRows(numF); kmin = kmin.bottomRows(numF); // kmax = dmax3.rowwise().norm(); // kmin = dmin3.rowwise().norm(); dmin3.rowwise().normalize(); dmax3.rowwise().normalize(); dmax.setZero(numF,2); dmin.setZero(numF,2); for (int i= 0; i (0, numF-1); igl::sparse(I, I, nonPlanarityMeasure, numF, numF, planarityWeight); } template IGL_INLINE void igl::ConjugateFFSolverData::precomputeConjugacyStuff() { H.resize(numF); UH.resize(numF); s.resize(numF); for (int i = 0; i Ht = H[i].transpose(); H[i] = .5*(H[i]+Ht); Eigen::EigenSolver > es(H[i]); s[i] = es.eigenvalues().real();//ok to do this because H symmetric //scale s[i] = s[i]/(s[i].cwiseAbs().minCoeff()); UH[i] = es.eigenvectors().real(); } } template IGL_INLINE void igl::ConjugateFFSolverData::computeLaplacians() { computeCoefficientLaplacian(2, DDA); computeCoefficientLaplacian(4, DDB); } template IGL_INLINE void igl::ConjugateFFSolverData:: precomputeInteriorEdges() { // Flag border edges numInteriorEdges = 0; isBorderEdge.setZero(numE,1); indFullToInterior = -1*Eigen::VectorXi::Ones(numE,1); for(unsigned i=0; i IGL_INLINE void igl::ConjugateFFSolverData:: computeCoefficientLaplacian(int n, Eigen::SparseMatrix > &D) { std::vector >> tripletList; // For every non-border edge for (unsigned eid=0; eid >(fid0, fid0, std::complex(1.))); tripletList.push_back(Eigen::Triplet >(fid1, fid1, std::complex(1.))); tripletList.push_back(Eigen::Triplet >(fid0, fid1, -1.*std::polar(1.,-1.*n*K[eid]))); tripletList.push_back(Eigen::Triplet >(fid1, fid0, -1.*std::polar(1.,1.*n*K[eid]))); } } D.resize(numF,numF); D.setFromTriplets(tripletList.begin(), tripletList.end()); } template IGL_INLINE void igl::ConjugateFFSolverData:: computek() { K.setZero(numE); // For every non-border edge for (unsigned eid=0; eid N0 = FN.row(fid0); Eigen::Matrix N1 = FN.row(fid1); // find common edge on triangle 0 and 1 int fid0_vc = -1; int fid1_vc = -1; for (unsigned i=0;i<3;++i) { if (F2E(fid0,i) == eid) fid0_vc = i; if (F2E(fid1,i) == eid) fid1_vc = i; } assert(fid0_vc != -1); assert(fid1_vc != -1); Eigen::Matrix common_edge = V.row(F(fid0,(fid0_vc+1)%3)) - V.row(F(fid0,fid0_vc)); common_edge.normalize(); // Map the two triangles in a new space where the common edge is the x axis and the N0 the z axis Eigen::Matrix P; Eigen::Matrix o = V.row(F(fid0,fid0_vc)); Eigen::Matrix tmp = -N0.cross(common_edge); P << common_edge, tmp, N0; // P.transposeInPlace(); Eigen::Matrix V0; V0.row(0) = V.row(F(fid0,0)) -o; V0.row(1) = V.row(F(fid0,1)) -o; V0.row(2) = V.row(F(fid0,2)) -o; V0 = (P*V0.transpose()).transpose(); Eigen::Matrix V1; V1.row(0) = V.row(F(fid1,0)) -o; V1.row(1) = V.row(F(fid1,1)) -o; V1.row(2) = V.row(F(fid1,2)) -o; V1 = (P*V1.transpose()).transpose(); // compute rotation R such that R * N1 = N0 // i.e. map both triangles to the same plane double alpha = -atan2(V1((fid1_vc+2)%3,2),V1((fid1_vc+2)%3,1)); Eigen::Matrix R; R << 1, 0, 0, 0, cos(alpha), -sin(alpha) , 0, sin(alpha), cos(alpha); V1 = (R*V1.transpose()).transpose(); // measure the angle between the reference frames // k_ij is the angle between the triangle on the left and the one on the right Eigen::Matrix ref0 = V0.row(1) - V0.row(0); Eigen::Matrix ref1 = V1.row(1) - V1.row(0); ref0.normalize(); ref1.normalize(); double ktemp = atan2(ref1(1),ref1(0)) - atan2(ref0(1),ref0(0)); // just to be sure, rotate ref0 using angle ktemp... Eigen::Matrix R2; R2 << cos(ktemp), -sin(ktemp), sin(ktemp), cos(ktemp); Eigen::Matrix tmp1 = R2*(ref0.head(2)).transpose(); K[eid] = ktemp; } } } template IGL_INLINE void igl::ConjugateFFSolverData:: evaluateConjugacy(const Eigen::Matrix &pvU, const Eigen::Matrix &pvV, Eigen::Matrix &conjValues) const { conjValues.resize(numF,1); for (int j =0; j x; x<