template<size_t K, class Scalar>
Rodin::Variational::H1Element class

Degree k H1-conforming Lagrange element.

Template parameters
K Polynomial degree

Continuous H1-conforming piecewise polynomial (degree k) scalar Lagrange element.

H1Element provides a high-order H1-conforming finite element implementation with arbitrary polynomial degree K. The element uses Lagrange nodal basis functions satisfying $ \phi_i(x_j) = \delta_{ij} $ .

The H1Element provides a k-th order finite element with Lagrange basis functions. This is a high-order-stable H1-conforming element family.

Mathematical Properties

  • DOF count: Depends on geometry and polynomial degree K
    • Segment: K+1
    • Triangle: (K+1)(K+2)/2
    • Quadrilateral: (K+1)²
    • Tetrahedron: (K+1)(K+2)(K+3)/6
    • Wedge: (K+1)·(K+1)(K+2)/2
  • Basis functions: Lagrange polynomials of degree K satisfying the Lagrange property: $ \phi_i(x_j) = \delta_{ij} $ where $ x_j $ are the Lagrange nodes
  • Gradient: Polynomial of degree K-1
  • Partition of unity: $ \sum_{i=1}^{n} \phi_i(x) = 1 $ for all $ x $
  • Continuity: Global C⁰ continuity across element interfaces (H¹-conforming)

Convergence Rates

For smooth solutions, H1 elements of degree k provide k-th order convergence:

  • L² norm: $ \|u - u_h\|_{L^2} = O(h^{k+1}) $
  • H¹ norm (energy): $ \|u - u_h\|_{H^1} = O(h^k) $

Basis Function Construction

  • 1D (Segment): Classical Lagrange interpolation $ L_i^K(x) = \prod_{j \neq i} \frac{x - x_j}{x_i - x_j} $
  • 2D/3D Simplices: Barycentric Lagrange polynomials $ L_n^K(\lambda) = \prod_{m=0}^{n-1} \frac{K\lambda - m}{m+1} $
  • Tensor Products: For quadrilaterals and wedges, tensor product of lower-dimensional bases

Usage Example

// H1 element of degree 3 on triangle - 10 DOFs
RealH1Element<3> h1_tri(Polytope::Type::Triangle);
std::cout << h1_tri.getCount() << std::endl;  // Output: 10

// Evaluate basis function at a point
Math::Vector<Real> pt{{0.25, 0.25}};
Real value = h1_tri.getBasis(0)(pt);

// Get gradient of basis function
auto grad = h1_tri.getBasis(0).getGradient()(pt);

Base classes

template<class Derived>
class FiniteElementBase<H1Element<K, Scalar>>
Base class for finite elements.

Public types

class BasisFunction
Represents a Lagrange basis function of degree K.
class LinearForm
Evaluates linear form at a given function/vector field.
using Parent = FiniteElementBase<H1Element<K, Scalar>>
Parent class.
using RangeType = ScalarType
Type of range.

Public static functions

static auto getNodes(Geometry::Polytope::Type g) -> const std::vector<Math::SpatialPoint>&
Builds the high-order stable nodes for the element geometry.

Constructors, destructors, conversion operators

H1Element()
Default constructor. Creates an H1 element on a Point geometry.
H1Element(Geometry::Polytope::Type geometry)
Constructs an H1 element for the specified geometry.
H1Element(const H1Element& other) constexpr
Copy constructor.
H1Element(H1Element&& other) constexpr
Move constructor.

Public functions

auto operator=(const H1Element& other) -> H1Element& constexpr
Copy assignment operator.
auto operator=(H1Element&& other) -> H1Element& constexpr
Move assignment operator.
auto getCount() const -> size_t constexpr
Gets the number of degrees of freedom in the finite element.
auto getNode(size_t i) const -> const Math::SpatialPoint& constexpr
Gets the spatial coordinates of the i-th Lagrange node.
auto getLinearForm(size_t i) const -> const LinearForm&
Gets the linear form (evaluation functional) for the i-th DOF.
auto getBasis(size_t i) const -> const BasisFunction&
Gets the i-th basis function.
auto getOrder() const -> size_t constexpr
Gets the total maximum polynomial order of the element.
template<class Archive>
void serialize(Archive& ar, const unsigned int version)
Serializes the element (for boost::serialization).

