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

Class for local quadrature based computations of sub element matrices corresponding to element matrices associated to a trace variable in the trial space. More...

#include <projects/dpg/loc_comp_dpg.h>

Inheritance diagram for projects::dpg::TraceElementMatrixProvider< SCALAR, 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

 TraceElementMatrixProvider (const TraceElementMatrixProvider &)=delete
 
 TraceElementMatrixProvider (TraceElementMatrixProvider &&) noexcept=default
 
TraceElementMatrixProvideroperator= (const TraceElementMatrixProvider &)=delete
 
TraceElementMatrixProvideroperator= (TraceElementMatrixProvider &&)=delete
 
 TraceElementMatrixProvider (std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_trial, std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_test, size_type trial_component, size_type test_component, COEFF beta)
 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...
 
 ~TraceElementMatrixProvider () 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

COEFF beta_
 functor providing the coefficient invoved in the bilinear form More...
 
std::array< std::vector< lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR > >, 5 > fe_precomp_trial_
 array containing the precomputed information for the trial space: fe_precomp_trial[i][j] contains the precomputed reference finite element of \( \hat{u} \) for ref_el i evaluated on the quadrature points transformed to edge j of the boundary. More...
 
std::array< std::vector< lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR > >, 5 > fe_precomp_test_
 array containing the precomputed information for the test space: fe_precomp_trial[i][j] contains the precomputed reference finite element of \( v \) for ref_el i evaluated on the quadrature points transformed to edge j of the boundary. More...
 
std::array< lf::quad::QuadRule, 5 > segment_qr_
 the original segment quad rule used to construct the quadrules on ref_el boundaries. segment_qr[i] contains the segment qr for ref_el i. More...
 
size_type trial_component_
 index of trial component \( \hat{u} \) More...
 
size_type test_component_
 index of test component v More...
 

Detailed Description

template<typename SCALAR, typename COEFF>
class projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >

Class for local quadrature based computations of sub element matrices corresponding to element matrices associated to a trace variable in the trial space.

Template Parameters
SCALARtype for the entries of the element matrices. Must be a field type such as double or std::complex<double>
COEFFa MeshFunction that defines a coefficient \( \mathbf{\beta} \). It should be vector 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

\[ (\hat u,v) \mapsto\int\limits_{\partial K} \hat u \mathbf{n}_K \cdot \boldsymbol{\beta}(\mathbf{x}) v\mathrm{d}\mathbf{x} \;, \]

\( \hat{u} \) is a trace component of the trial space and \( v \) a component of the test space.

Template parameter requirement

Definition at line 1053 of file loc_comp_dpg.h.

Member Typedef Documentation

◆ elem_mat_t

template<typename SCALAR , typename COEFF >
using projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::elem_mat_t = typename SubElementMatrixProvider<SCALAR>::elem_mat_t

inherited types for element matrices

Definition at line 1058 of file loc_comp_dpg.h.

◆ ElemMat

template<typename SCALAR , typename COEFF >
using projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::ElemMat = typename SubElementMatrixProvider<SCALAR>::ElemMat

Definition at line 1059 of file loc_comp_dpg.h.

Constructor & Destructor Documentation

◆ TraceElementMatrixProvider() [1/3]

template<typename SCALAR , typename COEFF >
projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::TraceElementMatrixProvider ( const TraceElementMatrixProvider< SCALAR, COEFF > &  )
delete

◆ TraceElementMatrixProvider() [2/3]

template<typename SCALAR , typename COEFF >
projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::TraceElementMatrixProvider ( TraceElementMatrixProvider< SCALAR, COEFF > &&  )
defaultnoexcept

◆ TraceElementMatrixProvider() [3/3]

template<typename SCALAR , typename COEFF >
projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::TraceElementMatrixProvider ( std::shared_ptr< ProductUniformFESpace< SCALAR > >  fe_space_trial,
std::shared_ptr< ProductUniformFESpace< SCALAR > >  fe_space_test,
size_type  trial_component,
size_type  test_component,
COEFF  beta 
)

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 \(\hat u\)
test_componentIndex of the test space component \(v\)
betamesh function for the vector valued coefficient involved in the bilinear form.

This constructur uses local quadrature rules. On each edge of the boundary a quadrature rule with degree of exactness chosen as the sum of the polynomial degree of the two components \( \hat{u} \) and \( v \) is used.

Definition at line 1150 of file loc_comp_dpg.h.

References projects::dpg::BoundaryQuadRule(), projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::fe_precomp_test_, projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::fe_precomp_trial_, lf::base::RefEl::kQuad(), lf::base::RefEl::kSegment(), lf::base::RefEl::kTria(), lf::quad::make_QuadRule(), projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::segment_qr_, projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::test_component_, and projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::trial_component_.

◆ ~TraceElementMatrixProvider()

template<typename SCALAR , typename COEFF >
projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::~TraceElementMatrixProvider ( )
overridedefault

Member Function Documentation

◆ Eval()

template<typename SCALAR , typename COEFF >
TraceElementMatrixProvider< SCALAR, COEFF >::ElemMat projects::dpg::TraceElementMatrixProvider< SCALAR, 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 1207 of file loc_comp_dpg.h.

References lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::base::RefEl::Id(), lf::base::RefEl::NumSubEntities(), projects::dpg::OuterNormals(), lf::mesh::Entity::RefEl(), and lf::geometry::Geometry::SubGeometry().

◆ operator=() [1/2]

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

◆ operator=() [2/2]

template<typename SCALAR , typename COEFF >
TraceElementMatrixProvider & projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::operator= ( TraceElementMatrixProvider< SCALAR, COEFF > &&  )
delete

◆ TestComponent()

template<typename SCALAR , typename COEFF >
size_type projects::dpg::TraceElementMatrixProvider< SCALAR, 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 1102 of file loc_comp_dpg.h.

References projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::test_component_.

◆ TrialComponent()

template<typename SCALAR , typename COEFF >
size_type projects::dpg::TraceElementMatrixProvider< SCALAR, 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 1098 of file loc_comp_dpg.h.

References projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::trial_component_.

Member Data Documentation

◆ beta_

template<typename SCALAR , typename COEFF >
COEFF projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::beta_
private

functor providing the coefficient invoved in the bilinear form

Definition at line 1110 of file loc_comp_dpg.h.

◆ fe_precomp_test_

template<typename SCALAR , typename COEFF >
std::array< std::vector<lf::uscalfe::PrecomputedScalarReferenceFiniteElement<SCALAR> >, 5> projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::fe_precomp_test_
private

array containing the precomputed information for the test space: fe_precomp_trial[i][j] contains the precomputed reference finite element of \( v \) for ref_el i evaluated on the quadrature points transformed to edge j of the boundary.

Definition at line 1129 of file loc_comp_dpg.h.

Referenced by projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::TraceElementMatrixProvider().

◆ fe_precomp_trial_

template<typename SCALAR , typename COEFF >
std::array< std::vector<lf::uscalfe::PrecomputedScalarReferenceFiniteElement<SCALAR> >, 5> projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::fe_precomp_trial_
private

array containing the precomputed information for the trial space: fe_precomp_trial[i][j] contains the precomputed reference finite element of \( \hat{u} \) for ref_el i evaluated on the quadrature points transformed to edge j of the boundary.

Definition at line 1120 of file loc_comp_dpg.h.

Referenced by projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::TraceElementMatrixProvider().

◆ segment_qr_

template<typename SCALAR , typename COEFF >
std::array<lf::quad::QuadRule, 5> projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::segment_qr_
private

the original segment quad rule used to construct the quadrules on ref_el boundaries. segment_qr[i] contains the segment qr for ref_el i.

Definition at line 1133 of file loc_comp_dpg.h.

Referenced by projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::TraceElementMatrixProvider().

◆ test_component_

template<typename SCALAR , typename COEFF >
size_type projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::test_component_
private

◆ trial_component_

template<typename SCALAR , typename COEFF >
size_type projects::dpg::TraceElementMatrixProvider< SCALAR, COEFF >::trial_component_
private

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