Rodin::QF::Centroid class final

Single-point centroid quadrature formula.

This class implements a quadrature formula that uses a single integration point located at the centroid (barycenter) of the reference polytope. The quadrature is exact for polynomials up to degree 1 (linear functions).

Mathematical Foundation

The centroid rule approximates integrals as:

\[ \int_K f(x) \, dx \approx w \cdot f(c) \]

where $ c $ is the centroid of the polytope $ K $ and $ w $ is the weight equal to the volume (measure) of the reference element.

Reference Element Centroids

The quadrature points are located at:

  • Point: $ (0) $
  • Segment: $ (0.5) $
  • Triangle: $ (1/3, 1/3) $
  • Quadrilateral: $ (0.5, 0.5) $
  • Tetrahedron: $ (0.25, 0.25, 0.25) $
  • Wedge: $ (1/3, 1/3, 0.5) $

Reference Element Weights

The weights correspond to the reference element measures:

  • Point: $ 1 $
  • Segment: $ 1 $
  • Triangle: $ 1/2 $
  • Quadrilateral: $ 1 $
  • Tetrahedron: $ 1/6 $
  • Wedge: $ 1/2 $

This is the simplest possible quadrature rule for each geometry type and is computationally efficient when low accuracy is acceptable.

Base classes

class QuadratureFormulaBase
Abstract base class for quadrature formulas.

Public types

using Parent = QuadratureFormulaBase
Parent class type.

Constructors, destructors, conversion operators

Centroid(Geometry::Polytope::Type g) constexpr
Constructs a centroid quadrature formula for the given geometry.

Public functions

auto getSize() const -> size_t override
Gets the number of quadrature points (always 1).
auto getPoint(size_t i) const -> const Math::SpatialVector<Real>& override
Gets the single quadrature point coordinates (the centroid).
auto getWeight(size_t i) const -> Real override
Gets the weight of the single quadrature point.
auto copy() const -> Centroid* override noexcept
Creates a copy of this quadrature formula.

Function documentation

Rodin::QF::Centroid::Centroid(Geometry::Polytope::Type g) constexpr

Constructs a centroid quadrature formula for the given geometry.

Parameters
g Geometry type (e.g., Triangle, Quadrilateral, Tetrahedron, etc.)

size_t Rodin::QF::Centroid::getSize() const override

Gets the number of quadrature points (always 1).

Returns Always returns 1 (single point quadrature)

const Math::SpatialVector<Real>& Rodin::QF::Centroid::getPoint(size_t i) const override

Gets the single quadrature point coordinates (the centroid).

Parameters
i Index of the quadrature point (must be 0)
Returns Reference to the centroid coordinates in reference space

Real Rodin::QF::Centroid::getWeight(size_t i) const override

Gets the weight of the single quadrature point.

Parameters
i Index of the quadrature point (must be 0)
Returns Weight equal to the reference element measure

Centroid* Rodin::QF::Centroid::copy() const override noexcept

Creates a copy of this quadrature formula.

Returns Pointer to a new Centroid instance (caller takes ownership)