LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Class for local quadrature based computation of element matrix for Lagrangian finite elements and a weighted \(L^2\) inner product. More...
#include <lf/fe/fe.h>
Public Types | |
using | Scalar = typename decltype(mesh::utils::MeshFunctionReturnType< REACTION_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 | |
MassElementMatrixProvider (std::shared_ptr< const ScalarFESpace< SCALAR > > fe_space, REACTION_COEFF gamma) | |
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... | |
~MassElementMatrixProvider ()=default | |
MassElementMatrixProvider (const MassElementMatrixProvider &)=delete | |
standard constructors More... | |
MassElementMatrixProvider (MassElementMatrixProvider &&) noexcept=default | |
MassElementMatrixProvider & | operator= (const MassElementMatrixProvider &)=delete |
MassElementMatrixProvider & | operator= (MassElementMatrixProvider &&)=delete |
Private Attributes | |
std::shared_ptr< const ScalarFESpace< SCALAR > > | fe_space_ |
quad::QuadRuleCache | qr_cache_ |
functors providing coefficient functions | |
REACTION_COEFF | gamma_ |
Class for local quadrature based computation of element matrix for Lagrangian finite elements and a weighted \(L^2\) inner product.
SCALAR | scalar type of the FiniteElementSpace. Must be a field type such as double or std::complex<double> |
REACTION_COEFF | a MeshFunction that defines the reaction coefficient \( \mathbf{\gamma} \). It should be either scalar- or matrix-valued. |
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 Example 2.8.3.29
double
Definition at line 238 of file loc_comp_ellbvp.h.
using lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> |
Definition at line 248 of file loc_comp_ellbvp.h.
using lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Scalar = typename decltype(mesh::utils::MeshFunctionReturnType<REACTION_COEFF>() * Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>())::Scalar |
type of returned element matrix
Definition at line 245 of file loc_comp_ellbvp.h.
|
delete |
standard constructors
|
defaultnoexcept |
lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::MassElementMatrixProvider | ( | std::shared_ptr< const ScalarFESpace< SCALAR > > | fe_space, |
REACTION_COEFF | gamma | ||
) |
Constructor: cell-independent precomputations.
fe_space | collection of specifications for scalar-valued parametric reference elements |
gamma | mesh function providing scalar-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 326 of file loc_comp_ellbvp.h.
|
default |
destructor
lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::ElemMat lf::fe::MassElementMatrixProvider< SCALAR, REACTION_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 o degree 2p.
Throws an assertion in case the finite element specification is missing for the type of the cell.
Definition at line 335 of file loc_comp_ellbvp.h.
References lf::geometry::Geometry::DimGlobal(), lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::geometry::Geometry::Global(), lf::geometry::Geometry::IntegrationElement(), lf::fe::mass_element_matrix_provider_logger, 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 279 of file loc_comp_ellbvp.h.
|
delete |
|
delete |
|
private |
FE Space
Definition at line 309 of file loc_comp_ellbvp.h.
|
private |
Reaction coefficient
Definition at line 305 of file loc_comp_ellbvp.h.
|
private |
Definition at line 311 of file loc_comp_ellbvp.h.