// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2013 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 "cotmatrix.h" #include // For error printing #include #include "cotmatrix_entries.h" // Bug in unsupported/Eigen/SparseExtra needs iostream first #include template IGL_INLINE void igl::cotmatrix( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, Eigen::SparseMatrix& L) { using namespace Eigen; using namespace std; L.resize(V.rows(),V.rows()); Matrix edges; int simplex_size = F.cols(); // 3 for triangles, 4 for tets assert(simplex_size == 3 || simplex_size == 4); if(simplex_size == 3) { // This is important! it could decrease the comptuation time by a factor of 2 // Laplacian for a closed 2d manifold mesh will have on average 7 entries per // row L.reserve(10*V.rows()); edges.resize(3,2); edges << 1,2, 2,0, 0,1; }else if(simplex_size == 4) { L.reserve(17*V.rows()); edges.resize(6,2); edges << 1,2, 2,0, 0,1, 3,0, 3,1, 3,2; }else { return; } // Gather cotangents Matrix C; cotmatrix_entries(V,F,C); vector > IJV; IJV.reserve(F.rows()*edges.rows()*4); // Loop over triangles for(int i = 0; i < F.rows(); i++) { // loop over edges of element for(int e = 0;e(source,dest,C(i,e))); IJV.push_back(Triplet(dest,source,C(i,e))); IJV.push_back(Triplet(source,source,-C(i,e))); IJV.push_back(Triplet(dest,dest,-C(i,e))); } } L.setFromTriplets(IJV.begin(),IJV.end()); } #include "massmatrix.h" #include "cotmatrix_entries.h" #include "massmatrix.h" #include #include template < typename DerivedV, typename DerivedI, typename DerivedC, typename Scalar> IGL_INLINE void igl::cotmatrix( const Eigen::MatrixBase & V, const Eigen::MatrixBase & I, const Eigen::MatrixBase & C, Eigen::SparseMatrix& L, Eigen::SparseMatrix& M, Eigen::SparseMatrix& P) { typedef Eigen::Matrix RowVector3S; typedef Eigen::Matrix MatrixXS; typedef Eigen::Matrix VectorXS; typedef Eigen::Index Index; // number of vertices const Index n = V.rows(); // number of polyfaces const Index m = C.size()-1; assert(V.cols() == 2 || V.cols() == 3); std::vector > Lfijv; std::vector > Mfijv; std::vector > Pijv; // loop over vertices; set identity for original vertices for(Index i = 0;i X = decltype(X)::Zero(np+1,3); for(Index i = 0;i(A).solve(b); X.row(np) = w.transpose()*X.topRows(np); // scatter w into new row of P for(Index i = 0;i M; igl::massmatrix(X,F,igl::MASSMATRIX_TYPE_DEFAULT,M); Mp = M.diagonal(); } // Scatter into fine Laplacian and mass matrices const auto J = [&n,&np,&p,&I,&C](Index i)->Index{return i==np?n+p:I(C(p)+i);}; // Should just build Mf as a vector... for(Index i = 0;i Lf(n+m,n+m); Lf.setFromTriplets(Lfijv.begin(),Lfijv.end()); Eigen::SparseMatrix Mf(n+m,n+m); Mf.setFromTriplets(Mfijv.begin(),Mfijv.end()); L = P.transpose() * Lf * P; // "unlumped" M const Eigen::SparseMatrix PTMP = P.transpose() * Mf * P; // Lump M const VectorXS Mdiag = PTMP * VectorXS::Ones(n,1); M = Eigen::SparseMatrix(Mdiag.asDiagonal()); MatrixXS Vf = P*V; Eigen::MatrixXi Ff(I.size(),3); { Index f = 0; for(Index p = 0;p, Eigen::Matrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::SparseMatrix&, Eigen::SparseMatrix&, Eigen::SparseMatrix&); // generated by autoexplicit.sh template void igl::cotmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::SparseMatrix&); template void igl::cotmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::SparseMatrix&); template void igl::cotmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::SparseMatrix&); template void igl::cotmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::SparseMatrix&); #endif