Function documentation

template<size_t K, class Scalar>
static const std::vector<Math::SpatialPoint>& Rodin::Variational::H1Element<K, Scalar>::getNodes(Geometry::Polytope::Type g)

Builds the high-order stable nodes for the element geometry.

Constructs the positions of DOF nodes based on the geometry type and polynomial degree. Uses:

  • Segment: Gauss-Lobatto-Legendre (GLL) nodes on [0,1]
  • Quadrilateral: Tensor product of 1D GLL nodes
  • Triangle: Fekete nodes on reference triangle
  • Tetrahedron: Fekete nodes on reference tetrahedron
  • Wedge: Tensor product of triangle Fekete nodes with 1D GLL nodes

template<size_t K, class Scalar>
Rodin::Variational::H1Element<K, Scalar>::H1Element(Geometry::Polytope::Type geometry)

Constructs an H1 element for the specified geometry.

Parameters
geometry Type of element geometry (Segment, Triangle, Quadrilateral, etc.)

Initializes the element and builds the Lagrange nodes for the specified geometry. The number and positions of nodes depend on both the geometry and polynomial degree K.

template<size_t K, class Scalar>
Rodin::Variational::H1Element<K, Scalar>::H1Element(const H1Element& other) constexpr

Copy constructor.

Parameters
other Element to copy from

template<size_t K, class Scalar>
Rodin::Variational::H1Element<K, Scalar>::H1Element(H1Element&& other) constexpr

Move constructor.

Parameters
other Element to move from

template<size_t K, class Scalar>
H1Element& Rodin::Variational::H1Element<K, Scalar>::operator=(const H1Element& other) constexpr

Copy assignment operator.

Parameters
other Element to copy from
Returns Reference to this element

template<size_t K, class Scalar>
H1Element& Rodin::Variational::H1Element<K, Scalar>::operator=(H1Element&& other) constexpr

Move assignment operator.

Parameters
other Element to move from
Returns Reference to this element

template<size_t K, class Scalar>
size_t Rodin::Variational::H1Element<K, Scalar>::getCount() const constexpr

Gets the number of degrees of freedom in the finite element.

Returns Number of degrees of freedom

The DOF count depends on geometry and polynomial degree K:

  • Segment: K+1
  • Triangle: (K+1)(K+2)/2
  • Quadrilateral: (K+1)²
  • Tetrahedron: (K+1)(K+2)(K+3)/6
  • Wedge: (K+1)·(K+1)(K+2)/2

template<size_t K, class Scalar>
const Math::SpatialPoint& Rodin::Variational::H1Element<K, Scalar>::getNode(size_t i) const constexpr

Gets the spatial coordinates of the i-th Lagrange node.

Parameters
i Node index (0 to getCount()-1)
Returns Reference coordinates of the node (typically in [0,1]^d)

Lagrange nodes are the interpolation points where the basis functions satisfy φ_i(x_j) = δ_ij. Their positions are determined by the polynomial degree and geometry type.

template<size_t K, class Scalar>
const LinearForm& Rodin::Variational::H1Element<K, Scalar>::getLinearForm(size_t i) const

Gets the linear form (evaluation functional) for the i-th DOF.

Parameters
i DOF index (0 to getCount()-1)
Returns Linear form object that evaluates functions at the i-th node

template<size_t K, class Scalar>
const BasisFunction& Rodin::Variational::H1Element<K, Scalar>::getBasis(size_t i) const

Gets the i-th basis function.

Parameters
i Basis function index (0 to getCount()-1)
Returns Basis function object

The basis functions satisfy the Lagrange property and form a partition of unity.

template<size_t K, class Scalar>
size_t Rodin::Variational::H1Element<K, Scalar>::getOrder() const constexpr

Gets the total maximum polynomial order of the element.

Returns Polynomial degree K

template<size_t K, class Scalar> template<class Archive>
void Rodin::Variational::H1Element<K, Scalar>::serialize(Archive& ar, const unsigned int version)

Serializes the element (for boost::serialization).

Parameters
ar Archive to serialize to/from
version Serialization version