PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
|
Finite Element Method (FEM) More...
Classes | |
struct | Mesh |
A generic stateful finite element mesh representation. More... | |
Concepts | |
concept | CElement |
Reference finite element. | |
concept | CMesh |
Finite element mesh. | |
Typedefs | |
template<int Order> | |
using | Hexahedron = typename detail::Hexahedron<Order> |
Hexahedron finite element. | |
template<int Order> | |
using | Line = typename detail::Line<Order> |
Line finite element. | |
template<int Order> | |
using | Quadrilateral = typename detail::Quadrilateral<Order> |
Quadrilateral finite element. | |
template<int Order> | |
using | Tetrahedron = typename detail::Tetrahedron<Order> |
Tetrahedron finite element. | |
template<int Order> | |
using | Triangle = typename detail::Triangle<Order> |
Triangle finite element. | |
Enumerations | |
enum class | EHyperElasticSpdCorrection : std::uint32_t { None , Projection , Absolute } |
Bit-flag enum for SPD projection type. | |
enum | EElementElasticityComputationFlags : int { Potential = 1 << 0 , Gradient = 1 << 1 , Hessian = 1 << 2 } |
Bit-flag enum for element elasticity computation flags. | |
Functions | |
template<CElement TElement, class TDerivedx, class TDerivedX> | |
auto | DeformationGradient (Eigen::MatrixBase< TDerivedx > const &x, Eigen::MatrixBase< TDerivedX > const &GP) -> Matrix< TDerivedx::RowsAtCompileTime, TElement::kDims > |
Computes the deformation gradient \( \frac{\partial \mathbf{x}(X)}{\partial X} \) of the deformation map \( \mathbf{x}(X) \). | |
template<CElement TElement, int Dims, math::linalg::mini::CMatrix TMatrixGF, math::linalg::mini::CMatrix TMatrixGP, class ScalarType = typename TMatrixGF::ScalarType> | |
auto | GradientSegmentWrtDofs (TMatrixGF const &GF, TMatrixGP const &GP, auto i) -> math::linalg::mini::SVector< ScalarType, Dims > |
Computes \( \frac{\partial \Psi}{\partial \mathbf{x}_i} \in \mathbb{R}^d \), i.e. the gradient of a scalar function \( \Psi \) w.r.t. the \( i^{\text{th}}
\) node's degrees of freedom. | |
template<CElement TElement, int Dims, math::linalg::mini::CMatrix TMatrixGF, math::linalg::mini::CMatrix TMatrixGP, class ScalarType = typename TMatrixGF::ScalarType> | |
auto | GradientWrtDofs (TMatrixGF const &GF, TMatrixGP const &GP) -> math::linalg::mini::SVector< ScalarType, TElement::kNodes *Dims > |
Computes gradient w.r.t. FEM degrees of freedom \( x \) of scalar function \(\Psi(\mathbf{F}) \), where \( F \) is the jacobian of \( x \), via chain rule. This is effectively a rank-3 to rank-1 tensor contraction. | |
template<CElement TElement, int Dims, math::linalg::mini::CMatrix TMatrixHF, math::linalg::mini::CMatrix TMatrixGP, class ScalarType = typename TMatrixHF::ScalarType> | |
auto | HessianBlockWrtDofs (TMatrixHF const &HF, TMatrixGP const &GP, auto i, auto j) -> math::linalg::mini::SMatrix< ScalarType, Dims, Dims > |
Computes \( \frac{\partial^2 \Psi}{\partial \mathbf{x}_i \partial \mathbf{x}_j} \), i.e. the hessian of a scalar function \( \Psi \) w.r.t. the \( i^{\text{th}} \) and \(j^{\text{th}} \) node's degrees of freedom. | |
template<CElement TElement, int Dims, math::linalg::mini::CMatrix TMatrixHF, math::linalg::mini::CMatrix TMatrixGP, class ScalarType = typename TMatrixHF::ScalarType> | |
auto | HessianWrtDofs (TMatrixHF const &HF, TMatrixGP const &GP) -> math::linalg::mini::SMatrix< ScalarType, TElement::kNodes *Dims, TElement::kNodes *Dims > |
Computes hessian w.r.t. FEM degrees of freedom \( x \) of scalar function \(\Psi(\mathbf{F}) \), where \( F \) is the jacobian of \( x \), via chain rule. This is effectively a rank-4 to rank-2 tensor contraction. | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedGNeg, class TDerivedIn, class TDerivedOut> | |
void | GemmGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Compute gradient matrix-vector multiply \( Y += \mathbf{G} X \). | |
template<CMesh TMesh, class TDerivedeg, class TDerivedGNeg, class TDerivedIn, class TDerivedOut> | |
void | GemmGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Compute gradient matrix-vector multiply \( Y += \mathbf{G} X \) using mesh. | |
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg, class TDerivedGNeg> | |
auto | GradientMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg) -> Eigen::SparseMatrix< typename TDerivedGNeg::Scalar, Options, typename TDerivedE::Scalar > |
Construct the gradient operator's sparse matrix representation. | |
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedeg, class TDerivedGNeg> | |
auto | GradientMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg) -> Eigen::SparseMatrix< typename TDerivedGNeg::Scalar, Options, typename TMesh::IndexType > |
Construct the gradient operator's sparse matrix representation using mesh. | |
math::linalg::EEigenvalueFilter | ToEigenvalueFilter (EHyperElasticSpdCorrection mode) |
Convert EHyperElasticSpdCorrection to EEigenvalueFilter. | |
template<CElement TElement, int Dims, physics::CHyperElasticEnergy THyperElasticEnergy, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx, class TDerivedUg, class TDerivedGg, class TDerivedHg> | |
void | ToElementElasticity (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x, Eigen::PlainObjectBase< TDerivedUg > &Ug, Eigen::PlainObjectBase< TDerivedGg > &Gg, Eigen::PlainObjectBase< TDerivedHg > &Hg, int eFlags=EElementElasticityComputationFlags::Potential, EHyperElasticSpdCorrection eSpdCorrection=EHyperElasticSpdCorrection::None) |
Compute element elasticity and its derivatives at the given shape. | |
template<physics::CHyperElasticEnergy THyperElasticEnergy, CMesh TMesh, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx, class TDerivedUg, class TDerivedGg, class TDerivedHg> | |
void | ToElementElasticity (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x, Eigen::PlainObjectBase< TDerivedUg > &Ug, Eigen::PlainObjectBase< TDerivedGg > &Gg, Eigen::PlainObjectBase< TDerivedHg > &Hg, int eFlags=EElementElasticityComputationFlags::Potential, EHyperElasticSpdCorrection eSpdCorrection=EHyperElasticSpdCorrection::None) |
Compute element elasticity using mesh. | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedHg, class TDerivedIn, class TDerivedOut> | |
void | GemmHyperElastic (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedHg > const &Hg, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Apply the hessian matrix as a linear operator \( Y += \mathbf{H} X \). | |
template<CMesh TMesh, class TDerivedeg, class TDerivedHg, class TDerivedIn, class TDerivedOut> | |
void | GemmHyperElastic (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedHg > const &Hg, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Apply the hessian matrix as a linear operator \( Y += \mathbf{H} X \) using mesh. | |
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg> | |
auto | ElasticHessianSparsity (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg) -> math::linalg::SparsityPattern< typename TDerivedE::Scalar, Options > |
Construct the hessian matrix's sparsity pattern. | |
template<Eigen::StorageOptions Options, CMesh TMesh> | |
auto | ElasticHessianSparsity (TMesh const &mesh) -> math::linalg::SparsityPattern< typename TMesh::IndexType, Options > |
Construct the hessian matrix's sparsity pattern. | |
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg, class TDerivedHg> | |
auto | HyperElasticHessian (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedHg > const &Hg) -> Eigen::SparseMatrix< typename TDerivedHg::Scalar, Options, typename TDerivedE::Scalar > |
Construct the hessian matrix's sparse representation. | |
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedeg, class TDerivedHg> | |
auto | HyperElasticHessian (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedHg > const &Hg) -> Eigen::SparseMatrix< typename TDerivedHg::Scalar, Options, typename TMesh::IndexType > |
Construct the hessian matrix using mesh. | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedHg, Eigen::StorageOptions Options, class TDerivedH> | |
void | ToHyperElasticHessian (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedHg > const &Hg, math::linalg::SparsityPattern< typename TDerivedE::Scalar, Options > const &sparsity, Eigen::SparseCompressedBase< TDerivedH > &H) |
Construct the hessian matrix using mesh with sparsity pattern. | |
template<CMesh TMesh, class TDerivedHg, Eigen::StorageOptions Options, class TDerivedH> | |
void | ToHyperElasticHessian (TMesh const &mesh, Eigen::DenseBase< TDerivedHg > const &Hg, math::linalg::SparsityPattern< typename TMesh::IndexType, Options > const &sparsity, Eigen::SparseCompressedBase< typename TDerivedHg::Scalar > &H) |
Construct the hessian matrix using mesh with sparsity pattern. | |
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedHg> | |
auto | HyperElasticHessian (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedHg > const &Hg, math::linalg::SparsityPattern< typename TDerivedE::Scalar, Options > const &sparsity) -> Eigen::SparseMatrix< typename TDerivedHg::Scalar, Options, typename TDerivedE::Scalar > |
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedHg> | |
auto | HyperElasticHessian (TMesh const &mesh, Eigen::DenseBase< TDerivedHg > const &Hg, math::linalg::SparsityPattern< typename TMesh::IndexType, Options > const &sparsity) -> Eigen::SparseMatrix< typename TDerivedHg::Scalar, Options, typename TMesh::IndexType > |
Construct the hessian matrix using mesh with sparsity pattern. | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedGg, class TDerivedG> | |
void | ToHyperElasticGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGg > const &Gg, Eigen::PlainObjectBase< TDerivedG > &G) |
Compute the gradient vector into existing output. | |
template<CElement TElement, int Dims, physics::CHyperElasticEnergy THyperElasticEnergy, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx, class TDerivedGg, class TDerivedG> | |
void | ToHyperElasticGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x, Eigen::MatrixBase< TDerivedGg > const &Gg, Eigen::PlainObjectBase< TDerivedG > &G) |
Compute the gradient vector into existing output. | |
template<CMesh TMesh, class TDerivedeg, class TDerivedGg, class TDerivedOut> | |
void | ToHyperElasticGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGg > const &Gg, Eigen::PlainObjectBase< TDerivedOut > &G) |
Compute the gradient vector using mesh (updates existing output) | |
template<CMesh TMesh, physics::CHyperElasticEnergy THyperElasticEnergy, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx, class TDerivedGg, class TDerivedG> | |
void | ToHyperElasticGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x, Eigen::MatrixBase< TDerivedGg > const &Gg, Eigen::PlainObjectBase< TDerivedG > &G) |
Compute the gradient vector into existing output. | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedGg> | |
auto | HyperElasticGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGg > const &Gg) -> Eigen::Vector< typename TDerivedGg::Scalar, Eigen::Dynamic > |
Compute the gradient vector (allocates output) | |
template<CMesh TMesh, class TDerivedeg, class TDerivedGg> | |
auto | HyperElasticGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGg > const &Gg) |
Compute the gradient vector (allocates output) | |
template<CElement TElement, int Dims, physics::CHyperElasticEnergy THyperElasticEnergy, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx, class TDerivedGg> | |
auto | HyperElasticGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x, Eigen::MatrixBase< TDerivedGg > const &Gg) |
Compute the gradient vector (allocates output) | |
template<CMesh TMesh, physics::CHyperElasticEnergy THyperElasticEnergy, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx, class TDerivedGg> | |
auto | HyperElasticGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x, Eigen::MatrixBase< TDerivedGg > const &Gg) |
Compute the gradient vector (allocates output) | |
template<physics::CHyperElasticEnergy THyperElasticEnergy, CMesh TMesh, class TDerivedeg, class TDerivedGg> | |
auto | HyperElasticGradient (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGg > const &Gg) |
Compute the gradient vector using mesh (allocates output) | |
template<class TDerivedUg> | |
auto | HyperElasticPotential (Eigen::DenseBase< TDerivedUg > const &Ug) -> typename TDerivedUg::Scalar |
Compute the total elastic potential. | |
template<CElement TElement, int Dims, physics::CHyperElasticEnergy THyperElasticEnergy, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedmug, class TDerivedlambdag, class TDerivedx> | |
auto | HyperElasticPotential (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, Eigen::DenseBase< TDerivedmug > const &mug, Eigen::DenseBase< TDerivedlambdag > const &lambdag, Eigen::MatrixBase< TDerivedx > const &x) |
Compute the total elastic potential from element data and quadrature. | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedGg, class TDerivedOut> | |
void | ToHyperElasticGradient (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedGg > const &Gg, Eigen::PlainObjectBase< TDerivedOut > &G) |
template<CElement TElement, class TDerivedX, class TDerivedx, common::CFloatingPoint TScalar = typename TDerivedX::Scalar> | |
auto | Jacobian (Eigen::MatrixBase< TDerivedX > const &X, Eigen::MatrixBase< TDerivedx > const &x) -> Eigen::Matrix< TScalar, TDerivedx::RowsAtCompileTime, TElement::kDims > |
Given a map \( x(\xi) = \sum_i x_i N_i(\xi) \), where \( \xi \in \Omega^{\text{ref}}
\) is in reference space coordinates, computes \( \nabla_\xi x \). | |
template<class TDerived> | |
auto | DeterminantOfJacobian (Eigen::MatrixBase< TDerived > const &J) -> typename TDerived::Scalar |
Computes the determinant of a (potentially non-square) Jacobian matrix. | |
template<CElement TElement, int QuadratureOrder, class TDerivedE, class TDerivedX> | |
auto | DeterminantOfJacobian (Eigen::DenseBase< TDerivedE > const &E, Eigen::MatrixBase< TDerivedX > const &X) -> Eigen::Matrix< typename TDerivedX::Scalar, TElement::template QuadratureType< QuadratureOrder, typename TDerivedX::Scalar >::kPoints, Eigen::Dynamic > |
Computes the determinant of the Jacobian matrix at element quadrature points. | |
template<int QuadratureOrder, CMesh TMesh> | |
auto | DeterminantOfJacobian (TMesh const &mesh) -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::template QuadratureType< QuadratureOrder, typename TMesh::ScalarType >::kPoints, Eigen::Dynamic > |
template<CElement TElement, class TDerivedEg, class TDerivedX, class TDerivedXi> | |
auto | DeterminantOfJacobianAt (Eigen::DenseBase< TDerivedEg > const &Eg, Eigen::MatrixBase< TDerivedX > const &X, Eigen::MatrixBase< TDerivedXi > const &Xi) -> Eigen::Vector< typename TDerivedX::Scalar, Eigen::Dynamic > |
Computes the determinant of the Jacobian matrix at element quadrature points. | |
template<CElement TElement, class TDerivedE, class TDerivedX, class TDerivedeg, class TDerivedXi> | |
auto | DeterminantOfJacobianAt (Eigen::DenseBase< TDerivedE > const &E, Eigen::MatrixBase< TDerivedX > const &X, Eigen::MatrixBase< TDerivedeg > const &eg, Eigen::MatrixBase< TDerivedXi > const &Xi) -> Eigen::Vector< typename TDerivedX::Scalar, Eigen::Dynamic > |
Computes the determinant of the Jacobian matrix at element quadrature points. | |
template<CElement TElement, class TDerivedX, class TDerivedx, common::CFloatingPoint TScalar = typename TDerivedX::Scalar> | |
auto | ReferencePosition (Eigen::MatrixBase< TDerivedX > const &X, Eigen::MatrixBase< TDerivedx > const &x, int const maxIterations=5, TScalar const eps=1e-10) -> Eigen::Vector< TScalar, TElement::kDims > |
Computes the reference position \( \xi \) such that \( x(\xi) = x \). This inverse problem is solved using Gauss-Newton iterations. | |
template<CElement TElement, class TDerivedEg, class TDerivedX, class TDerivedXg, common::CFloatingPoint TScalar = typename TDerivedXg::Scalar> | |
auto | ReferencePositions (Eigen::DenseBase< TDerivedEg > const &Eg, Eigen::MatrixBase< TDerivedX > const &X, Eigen::MatrixBase< TDerivedXg > const &Xg, int maxIterations=5, TScalar eps=1e-10) -> Eigen::Matrix< TScalar, TElement::kDims, Eigen::Dynamic > |
Computes reference positions \( \xi \) such that \( X(\xi) = X_n \) for every point in \( \mathbf{X} \in \mathbb{R}^{d \times n} \). | |
template<CElement TElement, class TDerivedE, class TDerivedX, class TDerivedEg, class TDerivedXg> | |
auto | ReferencePositions (Eigen::DenseBase< TDerivedE > const &E, Eigen::MatrixBase< TDerivedX > const &X, Eigen::DenseBase< TDerivedEg > const &eg, Eigen::MatrixBase< TDerivedXg > const &Xg, int maxIterations=5, typename TDerivedXg::Scalar eps=1e-10) -> Eigen::Matrix< typename TDerivedXg::Scalar, TElement::kDims, Eigen::Dynamic > |
template<CMesh TMesh, class TDerivedEg, class TDerivedXg> | |
auto | ReferencePositions (TMesh const &mesh, Eigen::DenseBase< TDerivedEg > const &eg, Eigen::MatrixBase< TDerivedXg > const &Xg, int maxIterations=5, typename TMesh::ScalarType eps=1e-10) -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::kDims, Eigen::Dynamic > |
Computes reference positions \( \xi \) such that \( X(\xi) = X_n \) for every point in \( \mathbf{X} \in \mathbb{R}^{d \times n} \). | |
template<CElement TElement, int Dims, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedIn, class TDerivedOut> | |
void | GemmLaplacian (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, int dims, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Compute Laplacian matrix-vector multiply \( Y += \mathbf{L} X \). | |
template<CMesh TMesh, class TDerivedeg, class TDerivedwg, class TDerivedGNeg, class TDerivedIn, class TDerivedOut> | |
void | GemmLaplacian (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, int dims, Eigen::MatrixBase< TDerivedIn > const &X, Eigen::DenseBase< TDerivedOut > &Y) |
Compute Laplacian matrix-vector multiply \( Y += \mathbf{L} X \) using mesh. | |
template<CElement TElement, int Dims, Eigen::StorageOptions Options, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedGNeg> | |
auto | LaplacianMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, int dims=1) -> Eigen::SparseMatrix< typename TDerivedGNeg::Scalar, Options, typename TDerivedE::Scalar > |
Construct the Laplacian operator's sparse matrix representation. | |
template<Eigen::StorageOptions Options, CMesh TMesh, class TDerivedeg, class TDerivedwg, class TDerivedGNeg> | |
auto | LaplacianMatrix (TMesh const &mesh, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedGNeg > const &GNeg, int dims=1) -> Eigen::SparseMatrix< typename TDerivedGNeg::Scalar, Options, typename TMesh::IndexType > |
Construct the Laplacian operator's sparse matrix representation using mesh. | |
template<CElement TElement, class TDerivedE, class TDerivedeg, class TDerivedwg, class TDerivedNeg, class TDerivedFeg> | |
auto | LoadVectors (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes, Eigen::DenseBase< TDerivedeg > const &eg, Eigen::DenseBase< TDerivedwg > const &wg, Eigen::MatrixBase< TDerivedNeg > const &Neg, Eigen::MatrixBase< TDerivedFeg > const &Feg) -> Eigen::Matrix< typename TDerivedFeg::Scalar, TDerivedFeg::RowsAtCompileTime, Eigen::Dynamic > |
Computes the load vector for a given FEM mesh. | |
template<typename TElement, typename TN> | |
auto | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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 | 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) | |
template<CElement TElement, int QuadratureOrder, class TDerivedDetJe> | |
void | ToMeshQuadratureWeights (Eigen::MatrixBase< TDerivedDetJe > &detJeThenWg) |
Computes the inner product weights \( \mathbf{w}_{ge} \in \mathbb{R}^{|G^e| \times |E|}
\) such that \( \int_\Omega \cdot d\Omega = \sum_e \sum_g w_{ge} \cdot \). | |
template<CElement TElement, int QuadratureOrder, class TDerivedE, class TDerivedX> | |
auto | MeshQuadratureWeights (Eigen::MatrixBase< TDerivedE > const &E, Eigen::MatrixBase< TDerivedX > const &X) -> Eigen::Matrix< typename TDerivedX::Scalar, TElement::template QuadratureType< QuadratureOrder, typename TDerivedX::Scalar >::kPoints, Eigen::Dynamic > |
Computes the inner product weights \( \mathbf{w}_{ge} \in \mathbb{R}^{|G^e| \times |E|}
\) such that \( \int_\Omega \cdot d\Omega = \sum_e \sum_g w_{ge} \cdot \). | |
template<int QuadratureOrder, CMesh TMesh> | |
auto | MeshQuadratureWeights (TMesh const &mesh) -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::template QuadratureType< QuadratureOrder, typename TMesh::ScalarType >::kPoints, Eigen::Dynamic > |
Computes the inner product weights \( \mathbf{w}_{ge} \in \mathbb{R}^{|G^e| \times |E|}
\) such that \( \int_\Omega \cdot d\Omega = \sum_e \sum_g w_{ge} \cdot \). | |
template<common::CIndex TIndex> | |
auto | MeshQuadratureElements (TIndex nElements, TIndex nQuadPtsPerElement) |
Computes the element quadrature points indices for each quadrature point, given the number of elements and the quadrature weights matrix. | |
template<common::CIndex TIndex, class TDerivedwg> | |
auto | MeshQuadratureElements (TIndex nElements, Eigen::DenseBase< TDerivedwg > const &wg) |
Computes the element quadrature points indices for each quadrature point, given the number of elements and the quadrature weights matrix. | |
template<class TDerivedE, class TDerivedwg> | |
auto | MeshQuadratureElements (Eigen::DenseBase< TDerivedE > const &E, Eigen::DenseBase< TDerivedwg > const &wg) |
Computes the element quadrature points indices for each quadrature point. | |
template<CElement TElement, int QuadratureOrder, common::CArithmetic TScalar = Scalar, common::CIndex TIndex = Index> | |
auto | MeshReferenceQuadraturePoints (TIndex nElements) |
Computes the element quadrature points in reference element space. | |
template<int QuadratureOrder, CMesh TMesh> | |
auto | MeshReferenceQuadraturePoints (TMesh const &mesh) |
Computes the element quadrature points in reference element space. | |
template<CElement TElement, int QuadratureOrder, common::CFloatingPoint TScalar = Scalar> | |
auto | ElementShapeFunctions () -> Eigen::Matrix< TScalar, TElement::kNodes, TElement::template QuadratureType< QuadratureOrder, TScalar >::kPoints > |
Computes shape functions at element quadrature points for a polynomial quadrature rule of order QuadratureOrder. | |
template<CElement TElement, int QuadratureOrder, common::CFloatingPoint TScalar, class TDerivedE> | |
auto | ShapeFunctionMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::Index nNodes) -> Eigen::SparseMatrix< TScalar, Eigen::RowMajor, typename TDerivedE::Scalar > |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the element quadrature points. | |
template<int QuadratureOrder, CMesh TMesh> | |
auto | ShapeFunctionMatrix (TMesh const &mesh) -> Eigen::SparseMatrix< typename TMesh::ScalarType, Eigen::RowMajor, typename TMesh::IndexType > |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the element quadrature points. | |
template<CElement TElement, class TDerivedE, class TDerivedXi, common::CIndex TIndex = typename TDerivedE::Scalar> | |
auto | ShapeFunctionMatrixAt (Eigen::DenseBase< TDerivedE > const &E, TIndex nNodes, Eigen::MatrixBase< TDerivedXi > const &Xi) -> Eigen::SparseMatrix< typename TDerivedXi::Scalar, Eigen::RowMajor, TIndex > |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the element quadrature points. | |
template<CMesh TMesh, class TDerivedEg, class TDerivedXi> | |
auto | ShapeFunctionMatrixAt (TMesh const &mesh, Eigen::DenseBase< TDerivedEg > const &eg, Eigen::MatrixBase< TDerivedXi > const &Xi) -> Eigen::SparseMatrix< typename TMesh::ScalarType, Eigen::RowMajor, typename TMesh::IndexType > |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the given evaluation points. | |
template<CElement TElement, class TDerivedXi> | |
auto | ShapeFunctionsAt (Eigen::DenseBase< TDerivedXi > const &Xi) -> Eigen::Matrix< typename TDerivedXi::Scalar, TElement::kNodes, Eigen::Dynamic > |
Compute shape functions at the given reference positions. | |
template<CElement TElement, class TDerivedXi, class TDerivedX> | |
auto | ElementShapeFunctionGradients (Eigen::MatrixBase< TDerivedXi > const &Xi, Eigen::MatrixBase< TDerivedX > const &X) -> Eigen::Matrix< typename TDerivedX::Scalar, TElement::kNodes, TDerivedX::RowsAtCompileTime > |
Computes gradients \( \mathbf{G}(X) = \nabla_X \mathbf{N}_e(X) \) of FEM basis functions. Given some FEM function. | |
template<CElement TElement, int Dims, int QuadratureOrder, class TDerivedE, class TDerivedX, common::CFloatingPoint TScalar = typename TDerivedX::Scalar, common::CIndex TIndex = typename TDerivedE::Scalar> | |
auto | ShapeFunctionGradients (Eigen::MatrixBase< TDerivedE > const &E, Eigen::MatrixBase< TDerivedX > const &X) -> Eigen::Matrix< TScalar, TElement::kNodes, Eigen::Dynamic > |
Computes nodal shape function gradients at each element quadrature points. | |
template<int QuadratureOrder, CMesh TMesh> | |
auto | ShapeFunctionGradients (TMesh const &mesh) -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::kNodes, Eigen::Dynamic > |
Computes nodal shape function gradients at each element quadrature points. | |
template<CElement TElement, int Dims, class TDerivedEg, class TDerivedX, class TDerivedXi> | |
auto | ShapeFunctionGradientsAt (Eigen::DenseBase< TDerivedEg > const &Eg, Eigen::MatrixBase< TDerivedX > const &X, Eigen::MatrixBase< TDerivedXi > const &Xi) -> Eigen::Matrix< typename TDerivedXi::Scalar, TElement::kNodes, Eigen::Dynamic > |
Computes nodal shape function gradients at evaluation points Xi. | |
template<CMesh TMesh, class TDerivedEg, class TDerivedXi> | |
auto | ShapeFunctionGradientsAt (TMesh const &mesh, Eigen::DenseBase< TDerivedEg > const &eg, Eigen::MatrixBase< TDerivedXi > const &Xi) -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::kNodes, Eigen::Dynamic > |
Computes nodal shape function gradients at evaluation points Xg. | |
Finite Element Method (FEM)
using pbat::fem::Hexahedron = typename detail::Hexahedron<Order> |
using pbat::fem::Line = typename detail::Line<Order> |
using pbat::fem::Quadrilateral = typename detail::Quadrilateral<Order> |
using pbat::fem::Tetrahedron = typename detail::Tetrahedron<Order> |
using pbat::fem::Triangle = typename detail::Triangle<Order> |
auto pbat::fem::DeformationGradient | ( | Eigen::MatrixBase< TDerivedx > const & | x, |
Eigen::MatrixBase< TDerivedX > const & | GP ) -> Matrix<TDerivedx::RowsAtCompileTime, TElement::kDims> |
Computes the deformation gradient \( \frac{\partial \mathbf{x}(X)}{\partial X} \) of the deformation map \( \mathbf{x}(X) \).
If the problem is discretized with displacement coefficients \( \mathbf{u} = \mathbf{x}(X) - X \), then simply feed this function with argument \( \mathbf{x} = X + \mathbf{u} \).
TDerivedU | Eigen matrix expression |
TDerivedX | Eigen matrix expression |
TElement | FEM element type |
x | Matrix of column-wise position nodal coefficients |
GP | Basis function gradients |
auto pbat::fem::DeterminantOfJacobian | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::MatrixBase< TDerivedX > const & | X ) -> Eigen::Matrix< typename TDerivedX::Scalar, TElement::template QuadratureType<QuadratureOrder, typename TDerivedX::Scalar>::kPoints, Eigen::Dynamic> |
Computes the determinant of the Jacobian matrix at element quadrature points.
TElement | FEM element type |
QuadratureOrder | Quadrature order |
TDerivedE | Eigen matrix expression for element matrix |
TDerivedX | Eigen matrix expression for nodal positions |
E | |# elem nodes| x |# elems| element matrix |
X | |# dims| x |# nodes| mesh nodal position matrix |
|# elem quad.pts.| x |# elems|
matrix of jacobian determinants at element quadrature points auto pbat::fem::DeterminantOfJacobian | ( | Eigen::MatrixBase< TDerived > const & | J | ) | -> typename TDerived::Scalar |
Computes the determinant of a (potentially non-square) Jacobian matrix.
If the Jacobian matrix is square, the determinant is computed directly, otherwise the singular values \( \sigma_i \) are computed and their product \( \Pi_i \sigma_i \) is returned.
TDerived | Eigen matrix expression |
J | Jacobian matrix |
auto pbat::fem::DeterminantOfJacobian | ( | TMesh const & | mesh | ) | -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::template QuadratureType<QuadratureOrder, typename TMesh::ScalarType>:: kPoints, Eigen::Dynamic> |
QuadratureOrder | |
TMesh |
mesh |
auto pbat::fem::DeterminantOfJacobianAt | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::MatrixBase< TDerivedX > const & | X, | ||
Eigen::MatrixBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedXi > const & | Xi ) -> Eigen::Vector<typename TDerivedX::Scalar, Eigen::Dynamic> |
Computes the determinant of the Jacobian matrix at element quadrature points.
TElement | FEM element type |
TDerivedE | Eigen matrix expression for element matrix |
TDerivedX | Eigen matrix expression for mesh nodal positions |
TDerivedeg | Eigen matrix expression for element indices at evaluation points |
TDerivedXi | Eigen matrix expression for evaluation points in reference space |
E | |# elem nodes| x |# elems| mesh element matrix |
X | |# dims| x |# nodes| mesh nodal position matrix |
eg | |# eval.pts.| element indices at evaluation points |
Xi | |# dims| x |# eval.pts.| evaluation points in reference space |
|# eval.pts.| x 1
vector of jacobian determinants at evaluation points auto pbat::fem::DeterminantOfJacobianAt | ( | Eigen::DenseBase< TDerivedEg > const & | Eg, |
Eigen::MatrixBase< TDerivedX > const & | X, | ||
Eigen::MatrixBase< TDerivedXi > const & | Xi ) -> Eigen::Vector<typename TDerivedX::Scalar, Eigen::Dynamic> |
Computes the determinant of the Jacobian matrix at element quadrature points.
TElement | FEM element type |
TDerivedEg | Eigen matrix expression for element indices at evaluation points |
TDerivedX | Eigen matrix expression for mesh nodal positions |
TDerivedXi | Eigen matrix expression for evaluation points in reference space |
Eg | |# elem nodes| x |# eval.pts.| element indices at evaluation points |
X | |# dims| x |# nodes| mesh nodal position matrix |
Xi | |# dims| x |# eval.pts.| evaluation points in reference space |
|# eval.pts.| x 1
vector of jacobian determinants at evaluation points auto pbat::fem::ElasticHessianSparsity | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg ) -> math::linalg::SparsityPattern<typename TDerivedE::Scalar, Options> |
Construct the hessian matrix's sparsity pattern.
TElement | Element type |
Dims | Number of spatial dimensions |
Options | Storage options for the matrix |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
auto pbat::fem::ElasticHessianSparsity | ( | TMesh const & | mesh | ) | -> math::linalg::SparsityPattern<typename TMesh::IndexType, Options> |
Construct the hessian matrix's sparsity pattern.
mesh | Mesh containing element connectivity and node positions |
|
inline |
Compute a matrix of horizontally stacked element mass matrices for all quadrature points.
TElement | Element type |
TN | Type of the shape function vector at the quadrature point |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
Neg | |# nodes per element| x |# quad.pts.| shape function matrix at all quadrature points |
wg | |# quad.pts.| x 1 quadrature weights (including Jacobian determinant) |
rhog | |# quad.pts.| x 1 mass density at quadrature points |
|# elem nodes| x |# elem nodes * # quad.pts.|
matrix of stacked element mass matrices
|
inline |
Compute the element mass matrix at a single quadrature point.
TElement | Element type |
TN | Type of the shape function vector at the quadrature point |
TScalar | Scalar type (e.g., double) |
N | Shape function vector at the quadrature point (size: |# nodes per element|) |
w | Quadrature weight (including Jacobian determinant) |
rho | Density at the quadrature point |
auto pbat::fem::ElementShapeFunctionGradients | ( | Eigen::MatrixBase< TDerivedXi > const & | Xi, |
Eigen::MatrixBase< TDerivedX > const & | X ) -> Eigen::Matrix<typename TDerivedX::Scalar, TElement::kNodes, TDerivedX::RowsAtCompileTime> |
Computes gradients \( \mathbf{G}(X) = \nabla_X \mathbf{N}_e(X) \) of FEM basis functions. Given some FEM function.
\[\begin{align*} f(X) &= \begin{bmatrix} \dots & \mathbf{f}_i & \dots \end{bmatrix} \begin{bmatrix} \vdots \\ \phi_i(X) \\ \vdots \end{bmatrix} \\ &= \mathbf{F}_e \mathbf{N}_e(X) \end{align*} \]
where \( \phi_i(X) = N_i(\xi(X)) \) are the nodal basis functions, we have that
\[\begin{align*} \nabla_X f(X) &= \mathbf{F}_e \nabla_X \mathbf{N}_e(X) \\ &= \mathbf{F}_e \nabla_\xi \mathbf{N}_e(\xi(X)) \frac{\partial \xi}{\partial X} \end{align*} \]
Recall that
\[\begin{align*} \frac{\partial \xi}{\partial X} &= (\frac{\partial X}{\partial \xi})^{-1} \\ &= \left( \mathbf{X}_e \nabla_\xi \mathbf{N}_e \right)^{-1} \end{align*} \]
where \( \mathbf{X}_e = \begin{bmatrix} \dots & \mathbf{X}_i & \dots \end{bmatrix} \) are the element nodal positions. We thus have that
\[\mathbf{G}(X)^T = \nabla_X \mathbf{N}_e(X)^T = ( \mathbf{X}_e \nabla_\xi \mathbf{N}_e )^{-T} \nabla_\xi \mathbf{N}_e^T \]
where the left-multiplication by \( ( \mathbf{X}_e \nabla_\xi \mathbf{N}_e )^{-T} \) is carried out via LU/SVD decomposition for a square/rectangular Jacobian.
TDerivedXi | Eigen dense expression type for reference position |
TDerivedX | Eigen dense expression type for element nodal positions |
TElement | Element type |
Xi | |# elem dims| x 1 point in reference element at which to evaluate the gradients |
X | |# dims| x |# nodes| element vertices, i.e. nodes of affine element |
|# nodes| x |# dims|
matrix of basis function gradients in rows auto pbat::fem::ElementShapeFunctions | ( | ) | -> Eigen::Matrix< TScalar, TElement::kNodes, TElement::template QuadratureType<QuadratureOrder, TScalar>::kPoints> |
Computes shape functions at element quadrature points for a polynomial quadrature rule of order QuadratureOrder.
TElement | Element type |
QuadratureOrder | Quadrature order |
TScalar | Floating point type, defaults to Scalar |
|# element nodes| x |# elem. quad.pts.|
matrix of nodal shape function values at quadrature points
|
inline |
Compute gradient matrix-vector multiply \( Y += \mathbf{G} X \).
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
X | |# nodes| x |# cols| input matrix |
Y | |# dims * # quad.pts.| x |# cols| output matrix |
|
inline |
Compute gradient matrix-vector multiply \( Y += \mathbf{G} X \) using mesh.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
X | |# nodes| x |# cols| input matrix |
Y | |# dims * # quad.pts.| x |# cols| output matrix |
|
inline |
Apply the hessian matrix as a linear operator \( Y += \mathbf{H} X \).
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedHg | Type of the element hessian matrices |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points |
X | |# dims * # nodes| x |# cols| input matrix |
Y | |# dims * # nodes| x |# cols| output matrix |
|
inline |
Apply the hessian matrix as a linear operator \( Y += \mathbf{H} X \) using mesh.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedHg | Type of the element hessian matrices |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points |
X | |# dims * # nodes| x |# cols| input matrix |
Y | |# dims * # nodes| x |# cols| output matrix |
|
inline |
Compute Laplacian matrix-vector multiply \( Y += \mathbf{L} X \).
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
dims | Dimensionality of the image of the FEM function space |
X | |# nodes * dims| x |# cols| input matrix |
Y | |# nodes * dims| x |# cols| output matrix |
|
inline |
Compute Laplacian matrix-vector multiply \( Y += \mathbf{L} X \) using mesh.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
dims | Dimensionality of the image of the FEM function space |
X | |# nodes * dims| x |# cols| input matrix |
Y | |# nodes * dims| x |# cols| output matrix |
|
inline |
Compute mass matrix-matrix multiply \( Y += \mathbf{M} X \) using precomputed element mass matrices.
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Me | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
X | |# nodes * dims| x |# cols| input matrix |
Y | |# nodes * dims| x |# cols| output matrix |
|
inline |
Compute mass matrix-matrix multiply \( Y += \mathbf{M} X \) using mesh and precomputed element mass matrices.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
TDerivedIn | Type of the input matrix |
TDerivedOut | Type of the output matrix |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Me | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
X | |# nodes * dims| x |# cols| input matrix |
Y | |# nodes * dims| x |# cols| output matrix |
auto pbat::fem::GradientMatrix | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg ) -> Eigen::SparseMatrix<typename TDerivedGNeg::Scalar, Options, typename TDerivedE::Scalar> |
Construct the gradient operator's sparse matrix representation.
TElement | Element type |
Dims | Number of spatial dimensions |
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
auto pbat::fem::GradientMatrix | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg ) -> Eigen::SparseMatrix<typename TDerivedGNeg::Scalar, Options, typename TMesh::IndexType> |
Construct the gradient operator's sparse matrix representation using mesh.
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
auto pbat::fem::GradientSegmentWrtDofs | ( | TMatrixGF const & | GF, |
TMatrixGP const & | GP, | ||
auto | i ) -> math::linalg::mini::SVector<ScalarType, Dims> |
Computes \( \frac{\partial \Psi}{\partial \mathbf{x}_i} \in \mathbb{R}^d \), i.e. the gradient of a scalar function \( \Psi \) w.r.t. the \( i^{\text{th}} \) node's degrees of freedom.
TElement | FEM element type |
Dims | Problem dimensionality |
TMatrixGF | Mini matrix type |
TMatrixGP | Mini matrix type |
Scalar | Coefficient type of returned gradient |
GF | \( d^2 \times 1 \) gradient of \( \Psi \) w.r.t. vectorized jacobian \(\text{vec}(\mathbf{F}) \) |
GP | \( * \times d \) basis function gradients \( \nabla \mathbf{N}_i \) |
i | Basis function index |
auto pbat::fem::GradientWrtDofs | ( | TMatrixGF const & | GF, |
TMatrixGP const & | GP ) -> math::linalg::mini::SVector<ScalarType, TElement::kNodes * Dims> |
Computes gradient w.r.t. FEM degrees of freedom \( x \) of scalar function \(\Psi(\mathbf{F}) \), where \( F \) is the jacobian of \( x \), via chain rule. This is effectively a rank-3 to rank-1 tensor contraction.
Let \( \mathbf{g}_i \) be the \( i^{\text{th}} \) basis function gradient, then
\[\frac{\partial \mathbf{F}}{\partial \mathbf{x}_i} = \frac{\partial}{\partial \mathbf{x}_i} \mathbf{x}_i \mathbf{g}_i^T, \mathbf{x}_i \in R^d, \mathbf{g}_i \in R^d \]
Thus,
\[\frac{\partial \text{vec}(\mathbf{F})}{\partial \mathbf{x}_i} = \mathbf{g}_i \otimes \mathbf{I}_{d x d} \in R^{d^2 \times d} \]
With \( \mathbf{G}_F = \nabla_{\text{vec}(\mathbf{F})} \Psi \in \mathbb{R}^{d^2} \), and \( \mathbf{G}_F^k \) the \( k^\text{th} \; d \times 1 \) block of \( \mathbf{G}_F \), then
\[\frac{\partial \Psi}{\partial \mathbf{x}_i} = \frac{\partial \Psi}{\partial \text{vec}(\mathbf{F})} \frac{\partial \text{vec}(\mathbf{F})}{\partial \mathbf{x}_i} = \sum_{k=1}^{d} \mathbf{G}_F^k \mathbf{g}_{ik} \]
TElement | FEM element type |
Dims | Problem dimensionality |
TMatrixGF | Mini matrix type |
TMatrixGP | Mini matrix type |
Scalar | Coefficient type of returned gradient |
GF | \( d^2 \times 1 \) gradient of \( \Psi \) w.r.t. vectorized jacobian \(\text{vec}(\mathbf{F}) \) |
GP | \( * \times d \) basis function gradients \( \nabla \mathbf{N}_i \) |
auto pbat::fem::HessianBlockWrtDofs | ( | TMatrixHF const & | HF, |
TMatrixGP const & | GP, | ||
auto | i, | ||
auto | j ) -> math::linalg::mini::SMatrix<ScalarType, Dims, Dims> |
Computes \( \frac{\partial^2 \Psi}{\partial \mathbf{x}_i \partial \mathbf{x}_j} \), i.e. the hessian of a scalar function \( \Psi \) w.r.t. the \( i^{\text{th}} \) and \(j^{\text{th}} \) node's degrees of freedom.
TElement | FEM element type |
Dims | Problem dimensionality |
TMatrixHF | Mini matrix type |
TMatrixGP | Mini matrix type |
Scalar | Coefficient type of returned hessian |
HF | \( d^2 \times d^2 \) hessian of \( \Psi \) w.r.t. vectorized jacobian \(\text{vec}(\mathbf{F}) \), i.e. \( \frac{\partial^2 \Psi}{\partial \text{vec}(\mathbf{F})^2} \) |
GP | \( * \times d \) basis function gradients \( \nabla \mathbf{N}_i \) |
i | Basis function index |
j | Basis function index |
auto pbat::fem::HessianWrtDofs | ( | TMatrixHF const & | HF, |
TMatrixGP const & | GP ) -> math::linalg::mini::SMatrix<ScalarType, TElement::kNodes * Dims, TElement::kNodes * Dims> |
Computes hessian w.r.t. FEM degrees of freedom \( x \) of scalar function \(\Psi(\mathbf{F}) \), where \( F \) is the jacobian of \( x \), via chain rule. This is effectively a rank-4 to rank-2 tensor contraction.
Let \( \mathbf{g}_i \) be the \( i^{\text{th}} \) basis function gradient, then
\[\frac{\partial \mathbf{F}}{\partial \mathbf{x}_i} = \frac{\partial}{\partial \mathbf{x}_i} \mathbf{x}_i \mathbf{g}_i^T, \mathbf{x}_i \in R^d, \mathbf{g}_i \in R^d \]
Thus,
\[\frac{\partial \text{vec}(\mathbf{F})}{\partial \mathbf{x}_i} = \mathbf{g}_i \otimes \mathbf{I}_{d x d} \in R^{d^2 \times d} \]
With \( \mathbf{H}_F = \nabla^2_{\text{vec}(\mathbf{F})} \Psi \in \mathbb{R}^{d^2 \times d^2} \), and \( \mathbf{H}_F^{uv} \) the \( (u,v)^\text{th} \; d \times d \) block of \(\mathbf{H}_F \), then
\[\frac{\partial^2 \Psi}{\partial \mathbf{x}_i \partial \mathbf{x}_j} = \frac{\partial^2 \Psi}{\partial \text{vec}(\mathbf{F})^2} \frac{\partial^2 \text{vec}(\mathbf{F})}{\partial \mathbf{x}_i \partial \mathbf{x}_j} = \sum_{u=1}^{d} \sum_{v=1}^{d} \mathbf{H}_F^{uv} \mathbf{g}_{iu} \mathbf{g}_{jv} \]
TElement | FEM element type |
Dims | Problem dimensionality |
TMatrixHF | Mini matrix type |
TMatrixGP | Mini matrix type |
Scalar | Coefficient type of returned hessian |
HF | \( d^2 \times d^2 \) hessian of \( \Psi \) w.r.t. vectorized jacobian \(\text{vec}(\mathbf{F}) \), i.e. \( \frac{\partial^2 \Psi}{\partial \text{vec}(\mathbf{F})^2} \) |
GP | \( * \times d \) basis function gradients \( \nabla \mathbf{N}_i \) |
auto pbat::fem::HyperElasticGradient | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg, | ||
Eigen::DenseBase< TDerivedmug > const & | mug, | ||
Eigen::DenseBase< TDerivedlambdag > const & | lambdag, | ||
Eigen::MatrixBase< TDerivedx > const & | x, | ||
Eigen::MatrixBase< TDerivedGg > const & | Gg ) |
Compute the gradient vector (allocates output)
TElement | Element type |
Dims | Number of spatial dimensions |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedmug | Type of the first Lame coefficients at quadrature points |
TDerivedlambdag | Type of the second Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
TDerivedGg | Type of the element gradient vectors at quadrature points |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element gradient vectors at quadrature points |
|# dims * # nodes| x 1
gradient vector auto pbat::fem::HyperElasticGradient | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedGg > const & | Gg ) -> Eigen::Vector<typename TDerivedGg::Scalar, Eigen::Dynamic> |
Compute the gradient vector (allocates output)
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGg | Type of the element gradient vectors |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
|# dims * # nodes| x 1
gradient vector auto pbat::fem::HyperElasticGradient | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg, | ||
Eigen::DenseBase< TDerivedmug > const & | mug, | ||
Eigen::DenseBase< TDerivedlambdag > const & | lambdag, | ||
Eigen::MatrixBase< TDerivedx > const & | x, | ||
Eigen::MatrixBase< TDerivedGg > const & | Gg ) |
Compute the gradient vector (allocates output)
TMesh | Type of the mesh |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedmug | Type of the first Lame coefficients at quadrature points |
TDerivedlambdag | Type of the second Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
TDerivedGg | Type of the element gradient vectors at quadrature points |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element gradient vectors at quadrature points |
|# dims| x |# nodes|
gradient vector auto pbat::fem::HyperElasticGradient | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedGg > const & | Gg ) |
Compute the gradient vector (allocates output)
TMesh | Type of the mesh |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGg | Type of the element gradient vectors |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
|# dims * # nodes| x 1
gradient vector auto pbat::fem::HyperElasticGradient | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedGg > const & | Gg ) |
Compute the gradient vector using mesh (allocates output)
THyperElasticEnergy | Hyper elastic energy type |
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGg | Type of the element gradient vectors |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
|# dims * # nodes| x 1
gradient vector auto pbat::fem::HyperElasticHessian | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedHg > const & | Hg ) -> Eigen::SparseMatrix<typename TDerivedHg::Scalar, Options, typename TDerivedE::Scalar> |
Construct the hessian matrix's sparse representation.
TElement | Element type |
Dims | Number of spatial dimensions |
Options | Storage options for the matrix |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedHg | Type of the element hessian matrices |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| array of element elastic hessian matrices at quadrature points |
auto pbat::fem::HyperElasticHessian | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedHg > const & | Hg ) -> Eigen::SparseMatrix<typename TDerivedHg::Scalar, Options, typename TMesh::IndexType> |
Construct the hessian matrix using mesh.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedHg | Type of the element hessian matrices |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points |
auto pbat::fem::HyperElasticHessian | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedHg > const & | Hg, | ||
math::linalg::SparsityPattern< typename TMesh::IndexType, Options > const & | sparsity ) -> Eigen::SparseMatrix<typename TDerivedHg::Scalar, Options, typename TMesh::IndexType> |
Construct the hessian matrix using mesh with sparsity pattern.
TMesh | Type of the mesh |
TDerivedHg | Type of the element hessian matrices |
mesh | Mesh containing element connectivity and node positions |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points |
sparsity | Sparsity pattern for the hessian matrix |
auto pbat::fem::HyperElasticPotential | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg, | ||
Eigen::DenseBase< TDerivedmug > const & | mug, | ||
Eigen::DenseBase< TDerivedlambdag > const & | lambdag, | ||
Eigen::MatrixBase< TDerivedx > const & | x ) |
Compute the total elastic potential from element data and quadrature.
TElement | Element type |
Dims | Number of spatial dimensions |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedmug | Type of the first Lame coefficients at quadrature points |
TDerivedlambdag | Type of the second Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
auto pbat::fem::HyperElasticPotential | ( | Eigen::DenseBase< TDerivedUg > const & | Ug | ) | -> typename TDerivedUg::Scalar |
Compute the total elastic potential.
TDerivedUg | Type of the elastic potentials at quadrature points |
Ug | |# quad.pts.| x 1 array of elastic potentials at quadrature points |
auto pbat::fem::Jacobian | ( | Eigen::MatrixBase< TDerivedX > const & | X, |
Eigen::MatrixBase< TDerivedx > const & | x ) -> Eigen::Matrix<TScalar, TDerivedx::RowsAtCompileTime, TElement::kDims> |
Given a map \( x(\xi) = \sum_i x_i N_i(\xi) \), where \( \xi \in \Omega^{\text{ref}} \) is in reference space coordinates, computes \( \nabla_\xi x \).
TElement | FEM element type |
TDerived | Eigen matrix expression for map coefficients |
X | Reference space coordinates \( \xi \) |
x | Map coefficients \( \mathbf{x} \) |
auto pbat::fem::LaplacianMatrix | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg, | ||
int | dims = 1 ) -> Eigen::SparseMatrix<typename TDerivedGNeg::Scalar, Options, typename TDerivedE::Scalar> |
Construct the Laplacian operator's sparse matrix representation.
TElement | Element type |
Dims | Number of spatial dimensions |
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
dims | Dimensionality of the image of the FEM function space |
auto pbat::fem::LaplacianMatrix | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg, | ||
int | dims = 1 ) -> Eigen::SparseMatrix<typename TDerivedGNeg::Scalar, Options, typename TMesh::IndexType> |
Construct the Laplacian operator's sparse matrix representation using mesh.
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients at quadrature points |
dims | Dimensionality of the image of the FEM function space |
auto pbat::fem::LoadVectors | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedNeg > const & | Neg, | ||
Eigen::MatrixBase< TDerivedFeg > const & | Feg ) -> Eigen::Matrix<typename TDerivedFeg::Scalar, TDerivedFeg::RowsAtCompileTime, Eigen::Dynamic> |
Computes the load vector for a given FEM mesh.
TElement | Element type |
TDerivedE | Eigen matrix expression for element matrix |
TDerivedeg | Eigen vector expression for elements at quadrature points |
TDerivedwg | Eigen vector expression for quadrature weights |
TDerivedNeg | Eigen matrix expression for nodal shape function values at quadrature points |
TDerivedFeg | Eigen matrix expression for external load values at quadrature points |
E | |# elem. nodes| x |# elements| element matrix |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 array of elements associated with quadrature points |
wg | |# quad.pts.| x 1 quadrature weights |
Neg | |# elem. nodes| x |# quad.pts.| nodal shape function values at quadrature points |
Feg | |# dims| x |# quad.pts.| load vector values at quadrature points |
|# dims| x |# nodes|
matrix s.t. each output dimension's load vectors are stored in rows 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)
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
TDerivedNeg | Type of the shape functions at quadrature points |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights (including Jacobian determinants) |
rhog | |# quad.pts.| x 1 vector of density at quadrature points |
Neg | |# nodes per element| x |# quad.pts.| shape functions at quadrature points |
dims | Dimensionality of the image of the FEM function space |
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)
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Meg | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
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.
TElement | Element type |
Dims | Number of spatial dimensions |
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
TDerivedNeg | Type of the shape functions at quadrature points |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights (including Jacobian determinants) |
rhog | |# quad.pts.| x 1 vector of density at quadrature points |
Neg | |# nodes per element| x |# quad.pts.| shape functions at quadrature points |
dims | Dimensionality of the image of the FEM function space |
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.
TElement | Element type |
Dims | Number of spatial dimensions |
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Meg | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
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.
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
TDerivedNeg | Type of the shape functions at quadrature points |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights (including Jacobian determinants) |
rhog | |# quad.pts.| x 1 vector of density at quadrature points |
Neg | |# nodes per element| x |# quad.pts.| shape functions at quadrature points |
dims | Dimensionality of the image of the FEM function space |
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.
Options | Storage options for the matrix (default: Eigen::ColMajor) |
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Meg | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
auto pbat::fem::MeshQuadratureElements | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::DenseBase< TDerivedwg > const & | wg ) |
Computes the element quadrature points indices for each quadrature point.
TDerivedE | Eigen matrix expression for element matrix |
TDerivedwg | Eigen matrix expression for quadrature weights |
E | |# elem nodes| x |# elems| mesh element matrix |
wg | |# quad.pts.| x |# elems| mesh quadrature weights matrix |
|# quad.pts.| x |# elems|
matrix of element indices at quadrature points auto pbat::fem::MeshQuadratureElements | ( | TIndex | nElements, |
Eigen::DenseBase< TDerivedwg > const & | wg ) |
Computes the element quadrature points indices for each quadrature point, given the number of elements and the quadrature weights matrix.
TIndex | Index type |
TDerivedwg | Eigen matrix expression for quadrature weights |
nElements | Number of elements |
wg | |# quad.pts.| x |# elems| mesh quadrature weights matrix |
|# quad.pts.| x |# elems|
matrix of element indices at quadrature points auto pbat::fem::MeshQuadratureElements | ( | TIndex | nElements, |
TIndex | nQuadPtsPerElement ) |
Computes the element quadrature points indices for each quadrature point, given the number of elements and the quadrature weights matrix.
TIndex | Index type |
TDerivedwg | Eigen matrix expression for quadrature weights |
nElements | Number of elements |
wg | |# quad.pts.| x |# elems| mesh quadrature weights matrix |
|# quad.pts.| x |# elems|
matrix of element indices at quadrature points auto pbat::fem::MeshQuadratureWeights | ( | Eigen::MatrixBase< TDerivedE > const & | E, |
Eigen::MatrixBase< TDerivedX > const & | X ) -> Eigen::Matrix< typename TDerivedX::Scalar, TElement::template QuadratureType<QuadratureOrder, typename TDerivedX::Scalar>::kPoints, Eigen::Dynamic> |
Computes the inner product weights \( \mathbf{w}_{ge} \in \mathbb{R}^{|G^e| \times |E|} \) such that \( \int_\Omega \cdot d\Omega = \sum_e \sum_g w_{ge} \cdot \).
In other words, \( w_{ge} = w_g \det(J^e_g) \) where \( J^e_g \) is the Jacobian of the element map at the \( g^\text{{th} \) quadrature point and \( w_g \) is the \( g^\text{th} \) quadrature weight.
TElement | FEM element type |
QuadratureOrder | Quadrature order |
TDerivedE | Eigen matrix expression for element matrix |
TDerivedX | Eigen matrix expression for mesh nodal positions |
E | |# elem nodes| x |# elems| mesh element matrix |
X | |# dims| x |# nodes| mesh nodal position matrix |
|# quad.pts.| x |# elements|
matrix of quadrature weights multiplied by jacobian determinants at element quadrature points auto pbat::fem::MeshQuadratureWeights | ( | TMesh const & | mesh | ) | -> Eigen::Matrix< typename TMesh::ScalarType, TMesh::ElementType::template QuadratureType<QuadratureOrder, typename TMesh::ScalarType>:: kPoints, Eigen::Dynamic> |
Computes the inner product weights \( \mathbf{w}_{ge} \in \mathbb{R}^{|G^e| \times |E|} \) such that \( \int_\Omega \cdot d\Omega = \sum_e \sum_g w_{ge} \cdot \).
In other words, \( w_{ge} = w_g \det(J^e_g) \) where \( J^e_g \) is the Jacobian of the element map at the \( g^\text{{th} \) quadrature point and \( w_g \) is the \( g^\text{th} \) quadrature weight.
QuadratureOrder | Quadrature order |
TMesh | Mesh type |
mesh | FEM mesh |
|# quad.pts.|x|# elements|
matrix of quadrature weights multiplied by jacobian determinants at element quadrature points auto pbat::fem::MeshReferenceQuadraturePoints | ( | TIndex | nElements | ) |
Computes the element quadrature points in reference element space.
TElement | FEM element type |
QuadratureOrder | Quadrature order |
TScalar | Scalar type |
TIndex | Index type |
nElements | Number of elements |
|# ref. dims| x |# quad.pts. * # elems|
matrix of quadrature points in reference element space auto pbat::fem::MeshReferenceQuadraturePoints | ( | TMesh const & | mesh | ) |
Computes the element quadrature points in reference element space.
QuadratureOrder | Quadrature order |
TMesh | Mesh type |
mesh | FEM mesh |
|# ref. dims| x |# quad.pts. * # elems|
matrix of quadrature points in reference element space auto pbat::fem::ReferencePosition | ( | Eigen::MatrixBase< TDerivedX > const & | X, |
Eigen::MatrixBase< TDerivedx > const & | x, | ||
int const | maxIterations = 5, | ||
TScalar const | eps = 1e-10 ) -> Eigen::Vector<TScalar, TElement::kDims> |
Computes the reference position \( \xi \) such that \( x(\xi) = x \). This inverse problem is solved using Gauss-Newton iterations.
We need to solve the inverse problem
\[\min_\xi f(\xi) = 1/2 ||x(\xi) - X||_2^2 \]
to find the reference position \( \xi \) that corresponds to the domain position \( X \) in the element whose vertices are \( x \).
We use Gauss-Newton iterations on \( \min f \). This gives the iteration
\[dx^{k+1} = [H(\xi^k)]^{-1} J_x(\xi^k)^T (x(\xi^k) - X) \]
where \( H(\xi^k) = J_x(\xi^k)^T J_x(\xi^k) \).
TElement | FEM element type |
TDerivedX | Eigen matrix expression for reference positions |
TDerivedx | Eigen matrix expression for domain positions |
X | Reference positions \( \xi \) |
x | Domain positions \( x \) |
maxIterations | Maximum number of Gauss-Newton iterations |
eps | Convergence tolerance |
auto pbat::fem::ReferencePositions | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::MatrixBase< TDerivedX > const & | X, | ||
Eigen::DenseBase< TDerivedEg > const & | eg, | ||
Eigen::MatrixBase< TDerivedXg > const & | Xg, | ||
int | maxIterations = 5, | ||
typename TDerivedXg::Scalar | eps = 1e-10 ) -> Eigen::Matrix<typename TDerivedXg::Scalar, TElement::kDims, Eigen::Dynamic> |
TElement | FEM element type |
TDerivedE | Eigen matrix expression for elements |
TDerivedX | Eigen matrix expression for nodal positions |
TDerivedEg | Eigen matrix expression for indices of elements at evaluation points |
TDerivedXg | Eigen matrix expression for evaluation points in domain space |
E | |# elem nodes| x |# elems| element matrix |
X | |# dims| x |# nodes| nodal position matrix |
eg | |# eval.pts.| x 1 Indices of elements at evaluation points |
Xg | |# dims| x |# eval.pts.| evaluation points in domain space |
maxIterations | Maximum number of Gauss-Newton iterations |
eps | Convergence tolerance |
|# element dims| x n
matrix of reference positions associated with domain points auto pbat::fem::ReferencePositions | ( | Eigen::DenseBase< TDerivedEg > const & | Eg, |
Eigen::MatrixBase< TDerivedX > const & | X, | ||
Eigen::MatrixBase< TDerivedXg > const & | Xg, | ||
int | maxIterations = 5, | ||
TScalar | eps = 1e-10 ) -> Eigen::Matrix<TScalar, TElement::kDims, Eigen::Dynamic> |
Computes reference positions \( \xi \) such that \( X(\xi) = X_n \) for every point in \( \mathbf{X} \in \mathbb{R}^{d \times n} \).
TElement | FEM element type |
TDerivedEg | Eigen matrix expression for element indices at evaluation points |
TDerivedX | Eigen matrix expression for mesh node positions |
TDerivedXg | Eigen matrix expression for evaluation points |
Eg | |# elem nodes| x |# eval. pts.| matrix of element indices at evaluation points |
X | |# dims| x |# nodes| matrix of mesh node positions |
Xg | |# dims| x |# eval. pts.| matrix of evaluation points in domain space |
maxIterations | Maximum number of Gauss-Newton iterations |
eps | Convergence tolerance |
|# element dims| x n
matrix of reference positions associated with domain points X in corresponding elements E auto pbat::fem::ReferencePositions | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedEg > const & | eg, | ||
Eigen::MatrixBase< TDerivedXg > const & | Xg, | ||
int | maxIterations = 5, | ||
typename TMesh::ScalarType | eps = 1e-10 ) -> Eigen::Matrix<typename TMesh::ScalarType, TMesh::ElementType::kDims, Eigen::Dynamic> |
Computes reference positions \( \xi \) such that \( X(\xi) = X_n \) for every point in \( \mathbf{X} \in \mathbb{R}^{d \times n} \).
TMesh | FEM mesh type |
TDerivedEg | Eigen matrix expression for indices of elements |
TDerivedXg | Eigen matrix expression for evaluation points |
mesh | FEM mesh |
eg | |# eval.pts.| x 1 indices of elements at evaluation points |
Xg | |# dims| x |# eval.pts.| evaluation points in domain space |
maxIterations | Maximum number of Gauss-Newton iterations |
eps | Convergence tolerance |
|# element dims| x |# eval.pts.|
matrix of reference positions associated with domain points X in corresponding elements E E.size() == X.cols()
auto pbat::fem::ShapeFunctionGradients | ( | Eigen::MatrixBase< TDerivedE > const & | E, |
Eigen::MatrixBase< TDerivedX > const & | X ) -> Eigen::Matrix<TScalar, TElement::kNodes, Eigen::Dynamic> |
Computes nodal shape function gradients at each element quadrature points.
TElement | Element type |
Dims | Number of dimensions of the element |
QuadratureOrder | Quadrature order |
TDerivedE | Eigen dense expression type for element indices |
TDerivedX | Eigen dense expression type for element nodal positions |
TScalar | Floating point type, defaults to Scalar |
TIndex | Index type, defaults to coefficient type of TDerivedE |
E | |# elem nodes| x |# elements| array of element node indices |
X | |# dims| x |# nodes| array of mesh nodal positions |
|# elem nodes| x |# dims * # elem quad.pts. * # elems|
matrix of nodal shape function gradients auto pbat::fem::ShapeFunctionGradients | ( | TMesh const & | mesh | ) | -> Eigen::Matrix<typename TMesh::ScalarType, TMesh::ElementType::kNodes, Eigen::Dynamic> |
Computes nodal shape function gradients at each element quadrature points.
Order | Quadrature order |
TMesh | Mesh type |
mesh | FEM mesh |
|# element nodes| x |# dims * # quad.pts. * # elements|
matrix of shape function gradients auto pbat::fem::ShapeFunctionGradientsAt | ( | Eigen::DenseBase< TDerivedEg > const & | Eg, |
Eigen::MatrixBase< TDerivedX > const & | X, | ||
Eigen::MatrixBase< TDerivedXi > const & | Xi ) -> Eigen::Matrix<typename TDerivedXi::Scalar, TElement::kNodes, Eigen::Dynamic> |
Computes nodal shape function gradients at evaluation points Xi.
TElement | Element type |
Dims | Number of dimensions of the element |
TDerivedE | Eigen dense expression type for element indices |
TDerivedX | Eigen dense expression type for element nodal positions |
TDerivedXi | Eigen dense expression type for evaluation points |
Eg | |# elem nodes| x |# eval.pts.| array of element node indices at evaluation points |
X | |# dims| x |# elem nodes| array of element nodal positions |
Xi | |# dims| x |# eval.pts.| evaluation points in reference space |
|# elem nodes| x |# dims * # eval.pts.|
matrix of nodal shape function gradients auto pbat::fem::ShapeFunctionGradientsAt | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedEg > const & | eg, | ||
Eigen::MatrixBase< TDerivedXi > const & | Xi ) -> Eigen::Matrix<typename TMesh::ScalarType, TMesh::ElementType::kNodes, Eigen::Dynamic> |
Computes nodal shape function gradients at evaluation points Xg.
TMesh | Mesh type |
TDerivedEg | Eigen dense expression type for element indices |
TDerivedXi | Eigen dense expression type for evaluation points |
mesh | FEM mesh |
E | |# eval.pts.| array of elements associated with evaluation points |
Xi | |# dims|x|# eval.pts.| evaluation points |
|# element nodes| x |eg.size() * mesh.dims|
nodal shape function gradients at evaluation points Xi auto pbat::fem::ShapeFunctionMatrix | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes ) -> Eigen::SparseMatrix<TScalar, Eigen::RowMajor, typename TDerivedE::Scalar> |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the element quadrature points.
TElement | Element type |
QuadratureOrder | Quadrature order |
TScalar | Floating point type |
TIndex | Index type, defaults to Index |
E | |# nodes|x|# elements| array of element node indices |
nNodes | Number of nodes in the mesh |
|# elements * # quad.pts.| x |# nodes|
shape function matrix auto pbat::fem::ShapeFunctionMatrix | ( | TMesh const & | mesh | ) | -> Eigen::SparseMatrix<typename TMesh::ScalarType, Eigen::RowMajor, typename TMesh::IndexType> |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the element quadrature points.
QuadratureOrder | Quadrature order |
TMesh | Mesh type |
mesh | FEM mesh |
|# elements * # quad.pts.| x |# nodes|
shape function matrix auto pbat::fem::ShapeFunctionMatrixAt | ( | Eigen::DenseBase< TDerivedE > const & | E, |
TIndex | nNodes, | ||
Eigen::MatrixBase< TDerivedXi > const & | Xi ) -> Eigen::SparseMatrix<typename TDerivedXi::Scalar, Eigen::RowMajor, TIndex> |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the element quadrature points.
TElement | Element type |
TDerivedE | Eigen dense expression type for element indices |
TDerivedXi | Eigen dense expression type for quadrature points |
TIndex | Index type, defaults to coefficient type of TDerivedE |
E | |# nodes| x |# quad.pts.| array of element node indices |
nNodes | Number of nodes in the mesh |
Xi | |# dims| x |# quad.pts.| array of quadrature points in reference space |
|# quad.pts.| x |# nodes|
shape function matrix auto pbat::fem::ShapeFunctionMatrixAt | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedEg > const & | eg, | ||
Eigen::MatrixBase< TDerivedXi > const & | Xi ) -> Eigen::SparseMatrix<typename TMesh::ScalarType, Eigen::RowMajor, typename TMesh::IndexType> |
Constructs a shape function matrix \( \mathbf{N} \) for a given mesh, i.e. at the given evaluation points.
TMesh | Mesh type |
TDerivedEg | Eigen type |
TDerivedXi | Eigen type |
mesh | FEM mesh |
eg | |# quad.pts.| array of elements associated with quadrature points |
Xi | |# dims|x|# quad.pts.| array of quadrature points in reference space |
|# quad.pts.| x |# nodes|
shape function matrix auto pbat::fem::ShapeFunctionsAt | ( | Eigen::DenseBase< TDerivedXi > const & | Xi | ) | -> Eigen::Matrix<typename TDerivedXi::Scalar, TElement::kNodes, Eigen::Dynamic> |
Compute shape functions at the given reference positions.
TElement | Element type |
TDerivedXi | Eigen dense expression type for reference positions |
Xi | |# dims|x|# quad.pts.| evaluation points |
|# element nodes| x |Xi.cols()|
matrix of nodal shape functions at reference points Xi PBAT_API math::linalg::EEigenvalueFilter pbat::fem::ToEigenvalueFilter | ( | EHyperElasticSpdCorrection | mode | ) |
Convert EHyperElasticSpdCorrection to EEigenvalueFilter.
mode | SPD correction mode |
void pbat::fem::ToElementElasticity | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::DenseBase< TDerivedwg > const & | wg, | ||
Eigen::MatrixBase< TDerivedGNeg > const & | GNeg, | ||
Eigen::DenseBase< TDerivedmug > const & | mug, | ||
Eigen::DenseBase< TDerivedlambdag > const & | lambdag, | ||
Eigen::MatrixBase< TDerivedx > const & | x, | ||
Eigen::PlainObjectBase< TDerivedUg > & | Ug, | ||
Eigen::PlainObjectBase< TDerivedGg > & | Gg, | ||
Eigen::PlainObjectBase< TDerivedHg > & | Hg, | ||
int | eFlags = EElementElasticityComputationFlags::Potential, | ||
EHyperElasticSpdCorrection | eSpdCorrection = EHyperElasticSpdCorrection::None ) |
Compute element elasticity and its derivatives at the given shape.
TElement | Element type |
Dims | Number of spatial dimensions |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedlameg | Type of the Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
TDerivedUg | Type of the elastic potentials at quadrature points |
TDerivedGg | Type of the element gradient vectors |
TDerivedHg | Type of the element hessian matrices |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
Ug | Optional |# quad.pts.| x 1 elastic potentials at quadrature points |
Gg | Optional |# dims * # elem nodes| x |# quad.pts.| element elastic potential gradients at quadrature points. |
Hg | Optional |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points. |
eFlags | Flags for the computation (see EElementElasticityComputationFlags) |
eSpdCorrection | SPD correction mode (see EHyperElasticSpdCorrection) |
|
inline |
Compute element elasticity using mesh.
THyperElasticEnergy | Hyper elastic energy type |
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedmug | Type of the first Lame coefficients at quadrature points |
TDerivedlambdag | Type of the second Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
TDerivedUg | Type of the elastic potentials at quadrature points |
TDerivedGg | Type of the element gradient vectors at quadrature points |
TDerivedHg | Type of the element hessian matrices at quadrature points |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
Ug | |# quad.pts.| x 1 elastic potentials at quadrature points |
Gg | Optional |# dims * # elem nodes| x |# quad.pts.| element elastic potential gradients at quadrature points. |
Hg | Optional |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points. |
eFlags | Flags for the computation (see EElementElasticityComputationFlags) |
eSpdCorrection | SPD correction mode (see EHyperElasticSpdCorrection) |
|
inline |
Compute a matrix of horizontally stacked element mass matrices for all quadrature points.
TElement | Element type |
TN | Type of the shape function vector at the quadrature point |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
TDerivedOut | Type of the output matrix |
Neg | |# nodes per element| x |# quad.pts.| shape function matrix at all quadrature points |
wg | |# quad.pts.| x 1 quadrature weights (including Jacobian determinant) |
rhog | |# quad.pts.| x 1 mass density at quadrature points |
M | |# elem nodes| x |# elem nodes * # quad.pts.| out matrix of stacked element mass matrices |
|
inline |
Compute the gradient vector into existing output.
TElement | Element type |
Dims | Number of spatial dimensions |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedmug | Type of the first Lame coefficients at quadrature points |
TDerivedlambdag | Type of the second Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
TDerivedGg | Type of the element gradient vectors at quadrature points |
TDerivedG | Type of the output gradient vector |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
G | |# dims * # nodes| output gradient vector |
void pbat::fem::ToHyperElasticGradient | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedeg > const & | eg, | ||
Eigen::MatrixBase< TDerivedGg > const & | Gg, | ||
Eigen::PlainObjectBase< TDerivedG > & | G ) |
Compute the gradient vector into existing output.
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGg | Type of the element gradient vectors |
TDerivedOut | Type of the output vector |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
G | |# dims * # nodes| output gradient vector |
|
inline |
Compute the gradient vector into existing output.
TMesh | Type of the mesh |
THyperElasticEnergy | Hyper elastic energy type |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedGNeg | Type of the shape function gradients at quadrature points |
TDerivedmug | Type of the first Lame coefficients at quadrature points |
TDerivedlambdag | Type of the second Lame coefficients at quadrature points |
TDerivedx | Type of the deformed nodal positions |
TDerivedGg | Type of the element gradient vectors at quadrature points |
TDerivedG | Type of the output gradient vector |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights |
GNeg | |# nodes per element| x |# dims * # quad.pts.| shape function gradients |
mug | |# quad.pts.| x 1 first Lame coefficients at quadrature points |
lambdag | |# quad.pts.| x 1 second Lame coefficients at quadrature points |
x | |# dims * # nodes| x 1 deformed nodal positions |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
G | |# dims * # nodes| output gradient vector |
|
inline |
Compute the gradient vector using mesh (updates existing output)
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedGg | Type of the element gradient vectors |
TDerivedOut | Type of the output vector |
mesh | Mesh containing element connectivity and node positions |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Gg | |# dims * # elem nodes| x |# quad.pts.| array of element elastic gradient vectors at quadrature points |
G | |# dims * # nodes| output gradient vector |
void pbat::fem::ToHyperElasticHessian | ( | Eigen::DenseBase< TDerivedE > const & | E, |
Eigen::Index | nNodes, | ||
Eigen::DenseBase< TDerivedHg > const & | Hg, | ||
math::linalg::SparsityPattern< typename TDerivedE::Scalar, Options > const & | sparsity, | ||
Eigen::SparseCompressedBase< TDerivedH > & | H ) |
Construct the hessian matrix using mesh with sparsity pattern.
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedHg | Type of the element hessian matrices |
Options | Storage options for the matrix |
TDerivedH | Type of the output sparse matrix |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points |
sparsity | Sparsity pattern for the hessian matrix |
H | Output sparse matrix for the hessian |
void pbat::fem::ToHyperElasticHessian | ( | TMesh const & | mesh, |
Eigen::DenseBase< TDerivedHg > const & | Hg, | ||
math::linalg::SparsityPattern< typename TMesh::IndexType, Options > const & | sparsity, | ||
Eigen::SparseCompressedBase< typename TDerivedHg::Scalar > & | H ) |
Construct the hessian matrix using mesh with sparsity pattern.
TMesh | Type of the mesh |
TDerivedHg | Type of the element hessian matrices |
Options | Storage options for the matrix |
TDerivedH | Type of the output sparse matrix |
mesh | Mesh containing element connectivity and node positions |
Hg | |# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| element elastic hessian matrices at quadrature points |
sparsity | Sparsity pattern for the hessian matrix |
H | Output sparse matrix for the hessian |
|
inline |
Compute lumped mass matrix's diagonal vector into existing output vector.
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
TDerivedNeg | Type of the shape functions at quadrature points |
TDerivedOut | Type of the output vector |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights (including Jacobian determinants) |
rhog | |# quad.pts.| x 1 vector of density at quadrature points |
Neg | |# nodes per element| x |# quad.pts.| shape functions at quadrature points |
dims | Dimensionality of the image of the FEM function space |
m | Output vector of lumped masses |# nodes * dims| x 1 |
|
inline |
Compute lumped mass vector from precomputed element mass matrices into existing output vector.
TElement | Element type |
Dims | Number of spatial dimensions |
TDerivedE | Type of the element matrix |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
TDerivedOut | Type of the output vector |
E | |# nodes per element| x |# elements| matrix of mesh elements |
nNodes | Number of mesh nodes |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Meg | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
m | Output vector of lumped masses |# nodes * dims| x 1 |
|
inline |
Compute lumped mass vector using mesh and precomputed element mass matrices into existing output vector.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedwg | Type of the quadrature weights |
TDerivedrhog | Type of the density at quadrature points |
TDerivedNeg | Type of the shape functions at quadrature points |
TDerivedOut | Type of the output vector |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
wg | |# quad.pts.| x 1 vector of quadrature weights (including Jacobian determinants) |
rhog | |# quad.pts.| x 1 vector of density at quadrature points |
Neg | |# nodes per element| x |# quad.pts.| shape functions at quadrature points |
dims | Dimensionality of the image of the FEM function space |
m | Output vector of lumped masses |# nodes * dims| x 1 |
|
inline |
Compute lumped mass vector using mesh and precomputed element mass matrices into existing output vector.
TMesh | Type of the mesh |
TDerivedeg | Type of the element indices at quadrature points |
TDerivedMe | Type of the precomputed element mass matrices |
TDerivedOut | Type of the output vector |
mesh | The finite element mesh |
eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
Meg | |# nodes per element| x |# nodes per element * # quad.pts.| precomputed element mass matrices |
dims | Dimensionality of the image of the FEM function space |
m | Output vector of lumped masses |# nodes * dims| x 1 |
void pbat::fem::ToMeshQuadratureWeights | ( | Eigen::MatrixBase< TDerivedDetJe > & | detJeThenWg | ) |
Computes the inner product weights \( \mathbf{w}_{ge} \in \mathbb{R}^{|G^e| \times |E|} \) such that \( \int_\Omega \cdot d\Omega = \sum_e \sum_g w_{ge} \cdot \).
In other words, \( w_{ge} = w_g \det(J^e_g) \) where \( J^e_g \) is the Jacobian of the element map at the \( g^\text{{th} \) quadrature point and \( w_g \) is the \( g^\text{th} \) quadrature weight.
QuadratureOrder | Quadrature order |
TMesh | Mesh type |
mesh | FEM mesh |
detJeThenWg | |# quad.pts.| x |# elements| matrix of jacobian determinants at element quadrature points |
|# quad.pts.|x|# elements|
matrix of quadrature weights multiplied by jacobian determinants at element quadrature points