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

cartesian/prodcut space of finite element functions on a hybrid 2D mesh. More...

#include <projects/dpg/product_fe_space.h>

Public Types

using scalar = SCALAR
 

Public Member Functions

 ProductUniformFESpace ()=default
 
 ProductUniformFESpace (const ProductUniformFESpace &)=delete
 
 ProductUniformFESpace (ProductUniformFESpace &&) noexcept=default
 
ProductUniformFESpaceoperator= (const ProductUniformFESpace &)=delete
 
ProductUniformFESpaceoperator= (ProductUniformFESpace &&) noexcept=default
 
 ProductUniformFESpace (std::shared_ptr< const lf::mesh::Mesh > mesh_p, std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_tria_v, std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_quad_v, std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_edge_v, std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_point_v)
 Main construcotr: sets up the local-to-global index mapping (dofhandler) More...
 
std::shared_ptr< const lf::mesh::MeshMesh () const
 acess to underlying mesh More...
 
const ProductUniformFEDofHandlerLocGlobMap () const
 access to associated local-to-global map More...
 
lf::fe::ScalarReferenceFiniteElement< SCALAR > const * ShapeFunctionLayout (lf::base::RefEl ref_el_type, size_type component) const
 access to shape function layout More...
 
size_type NumRefShapeFunctions (lf::base::RefEl ref_el_type, size_type component) const
 number of interior shape functions associated to entities of various types for a given component More...
 
size_type NumComponents () const
 number of components of the product space More...
 
std::shared_ptr< const lf::uscalfe::UniformScalarFESpace< SCALAR > > ComponentFESpace (size_type component) const
 Constructs a lf::uscalfe::UniformScalarFESpace that describes the finite element space associated to a certain component. More...
 
virtual ~ProductUniformFESpace ()=default
 no special destructor More...
 

Private Member Functions

void init ()
 

Private Attributes

std::shared_ptr< const lf::mesh::Meshmesh_p_
 
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_tria_v_
 
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_quad_v_
 
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_edge_v_
 
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > > rfs_point_v_
 
std::vector< size_typenum_rsf_node_
 
std::vector< size_typenum_rsf_edge_
 
std::vector< size_typenum_rsf_tria_
 
std::vector< size_typenum_rsf_quad_
 
std::unique_ptr< ProductUniformFEDofHandlerdofh_p_
 
size_type numComponents_ = 0
 

Detailed Description

template<typename SCALAR>
class projects::dpg::ProductUniformFESpace< SCALAR >

cartesian/prodcut space of finite element functions on a hybrid 2D mesh.

Template Parameters
SCALARunderlying scalar type for all components, usually 'double'

This class provides functionality similar to the class lf::uscalfe::UniformScalarFESpace for product/cartesian finite element spaces. These spaces are of the form

\[ U = U_0 \times U_1 \times \dots \times U_{n-1} \]

where each of the \( U_i \) is a function space of scalar valued functions (see the description in ProductUniformFEDofHandler).

To allow this, this class provides some further member functions that give extra information about the contained FE spaces. It furthermore uses a ProductUniformFEDofHandler internally resulting in a different dof ordering.

Furthermore, this class also weakens some constraints on the shape function layouts compared to the lf::uscalfe::UniformScalarFESpace in order to represent function spaces as they occur in the description of DPG methods:

Note
Some of the pointers may be NULL, which indicates that the shape function descriptions are missing.

Definition at line 58 of file product_fe_space.h.

Member Typedef Documentation

◆ scalar

template<typename SCALAR >
using projects::dpg::ProductUniformFESpace< SCALAR >::scalar = SCALAR

Definition at line 60 of file product_fe_space.h.

Constructor & Destructor Documentation

◆ ProductUniformFESpace() [1/4]

template<typename SCALAR >
projects::dpg::ProductUniformFESpace< SCALAR >::ProductUniformFESpace ( )
default

default constructor

Note
creates an invalid object that can not be used.

◆ ProductUniformFESpace() [2/4]

template<typename SCALAR >
projects::dpg::ProductUniformFESpace< SCALAR >::ProductUniformFESpace ( const ProductUniformFESpace< SCALAR > &  )
delete

◆ ProductUniformFESpace() [3/4]

template<typename SCALAR >
projects::dpg::ProductUniformFESpace< SCALAR >::ProductUniformFESpace ( ProductUniformFESpace< SCALAR > &&  )
defaultnoexcept

◆ ProductUniformFESpace() [4/4]

