LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Class for computing element matrices for general scalar-valued finite elements and homogeneous 2nd-order elliptic bilinear forms. More...
#include <lf/fe/fe.h>
Public Types | |
using | Scalar |
type of returned element matrix | |
using | ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> |
Public Member Functions | |
DiffusionElementMatrixProvider (std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space, DIFF_COEFF alpha) | |
Constructor. | |
bool | isActive (const lf::mesh::Entity &) const |
All cells are considered active. | |
ElemMat | Eval (const lf::mesh::Entity &cell) const |
Routine for the computation of element matrices. | |
~DiffusionElementMatrixProvider ()=default | |
DiffusionElementMatrixProvider (const DiffusionElementMatrixProvider &)=delete | |
standard constructors | |
DiffusionElementMatrixProvider (DiffusionElementMatrixProvider &&) noexcept=default | |
DiffusionElementMatrixProvider & | operator= (const DiffusionElementMatrixProvider &)=delete |
DiffusionElementMatrixProvider & | operator= (DiffusionElementMatrixProvider &&)=delete |
Private Attributes | |
std::shared_ptr< const ScalarFESpace< SCALAR > > | fe_space_ |
quad::QuadRuleCache | qr_cache_ |
functors providing coefficient functions | |
DIFF_COEFF | alpha_ |
Class for computing element matrices for general scalar-valued finite elements and homogeneous 2nd-order elliptic bilinear forms.
SCALAR | scalar type of the ScalarFESpace Must be a field type such as double or std::complex<double> |
DIFF_COEFF | a MeshFunction that defines the diffusion coefficient \( \mathbf{\alpha} \). It should be either scalar- or matrix-valued. |
The element matrix corresponds to the (local) bilinear form
\[ (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 2.8.3.29
double
lf::uscalfe
please consider switching to lf::uscalfe::ReactionDiffusionElementMatrixProvider
provided there, as it is specifically optimized for uniform FE Spaces.This class logs additional information to DiffusionElementMatrixProviderLogger(). See Loggers and Debug output for more information.
The following code snippet computes the solution of the BVP
\begin{align} - \Delta u &= 1 && \text{on }\Omega := [0,1]^2 \\ u &= 0 && \text{on }\partial \Omega \end{align}
Definition at line 76 of file loc_comp_ellbvp.h.
using lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> |
Definition at line 84 of file loc_comp_ellbvp.h.
using lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Scalar |
type of returned element matrix
Definition at line 81 of file loc_comp_ellbvp.h.
|
delete |
standard constructors
|
defaultnoexcept |
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider | ( | std::shared_ptr< const ScalarFESpace< SCALAR > > | fe_space, |
DIFF_COEFF | alpha ) |
Constructor.
fe_space | collection of specifications for scalar-valued parametric reference elements |
alpha | mesh function for the (possibly matrix-valued) diffusion coefficient |
fe_space
is from lf::uscalfe, please consider using the lf::uscalfe::DiffusionElementMatrixProvider provided there, as that is specifically optimized for uniform FE spaces. Definition at line 161 of file loc_comp_ellbvp.h.
|
default |
destructor
lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::ElemMat lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval | ( | const lf::mesh::Entity & | cell | ) | const |
Routine for the computation of element matrices.
cell | reference to the (triangular or quadrilateral) cell for which the element matrix should be computed. |
Evaluation of the integral over the basis functions is done using numerical quadrature with double the degree of exactness as the polynomial degree of the basis functions on the provided cell.
Throws an assertion in case the finite element specification is missing for the type of the cell.
Definition at line 171 of file loc_comp_ellbvp.h.
References lf::fe::DiffusionElementMatrixProviderLogger(), lf::geometry::Geometry::DimGlobal(), lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::geometry::Geometry::Global(), lf::geometry::Geometry::IntegrationElement(), lf::geometry::Geometry::JacobianInverseGramian(), lf::base::RefEl::NodeCoords(), lf::quad::QuadRule::NumPoints(), lf::quad::QuadRule::Points(), lf::mesh::Entity::RefEl(), and lf::quad::QuadRule::Weights().
|
inline |
All cells are considered active.
Definition at line 116 of file loc_comp_ellbvp.h.
|
delete |
|
delete |
|
private |
Diffusion coefficient
Definition at line 140 of file loc_comp_ellbvp.h.
|
private |
FE Space
Definition at line 144 of file loc_comp_ellbvp.h.
|
private |
Definition at line 146 of file loc_comp_ellbvp.h.