LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Public Types | Public Member Functions | Private Attributes | List of all members
projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF > Class Template Reference

Class for local quadrature based computations of sub element matrices corresponding to reaction element matrices. More...

#include <projects/dpg/loc_comp_dpg.h>

Inheritance diagram for projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >:
projects::dpg::SubElementMatrixProvider< SCALAR >

Public Types

using elem_mat_t = typename SubElementMatrixProvider< SCALAR >::elem_mat_t
 inherited types for element matrices More...
 
using ElemMat = typename SubElementMatrixProvider< SCALAR >::ElemMat
 
- Public Types inherited from projects::dpg::SubElementMatrixProvider< SCALAR >
using size_type = lf::uscalfe::size_type
 
using elem_mat_t = Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic >
 internal type for element matrices More...
 
using ElemMat = elem_mat_t
 

Public Member Functions

 ReactionElementMatrixProvider (const ReactionElementMatrixProvider &)=delete
 
 ReactionElementMatrixProvider (ReactionElementMatrixProvider &&) noexcept=default
 
ReactionElementMatrixProvideroperator= (const ReactionElementMatrixProvider &)=delete
 
ReactionElementMatrixProvideroperator= (ReactionElementMatrixProvider &&)=delete
 
 ReactionElementMatrixProvider (std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_trial, std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_test, size_type trial_component, size_type test_component, REACTION_COEFF gamma)
 Constructor: performs cell-independent precomputations. More...
 
ElemMat Eval (const lf::mesh::Entity &cell) override
 main routine for the computation of element matrices More...
 
size_type TrialComponent () const override
 returns the index of the trial space component \( u_{i_k} \) which is the trial space for the bilinear form \( b_k \) More...
 
size_type TestComponent () const override
 returns the index of the test space component \( v_{j_k} \) which is the test space for the bilinear form \( b_k \) More...
 
 ~ReactionElementMatrixProvider () override=default
 
- Public Member Functions inherited from projects::dpg::SubElementMatrixProvider< SCALAR >
 SubElementMatrixProvider ()=default
 
virtual ~SubElementMatrixProvider ()=default
 
 SubElementMatrixProvider (const SubElementMatrixProvider &)=delete
 standard constructors More...
 
 SubElementMatrixProvider (SubElementMatrixProvider &&) noexcept=default
 
SubElementMatrixProvideroperator= (const SubElementMatrixProvider &)=delete
 
SubElementMatrixProvideroperator= (SubElementMatrixProvider &&)=delete
 
virtual bool isActive (const lf::mesh::Entity &)
 All cells are considered active in the default implementation. More...
 
virtual ElemMat Eval (const lf::mesh::Entity &cell)=0
 main routine for the computation of (sub) element matrices More...
 
virtual size_type TrialComponent () const =0
 returns the index of the trial space component \( u_{i_k} \) which is the trial space for the bilinear form \( b_k \) More...
 
virtual size_type TestComponent () const =0
 returns the index of the test space component \( v_{j_k} \) which is the test space for the bilinear form \( b_k \) More...
 

Private Attributes

REACTION_COEFF gamma_
 functor providing reaction coefficient More...
 
std::array< lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, 5 > fe_precomp_trial_
 array containing the precomputed information for the trial space: fe_precomp_trial[i] contains the precomputed reference finite element of the trial space component u for ref_el i. More...
 
std::array< lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, 5 > fe_precomp_test_
 array containing the precomputed information for the test space: fe_precomp_trial[i] contains the precomputed reference finite element of the test space component v for ref_el i. More...
 
size_type trial_component_
 index of trial component u More...
 
size_type test_component_
 index of test component v More...
 

Detailed Description

template<typename SCALAR, typename REACTION_COEFF>
class projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >

Class for local quadrature based computations of sub element matrices corresponding to reaction element matrices.

Template Parameters
SCALARtype for the entries of the element matrices. Must be a field type such as double or std::complex<double>
REACTION_COEFFa MeshFunction that defines the reaction coefficient \( \gamma \). It should be scalar 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 corresponds to the (local) bilinear form

\[ (u,v) \mapsto\int\limits_{K} \gamma(\mathbf{x})u\,v\,\mathrm{d}\mathbf{x} \;, \]

with reaction coefficient \(\gamma \). \(u\) is a component of the trial space and \(v\) a component of the test space.

Template parameter requirement

Definition at line 303 of file loc_comp_dpg.h.

Member Typedef Documentation

◆ elem_mat_t

template<typename SCALAR , typename REACTION_COEFF >
using projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::elem_mat_t = typename SubElementMatrixProvider<SCALAR>::elem_mat_t

inherited types for element matrices

Definition at line 308 of file loc_comp_dpg.h.

◆ ElemMat

template<typename SCALAR , typename REACTION_COEFF >
using projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ElemMat = typename SubElementMatrixProvider<SCALAR>::ElemMat

Definition at line 309 of file loc_comp_dpg.h.

Constructor & Destructor Documentation

◆ ReactionElementMatrixProvider() [1/3]

template<typename SCALAR , typename REACTION_COEFF >
projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ReactionElementMatrixProvider ( const ReactionElementMatrixProvider< SCALAR, REACTION_COEFF > &  )
delete

