PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
Loading...
Searching...
No Matches
pbat::fem::CElement Concept Reference

Reference finite element. More...

#include <Concepts.h>

Concept definition

template<class T>
concept pbat::fem::CElement = requires(T t)
{
typename T::AffineBaseType;
{T::bHasConstantJacobian}->std::convertible_to<bool>;
typename T::template QuadratureType<1, Scalar>;
{T::kOrder}->std::convertible_to<int>;
{T::kDims}->std::convertible_to<int>;
{T::kNodes}->std::convertible_to<int>;
requires common::CContiguousIndexRange<decltype(T::Coordinates)>;
requires common::CContiguousIndexRange<decltype(T::Vertices)>;
{t.N(Vector<T::kDims>{})}->std::convertible_to<Vector<T::kNodes>>;
{t.GradN(Vector<T::kDims>{})}->std::convertible_to<Matrix<T::kNodes, T::kDims>>;
}
Reference finite element.
Definition Concepts.h:63
Eigen::Vector< Scalar, N > Vector
Fixed-size vector type.
Definition Aliases.h:24
Eigen::Matrix< Scalar, Rows, Cols > Matrix
Fixed-size matrix type.
Definition Aliases.h:31

Detailed Description

Reference finite element.

Example type TElement satisfying concept CElement

template <int Order>
struct TElement
{
using AffineBaseType = TElement<1>;
static int constexpr kOrder;
static int constexpr kDims;
static int constexpr kNodes;
static std::array<int, kNodes * kDims> constexpr Coordinates;
static std::array<int, AffineBaseType::kNodes> constexpr Vertices;
static bool constexpr bHasConstantJacobian;
template <int PolynomialOrder>
using QuadratureType = math::SymmetricSimplexPolynomialQuadratureRule<kDims,
PolynomialOrder>;
template <class TDerived, class TScalar = typename TDerived::Scalar>
static Eigen::Vector<TScalar, kNodes>
N(Eigen::DenseBase<TDerived> const& X);
static Matrix<kNodes, kDims> GradN(Vector<kDims> const& X);
};
typename detail::SymmetricSimplexPolynomialQuadratureRule< Dims, Order, TScalar > SymmetricSimplexPolynomialQuadratureRule
Symmetric quadrature rule for reference simplices in dimensions .
Definition SymmetricQuadratureRules.h:96
Note
  • Divide Coordinates by kOrder to obtain actual coordinates in the reference element
  • Vertices are indices into nodes [0,kNodes-1] revealing vertices of the element
  • bHasConstantJacobian is true if the Jacobian is constant over the element
  • QuadratureType proposes a quadrature rule suitable for integrating over the element
  • N computes the nodal shape functions evaluated at the given reference point X
  • GradN computes the gradient of the nodal shape functions evaluated at the given reference point X
Template Parameters
TElement type