template<typename SCALAR >
projects::dpg::ProductUniformFESpace< SCALAR >::ProductUniformFESpace ( std::shared_ptr< const lf::mesh::Mesh mesh_p,
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > >  rfs_tria_v,
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > >  rfs_quad_v,
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > >  rfs_edge_v,
std::vector< std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > >  rfs_point_v 
)
inline

Main construcotr: sets up the local-to-global index mapping (dofhandler)

Parameters
mesh_pshared pointer to underlying mesh (imutable)
rfs_tria_vvector of shared pointers to layout description of reference shape functions on triangular cells for each component.
rfs_quad_vvector of shared pointers to layout description of reference shape functions on quadrilateral cells for each component.
rfs_edge_vvector of shared pointers to layout description of reference shape functions on the edges for each component.

The schemes for local shape functions belonging to some component have to satisfy the following compatibility conditions.

  • The number of interior shape functions for edges of triangles and quadrilaterals must agree.
Note
none of the shape function layouts needs to be specified; just pass a null pointer. This will restrict the applicability of the resulting finite element space object.

Definition at line 92 of file product_fe_space.h.

References projects::dpg::ProductUniformFESpace< SCALAR >::init(), projects::dpg::ProductUniformFESpace< SCALAR >::numComponents_, and projects::dpg::ProductUniformFESpace< SCALAR >::rfs_tria_v_.

◆ ~ProductUniformFESpace()

template<typename SCALAR >
virtual projects::dpg::ProductUniformFESpace< SCALAR >::~ProductUniformFESpace ( )
virtualdefault

no special destructor

Member Function Documentation

◆ ComponentFESpace()

template<typename SCALAR >
std::shared_ptr< const lf::uscalfe::UniformScalarFESpace< SCALAR > > projects::dpg::ProductUniformFESpace< SCALAR >::ComponentFESpace ( size_type  component) const
inline

Constructs a lf::uscalfe::UniformScalarFESpace that describes the finite element space associated to a certain component.

Note
The resulting FESpace has no more information about the other components
Warning
For some components the construction of the lf::uscalfe::UniformScalarFESpace may be impossible, since some of the constraints regarding the shape function layouts were weakened.

Definition at line 163 of file product_fe_space.h.

References projects::dpg::ProductUniformFESpace< SCALAR >::mesh_p_, projects::dpg::ProductUniformFESpace< SCALAR >::rfs_edge_v_, projects::dpg::ProductUniformFESpace< SCALAR >::rfs_point_v_, projects::dpg::ProductUniformFESpace< SCALAR >::rfs_quad_v_, and projects::dpg::ProductUniformFESpace< SCALAR >::rfs_tria_v_.

Referenced by projects::dpg::InitEssentialConditionsFromFunctions().

◆ init()

template<typename SCALAR >
void projects::dpg::ProductUniformFESpace< SCALAR >::init
private

Initialization of class member variables and consistency checks

Definition at line 210 of file product_fe_space.h.

References lf::base::RefEl::kPoint(), lf::base::RefEl::kQuad(), lf::base::RefEl::kSegment(), and lf::base::RefEl::kTria().

Referenced by projects::dpg::ProductUniformFESpace< SCALAR >::ProductUniformFESpace().

◆ LocGlobMap()

template<typename SCALAR >
const ProductUniformFEDofHandler & projects::dpg::ProductUniformFESpace< SCALAR >::LocGlobMap ( ) const
inline

access to associated local-to-global map

Returns
a reference to the ProductUniformFEDofHandler object (immutable)

Definition at line 126 of file product_fe_space.h.

References projects::dpg::ProductUniformFESpace< SCALAR >::dofh_p_, and projects::dpg::ProductUniformFESpace< SCALAR >::mesh_p_.

Referenced by projects::dpg::InitEssentialConditionsFromFunctions().

◆ Mesh()

template<typename SCALAR >
std::shared_ptr< const lf::mesh::Mesh > projects::dpg::ProductUniformFESpace< SCALAR >::Mesh ( ) const
inline

acess to underlying mesh

Returns
a shared pointer to the mesh

Definition at line 118 of file product_fe_space.h.

References projects::dpg::ProductUniformFESpace< SCALAR >::mesh_p_.

◆ NumComponents()

template<typename SCALAR >
size_type projects::dpg::ProductUniformFESpace< SCALAR >::NumComponents ( ) const
inline

number of components of the product space

Definition at line 152 of file product_fe_space.h.

References projects::dpg::ProductUniformFESpace< SCALAR >::numComponents_.

◆ NumRefShapeFunctions()

