DUNE PDELab (2.8)

Dune::PDELab::LocalFiniteElementMapInterface< T, Imp > Class Template Reference

interface for a finite element map More...

#include <dune/pdelab/finiteelementmap/finiteelementmap.hh>

Public Types

typedef T Traits
 Export traits.
 

Public Member Functions

template<class EntityType >
const Traits::FiniteElementType & find (const EntityType &e) const =delete
 Return local basis for the given entity. More...
 

Static Public Member Functions

static constexpr bool hasDOFs (int codim)=delete
 return if FiniteElementMap has degrees of freedom for given codimension
 
static constexpr std::size_t maxLocalSize ()=delete
 compute an upper bound for the local number of DOFs. More...
 

Size calculation

The FiniteElementMap provides different methods to compute the size of the GridFunctionSpace (if possible) without iterating the grid. The approach is as follows (pseudo code):

computeNumberOfDofs(GridView, FEM):
if(FEM.fixedSize()):
sum(FEM.size(gt)*GridView.size(gt) for gt in GeometryTypes)
else
sum(FEM.find(E).basis().size() for E in GridView.entities<0>())
Grid view abstract base class.
Definition: gridview.hh:63
static constexpr bool fixedSize()=delete
a FiniteElementMap is fixedSize iif the size of the local functions space for each GeometryType is fi...
const Traits::FiniteElementType & find(const EntityType &e) const =delete
Return local basis for the given entity.
static constexpr std::size_t size(GeometryType gt)=delete
if the FiniteElementMap is fixedSize, the size methods computes the number of DOFs for given Geometry...
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:156
static constexpr int dimension = -1
 dimension of the domain this FEM is defined on.
 
static constexpr bool fixedSize ()=delete
 a FiniteElementMap is fixedSize iif the size of the local functions space for each GeometryType is fixed.
 
static constexpr std::size_t size (GeometryType gt)=delete
 if the FiniteElementMap is fixedSize, the size methods computes the number of DOFs for given GeometryType.
 

Detailed Description

template<class T, class Imp>
class Dune::PDELab::LocalFiniteElementMapInterface< T, Imp >

interface for a finite element map

Member Function Documentation

◆ find()

template<class T , class Imp >
template<class EntityType >
const Traits::FiniteElementType & Dune::PDELab::LocalFiniteElementMapInterface< T, Imp >::find ( const EntityType &  e) const
delete

Return local basis for the given entity.

The return value is a reference to Traits::LocalBasisType. If there is a different local basis for two elements then this type must be polymorphic.

◆ maxLocalSize()

template<class T , class Imp >
static constexpr std::size_t Dune::PDELab::LocalFiniteElementMapInterface< T, Imp >::maxLocalSize ( )
staticconstexprdelete

compute an upper bound for the local number of DOFs.

this upper bound is used to avoid reallocations in std::vectors used during the assembly.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)