LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Class for local quadrature based computations for Lagrangian finite elements and second-order scalar elliptic BVPs. More...
#include <lf/uscalfe/uscalfe.h>
Public Types | |
using | Scalar |
type of returned element matrix | |
using | ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> |
Public Member Functions | |
ReactionDiffusionElementMatrixProvider (std::shared_ptr< const UniformScalarFESpace< SCALAR > > fe_space, DIFF_COEFF alpha, REACTION_COEFF gamma) | |
Constructor: cell-independent precomputations. | |
ReactionDiffusionElementMatrixProvider (std::shared_ptr< const UniformScalarFESpace< SCALAR > > fe_space, DIFF_COEFF alpha, REACTION_COEFF gamma, quad_rule_collection_t qr_collection) | |
Constructor: cell-independent precomputations and custom quadrature rule. | |
virtual bool | isActive (const lf::mesh::Entity &) |
All cells are considered active in the default implementation. | |
ElemMat | Eval (const lf::mesh::Entity &cell) |
main routine for the computation of element matrices | |
virtual | ~ReactionDiffusionElementMatrixProvider ()=default |
ReactionDiffusionElementMatrixProvider (const ReactionDiffusionElementMatrixProvider &)=delete | |
standard constructors | |
ReactionDiffusionElementMatrixProvider (ReactionDiffusionElementMatrixProvider &&) noexcept=default | |
ReactionDiffusionElementMatrixProvider & | operator= (const ReactionDiffusionElementMatrixProvider &)=delete |
ReactionDiffusionElementMatrixProvider & | operator= (ReactionDiffusionElementMatrixProvider &&)=delete |
Private Attributes | |
std::array< PrecomputedScalarReferenceFiniteElement< SCALAR >, 5 > | fe_precomp_ |
functors providing coefficient functions | |
DIFF_COEFF | alpha_ |
REACTION_COEFF | gamma_ |
Class for local quadrature based computations for Lagrangian finite elements and second-order scalar elliptic BVPs.
SCALAR | scalar type of the UniformScalarFESpace. 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. |
REACTION_COEFF | a MeshFunction that defines the reaction coefficient \( \mathbf{\gamma} \). It should be either scalar- or matrix-valued. |
The element matrix is corresponds to the (local) bilinear form
\[ (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 Example 2.8.3.29
double
operator (const Entity &,ref_coord_t)
that returns either a scalar or a matrix type that is compatible with Eigen's matrices. Usually it will be an Eigen::Matrix either of variable of fixed size.This class logs additional information to ReactionDiffusionElementMatrixProviderLogger(). See Loggers and Debug output for more information.
Definition at line 87 of file loc_comp_ellbvp.h.
using lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ElemMat = Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> |
Definition at line 97 of file loc_comp_ellbvp.h.
using lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Scalar |
type of returned element matrix
Definition at line 92 of file loc_comp_ellbvp.h.
|
delete |
standard constructors
|
defaultnoexcept |
lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider | ( | std::shared_ptr< const UniformScalarFESpace< SCALAR > > | fe_space, |
DIFF_COEFF | alpha, | ||
REACTION_COEFF | gamma ) |
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 |
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 212 of file loc_comp_ellbvp.h.
References lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::fe_precomp_, lf::base::RefEl::kQuad(), lf::base::RefEl::kTria(), and lf::quad::make_QuadRule().
lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider | ( | std::shared_ptr< const UniformScalarFESpace< SCALAR > > | fe_space, |
DIFF_COEFF | alpha, | ||
REACTION_COEFF | gamma, | ||
quad_rule_collection_t | qr_collection ) |
Constructor: cell-independent precomputations and custom quadrature rule.
fe_space | collection of specifications for scalar-valued parametric reference elements |
alpha | mesh function for the (possibly matrix-valued) diffusion coefficient |
gamma | mesh function providing scalar-valued diffusion coefficient |
qr_collection | collection of quadrature rules. A quadrature rule is required for every cell type for which the finite element space provides local shape functions. If a quadrature rule is not specified for a cell type and the Eval() method is called for such a cell, exception will be thrown. |
Definition at line 236 of file loc_comp_ellbvp.h.
References lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::fe_precomp_, lf::base::RefEl::kQuad(), lf::base::RefEl::kTria(), and lf::quad::QuadRule::RefEl().
|
virtualdefault |
Virtual destructor
lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ElemMat lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval | ( | const lf::mesh::Entity & | cell | ) |
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 Lagrangian finite element spaces: for polynomial degree p a quadrature rule is chosen that is exact for polynomials o degree 2p.
base::LfException | in case the finite element specification is missing for the type of the cell or if there is no quadrature rule specified for the given cell type. |
Definition at line 271 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::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::isInitialized(), lf::geometry::Geometry::JacobianInverseGramian(), lf::quad::QuadRule::NumPoints(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::quad::QuadRule::Points(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompGradientsReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Qr(), lf::uscalfe::ReactionDiffusionElementMatrixProviderLogger(), lf::mesh::Entity::RefEl(), and lf::quad::QuadRule::Weights().
|
inlinevirtual |
All cells are considered active in the default implementation.
This method is meant to be overloaded if assembly should be restricted to a subset of cells.
Definition at line 155 of file loc_comp_ellbvp.h.
|
delete |
|
delete |
|
private |
Diffusion coefficient
Definition at line 182 of file loc_comp_ellbvp.h.
|
private |
Definition at line 188 of file loc_comp_ellbvp.h.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider(), and lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider().
|
private |
Reaction coefficient
Definition at line 184 of file loc_comp_ellbvp.h.