template<typename SCALAR >
size_type projects::dpg::ProductUniformFESpace< SCALAR >::NumRefShapeFunctions ( lf::base::RefEl  ref_el_type,
size_type  component 
) const

number of interior shape functions associated to entities of various types for a given component

Definition at line 380 of file product_fe_space.h.

References lf::base::RefEl::kPoint(), lf::base::RefEl::kQuad(), lf::base::RefEl::kSegment(), and lf::base::RefEl::kTria().

◆ operator=() [1/2]

template<typename SCALAR >
ProductUniformFESpace & projects::dpg::ProductUniformFESpace< SCALAR >::operator= ( const ProductUniformFESpace< SCALAR > &  )
delete

◆ operator=() [2/2]

template<typename SCALAR >
ProductUniformFESpace & projects::dpg::ProductUniformFESpace< SCALAR >::operator= ( ProductUniformFESpace< SCALAR > &&  )
defaultnoexcept

◆ ShapeFunctionLayout()

template<typename SCALAR >
lf::fe::ScalarReferenceFiniteElement< SCALAR > const * projects::dpg::ProductUniformFESpace< SCALAR >::ShapeFunctionLayout ( lf::base::RefEl  ref_el_type,
size_type  component 
) const

access to shape function layout

Parameters
ref_el_typetype of entit, can be anything except for lf::base::RefEl::kPoint()
componentindex of the component whose shape function layout is returned.
Warning
NULL pointers may be returned by this method in case a finite element specification was not given for a particular topological type of entity and e particular component.

Definition at line 348 of file product_fe_space.h.

References lf::base::RefEl::kPoint(), lf::base::RefEl::kQuad(), lf::base::RefEl::kSegment(), and lf::base::RefEl::kTria().

Member Data Documentation

◆ dofh_p_

template<typename SCALAR >
std::unique_ptr<ProductUniformFEDofHandler> projects::dpg::ProductUniformFESpace< SCALAR >::dofh_p_
private

local-to-global index map for the finite element space

Definition at line 199 of file product_fe_space.h.

Referenced by projects::dpg::ProductUniformFESpace< SCALAR >::LocGlobMap().

◆ mesh_p_

template<typename SCALAR >
std::shared_ptr<const lf::mesh::Mesh> projects::dpg::ProductUniformFESpace< SCALAR >::mesh_p_
private

◆ num_rsf_edge_

template<typename SCALAR >
std::vector<size_type> projects::dpg::ProductUniformFESpace< SCALAR >::num_rsf_edge_
private

Definition at line 194 of file product_fe_space.h.

◆ num_rsf_node_

template<typename SCALAR >
std::vector<size_type> projects::dpg::ProductUniformFESpace< SCALAR >::num_rsf_node_
private

numbers of local shape functions for all components on different types of entities.

Definition at line 193 of file product_fe_space.h.

◆ num_rsf_quad_

template<typename SCALAR >
std::vector<size_type> projects::dpg::ProductUniformFESpace< SCALAR >::num_rsf_quad_
private

Definition at line 196 of file product_fe_space.h.

◆ num_rsf_tria_

template<typename SCALAR >
std::vector<size_type> projects::dpg::ProductUniformFESpace< SCALAR >::num_rsf_tria_
private

Definition at line 195 of file product_fe_space.h.

◆ numComponents_

template<typename SCALAR >
size_type projects::dpg::ProductUniformFESpace< SCALAR >::numComponents_ = 0
private

◆ rfs_edge_v_

template<typename SCALAR >
std::vector< std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR> > > projects::dpg::ProductUniformFESpace< SCALAR >::rfs_edge_v_
private

◆ rfs_point_v_

template<typename SCALAR >
std::vector< std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR> > > projects::dpg::ProductUniformFESpace< SCALAR >::rfs_point_v_
private

◆ rfs_quad_v_

template<typename SCALAR >
std::vector< std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR> > > projects::dpg::ProductUniformFESpace< SCALAR >::rfs_quad_v_
private

◆ rfs_tria_v_

template<typename SCALAR >
std::vector< std::shared_ptr<const lf::fe::ScalarReferenceFiniteElement<SCALAR> > > projects::dpg::ProductUniformFESpace< SCALAR >::rfs_tria_v_
private

descritpions of reference shape functions for all components on different types of entities.

Definition at line 180 of file product_fe_space.h.

Referenced by projects::dpg::ProductUniformFESpace< SCALAR >::ComponentFESpace(), and projects::dpg::ProductUniformFESpace< SCALAR >::ProductUniformFESpace().


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