LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Class for local quadrature based computations of sub element matrices corresponding to diffusion element matrices. More...
#include <projects/dpg/loc_comp_dpg.h>
Public Types | |
using | elem_mat_t = typename SubElementMatrixProvider< SCALAR >::elem_mat_t |
inherited types for element matrices More... | |
using | ElemMat = typename SubElementMatrixProvider< SCALAR >::ElemMat |
![]() | |
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 | |
DiffusionElementMatrixProvider (const DiffusionElementMatrixProvider &)=delete | |
DiffusionElementMatrixProvider (DiffusionElementMatrixProvider &&) noexcept=default | |
DiffusionElementMatrixProvider & | operator= (const DiffusionElementMatrixProvider &)=delete |
DiffusionElementMatrixProvider & | operator= (DiffusionElementMatrixProvider &&)=delete |
DiffusionElementMatrixProvider (std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_trial, std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_test, size_type trial_component, size_type test_component, DIFF_COEFF alpha) | |
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... | |
~DiffusionElementMatrixProvider () override=default | |
![]() | |
SubElementMatrixProvider ()=default | |
virtual | ~SubElementMatrixProvider ()=default |
SubElementMatrixProvider (const SubElementMatrixProvider &)=delete | |
standard constructors More... | |
SubElementMatrixProvider (SubElementMatrixProvider &&) noexcept=default | |
SubElementMatrixProvider & | operator= (const SubElementMatrixProvider &)=delete |
SubElementMatrixProvider & | operator= (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 | |
DIFF_COEFF | alpha_ |
functor providing the diffusion 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... | |
Class for local quadrature based computations of sub element matrices corresponding to diffusion element matrices.
SCALAR | type for the entries of the element matrices. 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. |
The element matrix corresponds to the (local) bilinear form
\[ (u,v) \mapsto\int\limits_{K}\mathbf{grad}\,u \cdot \boldsymbol{\alpha}(\mathbf{x})\mathbf{grad}\,v \mathrm{d}\mathbf{x} \;, \]
with diffusion coefficient \(\mathbf{\alpha}\). \(u\) is a component of the trial space and \(v\) a component of the test space.
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 or fixed size. Definition at line 70 of file loc_comp_dpg.h.
using projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::elem_mat_t = typename SubElementMatrixProvider<SCALAR>::elem_mat_t |
inherited types for element matrices
Definition at line 75 of file loc_comp_dpg.h.
using projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::ElemMat = typename SubElementMatrixProvider<SCALAR>::ElemMat |
Definition at line 76 of file loc_comp_dpg.h.
|
delete |
|
defaultnoexcept |
projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider | ( | std::shared_ptr< ProductUniformFESpace< SCALAR > > | fe_space_trial, |
std::shared_ptr< ProductUniformFESpace< SCALAR > > | fe_space_test, | ||
size_type | trial_component, | ||
size_type | test_component, | ||
DIFF_COEFF | alpha | ||
) |
Constructor: performs cell-independent precomputations.
fe_space_trial | collection of specifications for the trial fe space |
fe_space_test | collection of specifications for the test fe space |
trial_component | Index of the trial space component \(u\) |
test_component | Index of the test space component \(v\) |
alpha | mesh function for the (matrix or scalar valued) diffusion 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 159 of file loc_comp_dpg.h.
References projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::fe_precomp_test_, projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::fe_precomp_trial_, lf::base::RefEl::kQuad(), lf::base::RefEl::kTria(), lf::quad::make_QuadRule(), projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::test_component_, and projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::trial_component_.
|
overridedefault |
|
overridevirtual |
main routine for the computation of element matrices
cell | reference to a (triangular or quadrilateral) cell for which the element matirx should be computed. |
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 196 of file loc_comp_dpg.h.
References lf::geometry::Geometry::DimGlobal(), lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::base::RefEl::Id(), 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 >::Qr(), lf::mesh::Entity::RefEl(), and lf::quad::QuadRule::Weights().
|
delete |
|
delete |
|
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 122 of file loc_comp_dpg.h.
References projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::test_component_.
|
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 118 of file loc_comp_dpg.h.
References projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::trial_component_.
|
private |
functor providing the diffusion coefficient
Definition at line 130 of file loc_comp_dpg.h.
|
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 141 of file loc_comp_dpg.h.
Referenced by projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider().
|
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 136 of file loc_comp_dpg.h.
Referenced by projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider().
|
private |
index of test component v
Definition at line 146 of file loc_comp_dpg.h.
Referenced by projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider(), and projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::TestComponent().
|
private |
index of trial component u
Definition at line 144 of file loc_comp_dpg.h.
Referenced by projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::DiffusionElementMatrixProvider(), and projects::dpg::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::TrialComponent().