PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
Loading...
Searching...
No Matches
Mass.h File Reference

FEM mass matrix operator. More...

#include "Concepts.h"
#include "ShapeFunctions.h"
#include "pbat/Aliases.h"
#include "pbat/common/Eigen.h"
#include "pbat/profiling/Profiling.h"
#include <exception>
#include <fmt/core.h>
#include <tbb/parallel_for.h>

Go to the source code of this file.

Namespaces

namespace  pbat
 The main namespace of the library.
 
namespace  pbat::fem
 Finite Element Method (FEM)
 

Functions

template<typename TElement, typename TN>
auto pbat::fem::ElementMassMatrix (Eigen::MatrixBase< TN > const &N, typename TN::Scalar w, typename TN::Scalar rho)
 Compute the element mass matrix at a single quadrature point.
 
template<typename TElement, typename TN, typename TDerivedwg, typename TDerivedrhog, typename TDerivedOut>
void pbat::fem::ToElementMassMatrices (Eigen::MatrixBase< TN > const &Neg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::PlainObjectBase< TDerivedOut > &Meg)
 Compute a matrix of horizontally stacked element mass matrices for all quadrature points.
 
template<typename TElement, typename TN, typename TDerivedwg, typename TDerivedrhog>
auto pbat::fem::ElementMassMatrices (Eigen::MatrixBase< TN > const &Neg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog)
 Compute a matrix of horizontally stacked element mass matrices for all quadrature points.
 
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedMe, class TDerivedIn, class TDerivedOut>
void pbat::fem::GemmMass (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Me, int dims, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y)
 Compute mass matrix-matrix multiply \( Y += \mathbf{M} X \) using precomputed element mass matrices.
 
template<CMesh TMesh, class TDerivedeg, class TDerivedMe, class TDerivedIn, class TDerivedOut>
void pbat::fem::GemmMass (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Me, int dims, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y)
 Compute mass matrix-matrix multiply \( Y += \mathbf{M} X \) using mesh and precomputed element mass matrices.
 
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedrhog, class TDerivedNeg>
auto pbat::fem::MassMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::MatrixBase< TDerivedNeg > const &Neg, int dims=1) -> Eigen::SparseMatrix< typename TDerivedNeg::Scalar, Options, typename TDerivedE::Scalar >
 Construct the mass matrix operator's sparse matrix representation.
 
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg, class TDerivedMe>
auto pbat::fem::MassMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Meg, int dims=1) -> Eigen::SparseMatrix< typename TDerivedMe::Scalar, Options, typename TDerivedE::Scalar >
 Construct the mass matrix operator's sparse matrix representation from precomputed element mass matrices.
 
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedeg, class TDerivedwg, class TDerivedrhog, class TDerivedNeg>
auto pbat::fem::MassMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::MatrixBase< TDerivedNeg > const &Neg, int dims=1) -> Eigen::SparseMatrix< typename TDerivedNeg::Scalar, Options, typename TMesh::IndexType >
 Construct the mass matrix operator's sparse matrix representation using mesh.
 
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedeg, class TDerivedMe>
auto pbat::fem::MassMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Meg, int dims=1) -> Eigen::SparseMatrix< typename TDerivedMe::Scalar, Options, typename TMesh::IndexType >
 Construct the mass matrix operator's sparse matrix representation using mesh and precomputed element mass matrices.
 
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedrhog, class TDerivedNeg, class TDerivedOut>
void pbat::fem::ToLumpedMassMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::MatrixBase< TDerivedNeg > const &Neg, int dims, Eigen::PlainObjectBase< TDerivedOut > &m)
 Compute lumped mass matrix's diagonal vector into existing output vector.
 
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedMe, class TDerivedOut>
void pbat::fem::ToLumpedMassMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Meg, int dims, Eigen::PlainObjectBase< TDerivedOut > &m)
 Compute lumped mass vector from precomputed element mass matrices into existing output vector.
 
template<CMesh TMesh, class TDerivedeg, class TDerivedwg, class TDerivedrhog, class TDerivedNeg, class TDerivedOut>
void pbat::fem::ToLumpedMassMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::MatrixBase< TDerivedNeg > const &Neg, int dims, Eigen::PlainObjectBase< TDerivedOut > &m)
 Compute lumped mass vector using mesh and precomputed element mass matrices into existing output vector.
 
template<CMesh TMesh, class TDerivedeg, class TDerivedMe, class TDerivedOut>
void pbat::fem::ToLumpedMassMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Meg, int dims, Eigen::PlainObjectBase< TDerivedOut > &m)
 Compute lumped mass vector using mesh and precomputed element mass matrices into existing output vector.
 
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedrhog, class TDerivedNeg>
auto pbat::fem::LumpedMassMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::MatrixBase< TDerivedNeg > const &Neg, int dims=1) -> Eigen::Vector< typename TDerivedNeg::Scalar, Eigen::Dynamic >
 Compute lumped mass matrix's diagonal vector (allocates output vector)
 
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedMe>
auto pbat::fem::LumpedMassMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Meg, int dims=1) -> Eigen::Vector< typename TDerivedMe::Scalar, Eigen::Dynamic >
 Compute lumped mass vector from precomputed element mass matrices (allocates output vector)
 
template<CMesh TMesh, class TDerivedeg, class TDerivedwg, class TDerivedrhog, class TDerivedNeg>
auto pbat::fem::LumpedMassMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::DenseBase< TDerivedrhog > const &rhog, Eigen::MatrixBase< TDerivedNeg > const &Neg, int dims=1)
 Compute lumped mass vector using mesh (allocates output vector)
 
template<CMesh TMesh, class TDerivedeg, class TDerivedMe>
auto pbat::fem::LumpedMassMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedMe > const &Meg, int dims=1)
 Compute lumped mass vector using mesh and precomputed element mass matrices (allocates output vector)
 

Detailed Description

FEM mass matrix operator.

Author
Quoc-Minh Ton-That (tonth.nosp@m.at.q.nosp@m.uocmi.nosp@m.nh@g.nosp@m.mail..nosp@m.com)
Date
2025-02-11

The mass matrix \( \mathbf{M}_{ij} = \int_\Omega \rho(X) \phi_i(X) \phi_j(X) \) of a finite element discretized function under Galerkin projection.

This file provides functions to compute mass matrix related quantities.