// This file is part of libigl, a simple c++ geometry processing library. // // Copyright (C) 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/. #include "biharmonic_coordinates.h" #include "cotmatrix.h" #include "massmatrix.h" #include "min_quad_with_fixed.h" #include "normal_derivative.h" #include "on_boundary.h" #include template < typename DerivedV, typename DerivedT, typename SType, typename DerivedW> IGL_INLINE bool igl::biharmonic_coordinates( const Eigen::PlainObjectBase & V, const Eigen::PlainObjectBase & T, const std::vector > & S, Eigen::PlainObjectBase & W) { return biharmonic_coordinates(V,T,S,2,W); } template < typename DerivedV, typename DerivedT, typename SType, typename DerivedW> IGL_INLINE bool igl::biharmonic_coordinates( const Eigen::PlainObjectBase & V, const Eigen::PlainObjectBase & T, const std::vector > & S, const int k, Eigen::PlainObjectBase & W) { using namespace Eigen; using namespace std; // This is not the most efficient way to build A, but follows "Linear // Subspace Design for Real-Time Shape Deformation" [Wang et al. 2015]. SparseMatrix A; { SparseMatrix N,Z,L,K,M; normal_derivative(V,T,N); Array I; Array C; on_boundary(T,I,C); { std::vector >ZIJV; for(int t =0;t Minv = ((VectorXd)M.diagonal().array().inverse()).asDiagonal(); switch(k) { default: assert(false && "unsupported"); case 2: // For C1 smoothness in 2D, one should use bi-harmonic A = K.transpose() * (Minv * K); break; case 3: // For C1 smoothness in 3D, one should use tri-harmonic A = K.transpose() * (Minv * (-L * (Minv * K))); break; } } // Vertices in point handles const size_t mp = count_if(S.begin(),S.end(),[](const vector & h){return h.size()==1;}); // number of region handles const size_t r = S.size()-mp; // Vertices in region handles size_t mr = 0; for(const auto & h : S) { if(h.size() > 1) { mr += h.size(); } } const size_t dim = T.cols()-1; // Might as well be dense... I think... MatrixXd J = MatrixXd::Zero(mp+mr,mp+r*(dim+1)); VectorXi b(mp+mr); MatrixXd H(mp+r*(dim+1),dim); { int v = 0; int c = 0; for(int h = 0;h= dim+1); for(int p = 0;p(),VectorXd(),true,W); } #ifdef IGL_STATIC_LIBRARY // Explicit template specialization template bool igl::biharmonic_coordinates, Eigen::Matrix, int, Eigen::Matrix >(Eigen::PlainObjectBase > const&, Eigen::PlainObjectBase > const&, std::vector >, std::allocator > > > const&, int, Eigen::PlainObjectBase >&); #endif