|
template<class TDerivedV, class TDerivedC> |
| Mesh (Eigen::MatrixBase< TDerivedV > const &V, Eigen::DenseBase< TDerivedC > const &C) |
| Constructs a finite element mesh given some input geometric mesh.
|
|
template<class TDerivedV, class TDerivedC> |
void | Construct (Eigen::MatrixBase< TDerivedV > const &V, Eigen::DenseBase< TDerivedC > const &C) |
| Constructs a finite element mesh given some input geometric mesh.
|
|
template<int QuadratureOrder> |
auto | QuadraturePoints () const -> Eigen::Matrix< ScalarType, Eigen::Dynamic, Eigen::Dynamic > |
| Compute quadrature points in domain space on this mesh.
|
|
template<int QuadratureOrder, int kPoints = TElement::template QuadratureType<QuadratureOrder, ScalarType>::kPoints> |
auto | QuadratureWeights () const -> Eigen::Vector< ScalarType, kPoints > |
| Obtain quadrature weights on the reference element of this mesh.
|
|
void | Serialize (io::Archive &archive) const |
| Serialize to HDF5 group.
|
|
void | Deserialize (io::Archive const &archive) |
| Deserialize from HDF5 group.
|
|
template<int QuadratureOrder> |
Eigen::Matrix< TScalar, Eigen::Dynamic, Eigen::Dynamic > | QuadraturePoints () const |
|
template<int QuadratureOrder, int kPoints> |
Eigen::Vector< TScalar, kPoints > | QuadratureWeights () const |
|
template<CElement TElement, int Dims, common::CFloatingPoint TScalar = Scalar, common::CIndex TIndex = Index>
struct pbat::fem::Mesh< TElement, Dims, TScalar, TIndex >
A generic stateful finite element mesh representation.
- Template Parameters
-
TElement | Type satisfying concept CElement |
Dims | Embedding dimensions of the mesh |
TScalar | Floating point type, defaults to Scalar |
TIndex | Index type, defaults to Index |
template<CElement TElement, int Dims, common::CFloatingPoint TScalar, common::CIndex TIndex>
template<class TDerivedV, class TDerivedC>
pbat::fem::Mesh< TElement, Dims, TScalar, TIndex >::Mesh |
( |
Eigen::MatrixBase< TDerivedV > const & | V, |
|
|
Eigen::DenseBase< TDerivedC > const & | C ) |
|
inline |
Constructs a finite element mesh given some input geometric mesh.
- Warning
- The cells of the input mesh should list its vertices in Lagrange order.
- Template Parameters
-
TDerivedV | Type of the vertex positions matrix |
TDerivedC | Type of the cell vertex indices matrix |
- Parameters
-
V | Dims x |# vertices| matrix of vertex positions |
C | |Element::AffineBase::Vertices| x |# cells| matrix of cell vertex indices into V |
template<CElement TElement, int Dims, common::CFloatingPoint TScalar, common::CIndex TIndex>
template<class TDerivedV, class TDerivedC>
void pbat::fem::Mesh< TElement, Dims, TScalar, TIndex >::Construct |
( |
Eigen::MatrixBase< TDerivedV > const & | V, |
|
|
Eigen::DenseBase< TDerivedC > const & | C ) |
|
inline |
Constructs a finite element mesh given some input geometric mesh.
- Warning
- The cells of the input mesh should list its vertices in Lagrange order.
- Template Parameters
-
TDerivedV | Type of the vertex positions matrix |
TDerivedC | Type of the cell vertex indices matrix |
- Parameters
-
V | Dims x |# vertices| matrix of vertex positions |
C | |Element::AffineBase::Vertices| x |# cells| matrix of cell vertex indices into V |
template<CElement TElement, int Dims, common::CFloatingPoint TScalar = Scalar, common::CIndex TIndex = Index>
template<int QuadratureOrder, int kPoints = TElement::template QuadratureType<QuadratureOrder, ScalarType>::kPoints>
Obtain quadrature weights on the reference element of this mesh.
- Template Parameters
-
QuadratureOrder | Quadrature order |
- Parameters
-
kPoints | Number of quadrature points, defaults to Element::template QuadratureType<QuadratureOrder>::kPoints . This template parameter should not be specified, it is there to enhance readability. |
- Returns
|# element quad.pts.|
vector of quadrature weights