template<class LHSDerived, class RHSDerived, class Range, class Mesh>
Rodin::Variational::QuadratureRule<ShapeFunctionBase<Dot<FunctionBase<LHSDerived>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<Range, Mesh>, TestSpace>, P1<Range, Mesh>, TestSpace>>, P1<Range, Mesh>, TestSpace>> class

Integration of the Dot product of some coefficient function and a P1 ShapeFunction.

This class represents the CTAD for the expression:

\[ \int f \cdot v \ dx \: , \]

where $ v \in \mathbb{P}_1 $ .

Judgement

The following judgement specifies that the expression is a well formed type of QuadratureRule.

\[ \dfrac {\vdash \int f \cdot v \ dx : \texttt{QuadratureRule}} {\vdash v : \mathbb{P}_1} \]

Base classes

template<class Number>
class LinearFormIntegratorBase<FormLanguage::Traits<ShapeFunctionBase<Dot<FunctionBase<LHSDerived>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<Range, Mesh>, TestSpace>>>>> ::ScalarType>
Base class for linear form integrators.

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 RHSDerived, class Range, class Mesh>
QuadratureRule* Rodin::Variational::QuadratureRule<ShapeFunctionBase<Dot<FunctionBase<LHSDerived>, ShapeFunctionBase<ShapeFunction<RHSDerived, P1<Range, Mesh>, TestSpace>, P1<Range, Mesh>, TestSpace>>, P1<Range, Mesh>, TestSpace>><LHSDerived, RHSDerived, Range, Mesh>::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.