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 = typename decltype(mesh::utils::MeshFunctionReturnType< DIFF_COEFF >() *Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >())::Scalar |
type of returned element matrix More... | |
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: cell-independent precomputations. More... | |
bool | isActive (const lf::mesh::Entity &) const |
All cells are considered active. More... | |
ElemMat | Eval (const lf::mesh::Entity &cell) const |
main routine for the computation of element matrices More... | |
~DiffusionElementMatrixProvider ()=default | |
DiffusionElementMatrixProvider (const DiffusionElementMatrixProvider &)=delete | |
standard constructors More... | |
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
Definition at line 59 of file loc_comp_ellbvp.h.
using lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> |
Definition at line 69 of file loc_comp_ellbvp.h.
using lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Scalar = typename decltype(mesh::utils::MeshFunctionReturnType<DIFF_COEFF>() * Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>())::Scalar |
type of returned element matrix
Definition at line 66 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: cell-independent precomputations.
fe_space | collection of specifications for scalar-valued parametric reference elements |
alpha | mesh function for the (possibly matrix-valued) diffusion coefficient |
This constructor uses local quadature rules with double the degree of exactness as the polynomial degree of the finite element space.
Definition at line 149 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 |
main routine for the computation of element matrices
cell | reference to the (triangular or quadrilateral) cell for which the element matrix should be computed. |
Actual computation of the element matrix based on numerical quadrature and mapping techniques. The order of the quadrature rule is tied to the polynomial degree of the underlying finite element spaces: for polynomial degree p a quadrature rule is chosen that is exact for polynomials of degree 2p.
Throws an assertion in case the finite element specification is missing for the type of the cell.
Definition at line 159 of file loc_comp_ellbvp.h.
References lf::fe::diffusion_element_matrix_provider_logger, 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 102 of file loc_comp_ellbvp.h.
|
delete |
|
delete |
|
private |
Diffusion coefficient
Definition at line 128 of file loc_comp_ellbvp.h.
|
private |
FE Space
Definition at line 132 of file loc_comp_ellbvp.h.
|
private |
Definition at line 134 of file loc_comp_ellbvp.h.