LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Class for local quadrature based computations of sub vectors corresponding to load vectors. More...
#include <projects/dpg/loc_comp_dpg.h>
Public Types | |
using | elem_vec_t = typename SubElementVectorProvider< SCALAR >::elem_vec_t |
inherited types for element vectors More... | |
using | ElemVec = typename SubElementVectorProvider< SCALAR >::ElemVec |
![]() | |
using | size_type = lf::uscalfe::size_type |
using | elem_vec_t = Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > |
internal type for element vectors More... | |
using | ElemVec = elem_vec_t |
Public Member Functions | |
LoadElementVectorProvider (const LoadElementVectorProvider &)=delete | |
LoadElementVectorProvider (LoadElementVectorProvider &&) noexcept=default | |
LoadElementVectorProvider & | operator= (const LoadElementVectorProvider &)=delete |
LoadElementVectorProvider & | operator= (LoadElementVectorProvider &&)=delete |
LoadElementVectorProvider (std::shared_ptr< ProductUniformFESpace< SCALAR > > fe_space_test, size_type test_component, FUNCTOR f) | |
Constructor: performs cell-independent precomputations. More... | |
ElemVec | Eval (const lf::mesh::Entity &cell) override |
main routine for the computation of element vectors More... | |
size_type | TestComponent () const override |
returns the index of the test space component \( v_{j_k} \) which is the test space for the linear form \( l_k \) More... | |
~LoadElementVectorProvider () override=default | |
![]() | |
SubElementVectorProvider ()=default | |
virtual | ~SubElementVectorProvider ()=default |
SubElementVectorProvider (const SubElementVectorProvider &)=delete | |
SubElementVectorProvider (SubElementVectorProvider &&) noexcept=default | |
SubElementVectorProvider & | operator= (const SubElementVectorProvider &)=delete |
SubElementVectorProvider & | operator= (SubElementVectorProvider &&)=delete |
virtual bool | isActive (const lf::mesh::Entity &) |
All cells are considered active in the default implementation. More... | |
virtual ElemVec | Eval (const lf::mesh::Entity &cell)=0 |
main routine for the computation of (sub) element vectors 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 linear form \( l_k \) More... | |
Private Attributes | |
FUNCTOR | f_ |
functor providing the source function More... | |
std::array< lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, 5 > | fe_precomp_test_ |
array containing precomputed information for the test space: fe_precomp_test_[i] contains the precomputed reference finite element for \( v\) on ref_el i. More... | |
size_type | test_component_ |
index of test component v More... | |
Class for local quadrature based computations of sub vectors corresponding to load vectors.
SCALAR | type for the entries of the element vectors. Must be a fiel type such as 'double'. |
FUNCTOR | a MeshFunction that defines the source function \( f \). It should be scalar-valued. |
The load vector corresponds to the local linear form
\[ v \mapsto \int_K f(\mathbf{x})\,v\,\mathrm{d}\mathbf{x}\;, \]
with source function \( f \). \(v \) is a component of the test space.
operator (const Entity &,ref_coord_t)
that returns a scalar Definition at line 1330 of file loc_comp_dpg.h.
using projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::elem_vec_t = typename SubElementVectorProvider<SCALAR>::elem_vec_t |
inherited types for element vectors
Definition at line 1335 of file loc_comp_dpg.h.
using projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::ElemVec = typename SubElementVectorProvider<SCALAR>::ElemVec |
Definition at line 1336 of file loc_comp_dpg.h.
|
delete |
|
defaultnoexcept |
projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::LoadElementVectorProvider | ( | std::shared_ptr< ProductUniformFESpace< SCALAR > > | fe_space_test, |
size_type | test_component, | ||
FUNCTOR | f | ||
) |
Constructor: performs cell-independent precomputations.
fe_space_test | collection of specifications for the fe space |
test_component | Index of the test space component \(v\) |
f | mesh function for the scalar valued source function |
This constructor uses local quadrature rules with twice the degree of exactnes of the polynomial degree of \( v\)
Definition at line 1393 of file loc_comp_dpg.h.
References projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::fe_precomp_test_, lf::base::RefEl::kQuad(), lf::base::RefEl::kTria(), lf::quad::make_QuadRule(), and projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::test_component_.
|
overridedefault |
|
overridevirtual |
main routine for the computation of element vectors
cell | reference to a (triangular or quadrilateral) cell for which the element vector should be computed |
Implements projects::dpg::SubElementVectorProvider< SCALAR >.
Definition at line 1417 of file loc_comp_dpg.h.
References lf::geometry::Geometry::DimLocal(), lf::mesh::Entity::Geometry(), lf::base::RefEl::Id(), lf::geometry::Geometry::IntegrationElement(), and lf::mesh::Entity::RefEl().
|
delete |
|
delete |
|
inlineoverridevirtual |
returns the index of the test space component \( v_{j_k} \) which is the test space for the linear form \( l_k \)
Implements projects::dpg::SubElementVectorProvider< SCALAR >.
Definition at line 1366 of file loc_comp_dpg.h.
References projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::test_component_.
|
private |
functor providing the source function
Definition at line 1374 of file loc_comp_dpg.h.
|
private |
array containing precomputed information for the test space: fe_precomp_test_[i] contains the precomputed reference finite element for \( v\) on ref_el i.
Definition at line 1379 of file loc_comp_dpg.h.
Referenced by projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::LoadElementVectorProvider().
|
private |
index of test component v
Definition at line 1382 of file loc_comp_dpg.h.
Referenced by projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::LoadElementVectorProvider(), and projects::dpg::LoadElementVectorProvider< SCALAR, FUNCTOR >::TestComponent().