◆ ReactionElementMatrixProvider() [2/3]

template<typename SCALAR , typename REACTION_COEFF >
projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ReactionElementMatrixProvider ( ReactionElementMatrixProvider< SCALAR, REACTION_COEFF > &&  )
defaultnoexcept

◆ ReactionElementMatrixProvider() [3/3]

template<typename SCALAR , typename REACTION_COEFF >
projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ReactionElementMatrixProvider ( std::shared_ptr< ProductUniformFESpace< SCALAR > >  fe_space_trial,
std::shared_ptr< ProductUniformFESpace< SCALAR > >  fe_space_test,
size_type  trial_component,
size_type  test_component,
REACTION_COEFF  gamma 
)

Constructor: performs cell-independent precomputations.

Parameters
fe_space_trialcollection of specifications for the trial fe space
fe_space_testcollection of specifications for the test fe space
trial_componentIndex of the trial space component \(u\)
test_componentIndex of the test space component \(v\)
gammamesh function for the (scalar valued) reaction coefficient.

This constructur uses local quadrature rules with the degree of exactness chosen as the sum of the polynomial degrees of the two components \( u\) and \( v \)

Definition at line 392 of file loc_comp_dpg.h.

References projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::fe_precomp_test_, projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::fe_precomp_trial_, lf::base::RefEl::kQuad(), lf::base::RefEl::kTria(), lf::quad::make_QuadRule(), projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::test_component_, and projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::trial_component_.

◆ ~ReactionElementMatrixProvider()

template<typename SCALAR , typename REACTION_COEFF >
projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::~ReactionElementMatrixProvider ( )
overridedefault

Member Function Documentation

◆ Eval()

template<typename SCALAR , typename REACTION_COEFF >
ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ElemMat projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval ( const lf::mesh::Entity cell)
overridevirtual

main routine for the computation of element matrices

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

Actual computation is based on numerical quadrature and mapping techniques.

Throws an assertion, in case any specification is missing for the type of cell.

Implements projects::dpg::SubElementMatrixProvider< SCALAR >.

Definition at line 430 of file loc_comp_dpg.h.

References lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::base::RefEl::Id(), lf::geometry::Geometry::IntegrationElement(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::isInitialized(), lf::quad::QuadRule::NumPoints(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::quad::QuadRule::Points(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Qr(), lf::mesh::Entity::RefEl(), and lf::quad::QuadRule::Weights().

◆ operator=() [1/2]

template<typename SCALAR , typename REACTION_COEFF >
ReactionElementMatrixProvider & projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::operator= ( const ReactionElementMatrixProvider< SCALAR, REACTION_COEFF > &  )
delete

◆ operator=() [2/2]

template<typename SCALAR , typename REACTION_COEFF >
ReactionElementMatrixProvider & projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::operator= ( ReactionElementMatrixProvider< SCALAR, REACTION_COEFF > &&  )
delete

◆ TestComponent()

template<typename SCALAR , typename REACTION_COEFF >
size_type projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::TestComponent ( ) const
inlineoverridevirtual

returns the index of the test space component \( v_{j_k} \) which is the test space for the bilinear form \( b_k \)

Implements projects::dpg::SubElementMatrixProvider< SCALAR >.

Definition at line 354 of file loc_comp_dpg.h.

References projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::test_component_.

◆ TrialComponent()

template<typename SCALAR , typename REACTION_COEFF >
size_type projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::TrialComponent ( ) const
inlineoverridevirtual

returns the index of the trial space component \( u_{i_k} \) which is the trial space for the bilinear form \( b_k \)

Implements projects::dpg::SubElementMatrixProvider< SCALAR >.

Definition at line 350 of file loc_comp_dpg.h.

References projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::trial_component_.

Member Data Documentation

◆ fe_precomp_test_

template<typename SCALAR , typename REACTION_COEFF >
std::array<lf::uscalfe::PrecomputedScalarReferenceFiniteElement<SCALAR>, 5> projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::fe_precomp_test_
private

array containing the precomputed information for the test space: fe_precomp_trial[i] contains the precomputed reference finite element of the test space component v for ref_el i.

Definition at line 373 of file loc_comp_dpg.h.

Referenced by projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ReactionElementMatrixProvider().

◆ fe_precomp_trial_

template<typename SCALAR , typename REACTION_COEFF >
std::array<lf::uscalfe::PrecomputedScalarReferenceFiniteElement<SCALAR>, 5> projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::fe_precomp_trial_
private

array containing the precomputed information for the trial space: fe_precomp_trial[i] contains the precomputed reference finite element of the trial space component u for ref_el i.

Definition at line 368 of file loc_comp_dpg.h.

Referenced by projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::ReactionElementMatrixProvider().

◆ gamma_

template<typename SCALAR , typename REACTION_COEFF >
REACTION_COEFF projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::gamma_
private

functor providing reaction coefficient

Definition at line 362 of file loc_comp_dpg.h.

◆ test_component_

template<typename SCALAR , typename REACTION_COEFF >
size_type projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::test_component_
private

◆ trial_component_

template<typename SCALAR , typename REACTION_COEFF >
size_type projects::dpg::ReactionElementMatrixProvider< SCALAR, REACTION_COEFF >::trial_component_
private

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