template<class TrialFES, class TestFES, class Operator>
Rodin::Variational::BilinearForm class final

Speciallization of BilinearForm for a matrix type.

Template parameters
TrialFES Type of trial finite element space
TestFES Type of test finite element space

Represents a bilinear form on a trial and test space.

This specialization aids in the construction of a $ m \times n $ matrix $ A $ , which is associated to the bilinear form. Here, $ n $ represents the size (total number of degrees-of-freedom) of the trial space, and $ m $ represents the size of the test space.

An object of type BilinearForm represents a bilinear map:

\[ \begin{aligned} a : U \times V &\rightarrow \mathbb{R}\\ (u, v) &\mapsto a(u, v) \end{aligned} \]

where $ U $ and $ V $ are finite element spaces.

Base classes

template<class Operator>
class BilinearFormBase<Operator>
Abstract base class for objects of type BilinearForm.

Public types

using OperatorType = Operator
Type of operator associated to the bilinear form.
using Parent = BilinearFormBase<OperatorType>
Parent class.

Constructors, destructors, conversion operators

BilinearForm(const TrialFunction<TrialFES>& u, const TestFunction<TestFES>& v) constexpr
Constructs a LinearForm with a reference to a TestFunction and a default constructed vector owned by the LinearForm instance.
BilinearForm(const TrialFunction<TrialFES>& u, const TestFunction<TestFES>& v, OperatorType& op) constexpr
Constructs a LinearForm with a reference to a TestFunction and an non-owned vector.
BilinearForm(const TrialFunction<TrialFES>& u, const TestFunction<TestFES>& v, Operator&& op) constexpr
Constructs a LinearForm with a references to a TrialFunction and a TestFunction, and an owned operator.

Public functions

auto operator()(const GridFunction<TrialFES>& u, const GridFunction<TestFES>& v) const -> ScalarType constexpr
Evaluates the linear form at the functions $ u $ and $ v $ .
void assemble() override
Assembles the bilinear form.
auto getTrialFunction() const -> const TrialFunction<TrialFES>& override
Gets the reference to the associated TrialFunction object.
auto getTestFunction() const -> const TestFunction<TestFES>& override
Gets the reference to the associated TestFunction object.
auto operator=(const FormLanguage::List<LocalBilinearFormIntegratorBaseType>& bfis) -> BilinearForm&
auto operator=(const FormLanguage::List<GlobalBilinearFormIntegratorBaseType>& bfis) -> BilinearForm&
auto from(const LocalBilinearFormIntegratorBaseType& bfi) -> BilinearForm&
Builds the bilinear form the given bilinear integrator.
auto add(const LocalBilinearFormIntegratorBaseType& bfi) -> BilinearForm&
Adds a bilinear integrator to the bilinear form.
auto from(const GlobalBilinearFormIntegratorBaseType& bfi) -> BilinearForm&
Builds the bilinear form the given bilinear integrator.
auto add(const GlobalBilinearFormIntegratorBaseType& bfi) -> BilinearForm&
Adds a bilinear integrator to the bilinear form.
auto copy() const -> BilinearForm* override noexcept
Copies the object and returns a non-owning pointer to the copied object.

Function documentation

template<class TrialFES, class TestFES, class Operator>
Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::BilinearForm(const TrialFunction<TrialFES>& u, const TestFunction<TestFES>& v) constexpr

Constructs a LinearForm with a reference to a TestFunction and a default constructed vector owned by the LinearForm instance.

Parameters
u
in Reference to a TestFunction

template<class TrialFES, class TestFES, class Operator>
Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::BilinearForm(const TrialFunction<TrialFES>& u, const TestFunction<TestFES>& v, OperatorType& op) constexpr

Constructs a LinearForm with a reference to a TestFunction and an non-owned vector.

Parameters
u
in Reference to a TestFunction
op

template<class TrialFES, class TestFES, class Operator>
ScalarType Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::operator()(const GridFunction<TrialFES>& u, const GridFunction<TestFES>& v) const constexpr

Evaluates the linear form at the functions $ u $ and $ v $ .

Returns The action $ a(u, v) $ which the bilinear form takes at $ ( u, v ) $ .

Given grid functions $ u $ and $ v $ , this function will compute the action of the bilinear mapping $ a(u, v) $ .

template<class TrialFES, class TestFES, class Operator>
void Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::assemble() override

Assembles the bilinear form.

This method will assemble the underlying sparse matrix associated the bilinear form.

template<class TrialFES, class TestFES, class Operator>
const TrialFunction<TrialFES>& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::getTrialFunction() const override

Gets the reference to the associated TrialFunction object.

Returns Reference to this (for method chaining)

template<class TrialFES, class TestFES, class Operator>
const TestFunction<TestFES>& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::getTestFunction() const override

Gets the reference to the associated TestFunction object.

Returns Reference to this (for method chaining)

template<class TrialFES, class TestFES, class Operator>
BilinearForm& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::operator=(const FormLanguage::List<LocalBilinearFormIntegratorBaseType>& bfis)

template<class TrialFES, class TestFES, class Operator>
BilinearForm& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::operator=(const FormLanguage::List<GlobalBilinearFormIntegratorBaseType>& bfis)

template<class TrialFES, class TestFES, class Operator>
BilinearForm& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::from(const LocalBilinearFormIntegratorBaseType& bfi)

Builds the bilinear form the given bilinear integrator.

Parameters
bfi in Bilinear integrator which will be used to build the bilinear form.
Returns Reference to this (for method chaining)

template<class TrialFES, class TestFES, class Operator>
BilinearForm& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::add(const LocalBilinearFormIntegratorBaseType& bfi)

Adds a bilinear integrator to the bilinear form.

Returns Reference to this (for method chaining)

template<class TrialFES, class TestFES, class Operator>
BilinearForm& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::from(const GlobalBilinearFormIntegratorBaseType& bfi)

Builds the bilinear form the given bilinear integrator.

Parameters
bfi in Bilinear integrator which will be used to build the bilinear form.
Returns Reference to this (for method chaining)

template<class TrialFES, class TestFES, class Operator>
BilinearForm& Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::add(const GlobalBilinearFormIntegratorBaseType& bfi)

Adds a bilinear integrator to the bilinear form.

Returns Reference to this (for method chaining)

template<class TrialFES, class TestFES, class Operator>
BilinearForm* Rodin::Variational::BilinearForm<TrialFES, TestFES, Operator>::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.