template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
Rodin::Variational::QuadratureRule<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>> class

Approximation of the integral of the the dot product between a trial shape function and a test shape function.

Represents the quadrature rule approximation of an integral:

\[ \int_{\mathcal{R}_h} \mathrm{IntegrandType} \ dx \approx \sum_{i = 1}^{n} w_i \ \mathrm{IntegrandType} (x_i) \]

where $ \mathcal{R}_h $ is some region of the mesh $ \mathcal{T}_h $ , the quadrature point $ x_i $ has an associated weight $ w_i $ and $ \mathrm{IntegrandType}(x_i) $ is the value of the integrand at the quadrature point.

Derived classes

template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
class BoundaryIntegral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>> final
Integration of the dot product of a trial and test operators.
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
class FaceIntegral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>> final
Integration of the dot product of a trial and test operators.
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
class Integral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>> final
Integration of the dot product of a trial and test operators.
template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
class InterfaceIntegral<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>> final
Integration of the dot product of a trial and test operators.

Public functions

auto copy() const -> QuadratureRule* override noexcept
Copies the object and returns a non-owning pointer to the copied object.

Function documentation

template<class LHSDerived, class TrialFES, class RHSDerived, class TestFES>
QuadratureRule* Rodin::Variational::QuadratureRule<Dot<ShapeFunctionBase<LHSDerived, TrialFES, TrialSpace>, ShapeFunctionBase<RHSDerived, TestFES, TestSpace>>><LHSDerived, TrialFES, RHSDerived, TestFES>::copy() const override noexcept

Copies the object and returns a non-owning pointer to the copied object.

Returns Non-owning pointer to the copied object.