// 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 "massmatrix.h" #include "massmatrix_intrinsic.h" #include "edge_lengths.h" #include "sparse.h" #include "doublearea.h" #include "volume.h" #include "voronoi_mass.h" #include "repmat.h" #include #include // For std::enable_if #include namespace igl { template struct MassMatrixHelper; // This would be easier with C++17 if constexpr // Specialization for triangles template struct MassMatrixHelper { static void compute( const Eigen::MatrixBase& V, const Eigen::MatrixBase& F, const MassMatrixType type, Eigen::SparseMatrix& M) { MassMatrixType eff_type = type == MASSMATRIX_TYPE_DEFAULT? MASSMATRIX_TYPE_VORONOI : type; Eigen::Matrix l; igl::edge_lengths(V, F, l); return massmatrix_intrinsic(l, F, eff_type, M); } }; // Specialization for tetrahedra template struct MassMatrixHelper { static void compute( const Eigen::MatrixBase& V, const Eigen::MatrixBase& F, const MassMatrixType type, Eigen::SparseMatrix& M) { const int n = V.rows(); const int m = F.rows(); using Eigen::Matrix; using Eigen::Dynamic; MassMatrixType eff_type = type == MASSMATRIX_TYPE_DEFAULT? MASSMATRIX_TYPE_BARYCENTRIC: type; Eigen::Matrix vol; volume(V, F, vol); vol = vol.array().abs(); Eigen::Matrix MI; Eigen::Matrix MJ; Eigen::Matrix MV; switch (eff_type) { case MASSMATRIX_TYPE_BARYCENTRIC: MI.resize(m*4,1); MJ.resize(m*4,1); MV.resize(m*4,1); MI.block(0*m,0,m,1) = F.col(0); MI.block(1*m,0,m,1) = F.col(1); MI.block(2*m,0,m,1) = F.col(2); MI.block(3*m,0,m,1) = F.col(3); MJ = MI; repmat(vol,4,1,MV); assert(MV.rows()==m*4&&MV.cols()==1); MV.array() /= 4.; break; case MASSMATRIX_TYPE_VORONOI: { MI = decltype(MI)::LinSpaced(n,0,n-1); MJ = MI; voronoi_mass(V,F,MV); break; } case MASSMATRIX_TYPE_FULL: MI.resize(m*16,1); MJ.resize(m*16,1); MV.resize(m*16,1); // indicies and values of the element mass matrix entries in the order // (1,0),(2,0),(3,0),(2,1),(3,1),(0,1),(3,2),(0,2),(1,2),(0,3),(1,3),(2,3),(0,0),(1,1),(2,2),(3,3); MI< struct MassMatrixHelper { static void compute( const Eigen::MatrixBase& V, const Eigen::MatrixBase& F, const MassMatrixType type, Eigen::SparseMatrix& M) { if (F.cols() == 3) { MassMatrixHelper::compute(V, F, type, M); } else if (F.cols() == 4) { MassMatrixHelper::compute(V, F, type, M); } else { // Handle unsupported simplex size at runtime assert(false && "Unsupported simplex size"); } } }; } template IGL_INLINE void igl::massmatrix( const Eigen::MatrixBase & V, const Eigen::MatrixBase & F, const MassMatrixType type, Eigen::SparseMatrix& M) { MassMatrixHelper::compute(V, F, type, M); } #ifdef IGL_STATIC_LIBRARY // Explicit template instantiation // generated by autoexplicit.sh template void igl::massmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, igl::MassMatrixType, Eigen::SparseMatrix&); // generated by autoexplicit.sh template void igl::massmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, igl::MassMatrixType, Eigen::SparseMatrix&); // generated by autoexplicit.sh template void igl::massmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, igl::MassMatrixType, Eigen::SparseMatrix&); // generated by autoexplicit.sh template void igl::massmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, igl::MassMatrixType, Eigen::SparseMatrix&); template void igl::massmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, igl::MassMatrixType, Eigen::SparseMatrix&); template void igl::massmatrix, Eigen::Matrix, double>(Eigen::MatrixBase > const&, Eigen::MatrixBase > const&, igl::MassMatrixType, Eigen::SparseMatrix&); #endif