// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 2014 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 "harmonic.h" #include "adjacency_matrix.h" #include "cotmatrix.h" #include "diag.h" #include "invert_diag.h" #include "isdiag.h" #include "massmatrix.h" #include "min_quad_with_fixed.h" #include "speye.h" #include "sum.h" #include template < typename DerivedV, typename DerivedF, typename Derivedb, typename Derivedbc, typename DerivedW> IGL_INLINE bool igl::harmonic( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, const Eigen::MatrixBase & b, const Eigen::MatrixBase & bc, const int k, Eigen::PlainObjectBase & W) { using namespace Eigen; typedef typename DerivedV::Scalar Scalar; SparseMatrix L,M; cotmatrix(V,F,L); if(k>1) { massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M); } return harmonic(L,M,b,bc,k,W); } template < typename DerivedF, typename Derivedb, typename Derivedbc, typename DerivedW> IGL_INLINE bool igl::harmonic( const Eigen::MatrixBase & F, const Eigen::MatrixBase & b, const Eigen::MatrixBase & bc, const int k, Eigen::PlainObjectBase & W) { using namespace Eigen; typedef typename Derivedbc::Scalar Scalar; SparseMatrix A; adjacency_matrix(F,A); // sum each row SparseVector Asum; sum(A,1,Asum); // Convert row sums into diagonal of sparse matrix SparseMatrix Adiag; diag(Asum,Adiag); SparseMatrix L = A-Adiag; SparseMatrix M; speye(L.rows(),M); return harmonic(L,M,b,bc,k,W); } template < typename DerivedL, typename DerivedM, typename Derivedb, typename Derivedbc, typename DerivedW> IGL_INLINE bool igl::harmonic( const Eigen::SparseMatrix & L, const Eigen::SparseMatrix & M, const Eigen::MatrixBase & b, const Eigen::MatrixBase & bc, const int k, Eigen::PlainObjectBase & W) { const int n = L.rows(); assert(n == L.cols() && "L must be square"); assert((k==1 || n == M.cols() ) && "M must be same size as L"); assert((k==1 || n == M.rows() ) && "M must be square"); assert((k==1 || igl::isdiag(M)) && "Mass matrix should be diagonal"); Eigen::SparseMatrix Q; igl::harmonic(L,M,k,Q); typedef DerivedL Scalar; min_quad_with_fixed_data data; min_quad_with_fixed_precompute(Q,b,Eigen::SparseMatrix(),true,data); W.resize(n,bc.cols()); typedef Eigen::Matrix VectorXS; const VectorXS B = VectorXS::Zero(n,1); for(int w = 0;w IGL_INLINE void igl::harmonic( const Eigen::SparseMatrix & L, const Eigen::SparseMatrix & M, const int k, Eigen::SparseMatrix & Q) { assert(L.rows() == L.cols()&&"L should be square"); Q = -L; if(k == 1) return; assert(L.rows() == M.rows()&&"L should match M's dimensions"); assert(M.rows() == M.cols()&&"M should be square"); Eigen::SparseMatrix Mi; invert_diag(M,Mi); // This is **not** robust for k>2. See KKT system in [Jacobson et al. 2010] // of the kharmonic function in gptoolbox for(int p = 1;p IGL_INLINE void igl::harmonic( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, const int k, Eigen::SparseMatrix & Q) { Eigen::SparseMatrix L,M; cotmatrix(V,F,L); if(k>1) { massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M); } return harmonic(L,M,k,Q); } #ifdef IGL_STATIC_LIBRARY // Explicit template instantiation // generated by autoexplicit.sh template bool igl::harmonic, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::PlainObjectBase >&); // generated by autoexplicit.sh template bool igl::harmonic, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::PlainObjectBase >&); // generated by autoexplicit.sh template bool igl::harmonic, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::PlainObjectBase >&); // generated by autoexplicit.sh template void igl::harmonic, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::SparseMatrix&); // generated by autoexplicit.sh template bool igl::harmonic, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::PlainObjectBase >&); // generated by autoexplicit.sh template bool igl::harmonic, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::PlainObjectBase >&); template bool igl::harmonic, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix, Eigen::Matrix >(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, int, Eigen::PlainObjectBase >&); #endif