PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
|
A matrix-free representation of a finite element load vector \( \mathbf{f}_i = \int_\Omega \mathbf{f}(X) \phi_i(X) \) under Galerkin projection. More...
#include <LoadVector.h>
Public Types | |
using | SelfType = LoadVector<TMesh, QuadratureOrder> |
Self type. | |
using | MeshType = TMesh |
Mesh type. | |
using | ElementType = typename TMesh::ElementType |
Element type. | |
using | QuadratureRuleType |
Quadrature rule type. | |
Public Member Functions | |
template<class TDerived> | |
LoadVector (MeshType const &mesh, Eigen::Ref< MatrixX const > const &detJe, Eigen::DenseBase< TDerived > const &fe, int dims=1) | |
Construct a new LoadVector object. | |
SelfType & | operator= (SelfType const &)=delete |
VectorX | ToVector () const |
Transforms this element-wise load vector representation into a global load vector. | |
template<class TDerived> | |
void | SetLoad (Eigen::DenseBase< TDerived > const &fe) |
Set the external loading. | |
void | CheckValidState () const |
Check if the state of this load vector is valid. | |
Public Attributes | |
MeshType const & | mesh |
The finite element mesh. | |
MatrixX | fe |
|dims|x|# elements| piecewise element constant load | |
MatrixX | N |
Eigen::Ref< MatrixX const > | detJe |
int | dims |
Dimensionality of external loading. Should have dims >= 1 . | |
Static Public Attributes | |
static int constexpr | kOrder = ElementType::kOrder |
Polynomial order of the load vector. | |
static int constexpr | kQuadratureOrder = QuadratureOrder |
Quadrature order. | |
A matrix-free representation of a finite element load vector \( \mathbf{f}_i = \int_\Omega \mathbf{f}(X) \phi_i(X) \) under Galerkin projection.
TMesh | Type satisfying concept CMesh |
QuadratureOrder | Quadrature order for integrating the load vector |
using pbat::fem::LoadVector< TMesh, QuadratureOrder >::QuadratureRuleType |
Quadrature rule type.
|
inline |
Construct a new LoadVector object.
TDerived | Eigen dense expression type |
mesh | Finite element mesh |
detJe | |# quad.pts.|x|# elements| affine element jacobian determinants at quadrature points |
fe | |dims|x|# elements| piecewise element constant load, or |dims|x1 constant load |
dims | Dimensionality of image of FEM function space |
dims >= 1
and fe.cols() == 1 || fe.cols() == mesh.E.cols()
and fe.rows() == dims
|
inline |
Set the external loading.
TDerived | Eigen dense expression type |
fe | |dims|x|# elements| piecewise element constant load, or |dims|x1 constant load |
|
inline |
Transforms this element-wise load vector representation into a global load vector.
Eigen::Ref<MatrixX const> pbat::fem::LoadVector< TMesh, QuadratureOrder >::detJe |
|# element quadrature points|x|# elements|
matrix of jacobian determinants at element quadrature points
MatrixX pbat::fem::LoadVector< TMesh, QuadratureOrder >::N |
|ElementType::kNodes|x|# elements|
integrated element shape functions. See (IntegratedShapeFunctions()).