![]() |
LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
LehrFEM++ includes built-in lf::assemble::EntityMatrixProvider
(EMPs) and lf::assemble::EntityVectorProvider
(EVPs) for common PDEs. Custom EMPs/EVPs can also be defined. These providers compute local element matrices/vectors for finite element methods and are typically passed to lf::assemble::AssembleMatrixLocally
or lf::assemble::AssembleVectorLocally
(see Assembly).
A custom EMP must implement the lf::assemble::EntityMatrixProvider
concept. Minimal definitions:
Below is a summary of key built-in EMPs:
Equation | Name (follow link for details) | Description |
---|---|---|
\( \int\limits_{K}\boldsymbol{\alpha}(\mathbf{x})\mathbf{grad}\,u \cdot\mathbf{grad}\,v\,\mathrm{d}\mathbf{x} \) | fe::DiffusionElementMatrixProvider | stiffness matrix for the diffusion equation. |
\( \int\limits_{K}\gamma(\mathbf{x})u\,\overline{v}\,\mathrm{d}\mathbf{x} \) | fe::MassElementMatrixProvider | mass matrix for the reaction-diffusion equation. |
\( \int\limits_e\gamma(x)u(x)\overline{v(x)}\,\mathrm{d}S(x) \) | fe::MassEdgeMatrixProvider | mass matrix for the reaction-diffusion equation. |
\( \int\limits_{K}\mathbf{grad}\,u\cdot\mathbf{grad}\,v\,\mathrm{d}\mathbf{x} \) | uscalfe::LinearFELaplaceElementMatrix | stiffness matrix for the Laplace equation. |
\( \int\limits_{K}\boldsymbol{\alpha}(\mathbf{x})\mathbf{grad}\,u\cdot\mathbf{grad}\,v\,\mathrm{d}\mathbf{x} \) \( + \int\limits_{K}\gamma(\mathbf{x})u\,\overline{v}\,\mathrm{d}\mathbf{x} \) | uscalfe::ReactionDiffusionElementMatrixProvider | stiffness matrix for the reaction-diffusion equation. |
\( \int\limits_e\gamma(x)u(x)\overline{v(x)}\,\mathrm{d}S(x) \) | uscale::MassEdgeMatrixProvider | mass matrix for the reaction-diffusion equation. |
lf::fe::DiffusionElementMatrixProvider
\[ (u,v) \mapsto\int\limits_{K}\boldsymbol{\alpha}(\mathbf{x})\mathbf{grad}\,u \cdot\mathbf{grad}\,v\,\mathrm{d}\mathbf{x} \; \]
with diffusion coefficient \(\mathbf{\alpha}\), see also Example Lecture Document Example 2.8.3.29
lf::fe::MassElementMatrixProvider
The element matrix corresponds to the (local) bilinear form
\[ (u,v) \mapsto\int\limits_{K}\gamma(\mathbf{x})u\,\overline{v}\,\mathrm{d}\mathbf{x} \;, \]
with reaction coefficient \(\gamma\), see also Lecture Document Example 2.8.3.29
lf::fe::MassEdgeMatrixProvider
The edge matrix corresponds to the (local) bilinear form
\[ (u,v) \mapsto \int\limits_e\gamma(x)u(x)\overline{v(x)}\,\mathrm{d}S(x)\;, \]
where \(e\) is an edge of the mesh, and \(\gamma\) a scalar-valued coefficient function.
lf::uscalfe::LinearFELaplaceElementMatrix
The element matrix corresponds to the (local) bilinear form
\[ (u,v) \mapsto\int\limits_{K}\mathbf{grad}\,u\cdot\mathbf{grad}\,v\,\mathrm{d}\mathbf{x} \;, \]
lf::uscalfe::ReactionDiffusionElementMatrixProvider
\[ (u,v) \mapsto\int\limits_{K}\boldsymbol{\alpha}(\mathbf{x})\mathbf{grad}\,u\cdot\mathbf{grad}\,v + \gamma(\mathbf{x})u\,\overline{v}\,\mathrm{d}\mathbf{x} \;, \]
with diffusion coefficient \(\mathbf{\alpha}\) and reaction coefficient \(\gamma\). See also Lecture Document Example 2.8.3.29.
lf::uscalfe::MassEdgeMatrixProvider
\[ (u,v) \mapsto \int\limits_e\gamma(x)u(x)\overline{v(x)}\,\mathrm{d}S(x)\;, \]
where \(e\) is an edge of the mesh, and \(\gamma\) a scalar-valued coefficient function. See Lecture Document Example 2.7.4.33 for a similar example on assembly of boundary contributions.
Entity Vector Providers (EVPs) in LehrFEM++ are used to compute local element vectors for finite element methods. Below is a table summarizing the key built-in EVPs in LehrFEM++.
Equation | Name (follow link for details) | Description |
---|---|---|
\( \int\limits_{K}f(\mathbf{x})v\,\mathrm{d}\mathbf{x} \) | fe::ScalarLoadElementVectorProvider | load vector for the scalar load on an element. |
\( \int\limits_{e}g(x)\overline{v(x)}\,\mathrm{d}S(x) \) | fe::ScalarLoadEdgeVectorProvider | load vector for the scalar load on an edge \(e\). |
\( \int\limits_{K}f(\mathbf{x})\overline{v(x)}\,\mathrm{d}\mathbf{x} \) | uscalfe::ScalarLoadElementVectorProvider | load vector for scalar finite elements; volume contributions only. |
\( \int\limits_{e}g(x)\overline{v(x)}\,\mathrm{d}S(x) \) | uscalfe::ScalarLoadEdgeVectorProvider | load vector for scalar finite elements |
\( \int\limits_{K}f(\mathbf{x})\overline{v(x)}\,\mathrm{d}\mathbf{x} \) | uscalfe::LinearFELocalLoadVector | load vector for scalar linear finite elements. |
lf::fe::ScalarLoadElementVectorProvider computes the local element vector corresponding to the (local) linear form
\[ v \mapsto\int\limits_{K}f(\mathbf{x})v\,\mathrm{d}\mathbf{x} \;, \]
where \(f\) is a scalar-valued function.
lf::fe::ScalarLoadEdgeVectorProvider computes the local edge contributions corresponding to the (local) linear form
\[ v \mapsto\int\limits_{e}g(x)\overline{v(x)}\,\mathrm{d}S(x) \;, \]
where \(e\) is an edge of the mesh, and \(g\) a scalar-valued function.
lf::uscalfe::ScalarLoadElementVectorProvider computes the local element vector corresponding to the (local) linear form
\[ v \mapsto\int\limits_{K}f(\mathbf{x})\overline{v(x)}\,\mathrm{d}\mathbf{x} \;, \]
where \(f\) is a scalar-valued function. Note: This provider only computes the volume contributions on \(K\).
lf::uscalfe::ScalarLoadEdgeVectorProvider computes the local edge contributions corresponding to the (local) linear form
\[ v \mapsto\int\limits_{e}g(x)\overline{v(x)}\,\mathrm{d}S(x) \;, \]
where \(e\) is an edge of the mesh, and \(g\) a scalar-valued function.
lf::uscalfe::LinearFELocalLoadVector computes the local element vector for linear FE corresponding to the (local) linear form
\[ v \mapsto\int\limits_{K}f(\mathbf{x})\overline{v(x)}\,\mathrm{d}\mathbf{x} \;, \]
where \(f\) is a scalar-valued function using edge midpoint quadrature.