LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF > Class Template Reference

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
 
ReactionDiffusionElementMatrixProvideroperator= (const ReactionDiffusionElementMatrixProvider &)=delete
 
ReactionDiffusionElementMatrixProvideroperator= (ReactionDiffusionElementMatrixProvider &&)=delete
 

Private Attributes

std::array< PrecomputedScalarReferenceFiniteElement< SCALAR >, 5 > fe_precomp_
 
functors providing coefficient functions
DIFF_COEFF alpha_
 
REACTION_COEFF gamma_
 

Detailed Description

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
class lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >

Class for local quadrature based computations for Lagrangian finite elements and second-order scalar elliptic BVPs.

Template Parameters
SCALARscalar type of the UniformScalarFESpace. Must be a field type such as double or std::complex<double>
DIFF_COEFFa MeshFunction that defines the diffusion coefficient \( \mathbf{\alpha} \). It should be either scalar- or matrix-valued.
REACTION_COEFFa MeshFunction that defines the reaction coefficient \( \mathbf{\gamma} \). It should be either scalar- or matrix-valued.
Note
This class complies with the type requirements for the template argument ENTITY_MATRIX_PROVIDER of the function lf::assemble::AssembleMatrixLocally().

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

Template parameter requirement

Note
The constructors of this class want an object of type UniformScalarFESpace, which holds a pointer to a mesh. However, for local builder classes global information about the mesh is irrelevant, and, therefore this object is used only to obtain information about the local shape functions. A revised implementation should directly pass this information to the constructor.

Logger

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.

Member Typedef Documentation

◆ ElemMat

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
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.

◆ Scalar

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
using lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Scalar
Initial value:
typename decltype(mesh::utils::MeshFunctionReturnType<DIFF_COEFF>() *
Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>() +
mesh::utils::MeshFunctionReturnType<REACTION_COEFF>() *
Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>())::Scalar
typename decltype(mesh::utils::MeshFunctionReturnType< DIFF_COEFF >() * Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >()+ mesh::utils::MeshFunctionReturnType< REACTION_COEFF >() * Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 >())::Scalar Scalar
type of returned element matrix

type of returned element matrix

Definition at line 92 of file loc_comp_ellbvp.h.

Constructor & Destructor Documentation

◆ ReactionDiffusionElementMatrixProvider() [1/4]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider ( const ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF > & )
delete

standard constructors

◆ ReactionDiffusionElementMatrixProvider() [2/4]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider ( ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF > && )
defaultnoexcept

◆ ReactionDiffusionElementMatrixProvider() [3/4]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
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.

Parameters
fe_spacecollection of specifications for scalar-valued parametric reference elements
alphamesh function for the (possibly matrix-valued) diffusion coefficient
gammamesh function providing scalar-valued diffusion coefficient
See also
LocCompLagrFEPreprocessor::LocCompLagrFEPreprocessor()

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().

◆ ReactionDiffusionElementMatrixProvider() [4/4]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
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.

Parameters
fe_spacecollection of specifications for scalar-valued parametric reference elements
alphamesh function for the (possibly matrix-valued) diffusion coefficient
gammamesh function providing scalar-valued diffusion coefficient
qr_collectioncollection 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.
See also
LocCompLagrFEPreprocessor::LocCompLagrFEPreprocessor()

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().

◆ ~ReactionDiffusionElementMatrixProvider()

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
virtual lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::~ReactionDiffusionElementMatrixProvider ( )
virtualdefault

Virtual destructor

Member Function Documentation

◆ Eval()

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
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

Parameters
cellreference to the (triangular or quadrilateral) cell for which the element matrix should be computed.
Returns
a small dense, containing the element matrix.

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.

Exceptions
base::LfExceptionin 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().

◆ isActive()

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
virtual bool lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::isActive ( const lf::mesh::Entity & )
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.

◆ operator=() [1/2]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
ReactionDiffusionElementMatrixProvider & lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::operator= ( const ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF > & )
delete

◆ operator=() [2/2]

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
ReactionDiffusionElementMatrixProvider & lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::operator= ( ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF > && )
delete

Member Data Documentation

◆ alpha_

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
DIFF_COEFF lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::alpha_
private

Diffusion coefficient

Definition at line 182 of file loc_comp_ellbvp.h.

◆ fe_precomp_

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
std::array<PrecomputedScalarReferenceFiniteElement<SCALAR>, 5> lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::fe_precomp_
private

◆ gamma_

template<base::Scalar SCALAR, mesh::utils::MeshFunction DIFF_COEFF, mesh::utils::MeshFunction REACTION_COEFF>
REACTION_COEFF lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::gamma_
private

Reaction coefficient

Definition at line 184 of file loc_comp_ellbvp.h.


The documentation for this class was generated from the following file: