Rodin/Variational/QuadratureRule.h file

Quadrature rule classes for numerical integration.

This file defines specializations of Rodin::Variational::QuadratureRule for integrating functions and shape-function expressions over mesh polytopes.

Mathematical foundation

A quadrature rule approximates an integral by a weighted sum:

\[ \int_K f(x)\,dx \approx \sum_{q=1}^{n_q} w_q\,f(x_q), \]

where:

  • $ x_q $ are quadrature points,
  • $ w_q $ are quadrature weights,
  • $ n_q $ is the number of quadrature points.

In Rodin, the reference quadrature formula lives in the Rodin::QF module, while the mapped geometric points live in Rodin::Geometry::PolytopeQuadrature. This separation allows:

  • canonical reuse of reference quadrature formulas,
  • mesh-owned caching of mapped geometric quadrature points,
  • simpler and more efficient integrator implementations.

Usage

Typical usage is:

auto qr = QuadratureRule(integrand);
qr.setPolytope(cell);
const auto value = qr.compute();

For local element integrators, setPolytope binds the integrator to a concrete polytope, selects an appropriate quadrature formula, and retrieves the corresponding cached Rodin::Geometry::PolytopeQuadrature from the mesh.

Namespaces

namespace Rodin
The Rodin library for finite element methods and shape optimization.
namespace Rodin::Variational
Module which provides the necessary tools for constructing and solving variational problems.

Classes

template<class FunctionDerived>
class Rodin::Variational::QuadratureRule<FunctionBase<FunctionDerived>>
Quadrature rule for integrating general functions on mesh polytopes.
template<class FES, class Data>
class Rodin::Variational::QuadratureRule<GridFunction<FES, Data>>
Integration of a scalar grid function over a mesh region.
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
class Rodin::Variational::QuadratureRule<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>>
Local bilinear-form quadrature for the dot product of a trial and a test shape function.
template<class NestedDerived, class FES>
class Rodin::Variational::QuadratureRule<ShapeFunctionBase<NestedDerived, FES, TestSpace>>
Local linear-form quadrature for a test shape function expression.