DUNE PDELab (unstable)
single component local function space More...
#include <dune/pdelab/gridfunctionspace/localfunctionspace.hh>
Public Types | |
typedef LeafNodeTag | NodeTag |
The type tag that describes a LeafNode. | |
Public Member Functions | |
template<typename Transformation > | |
LeafLocalFunctionSpaceNode (std::shared_ptr< const GFS > gfs, const Transformation &t) | |
initialize with grid function space | |
const Traits::FiniteElementType & | finiteElement () const |
get finite element | |
const Traits::ConstraintsType & | constraints () const |
get constraints engine | |
template<typename Entity , typename DOFIndexIterator , bool fast> | |
void | dofIndices (const Entity &e, DOFIndexIterator it, DOFIndexIterator endit, std::integral_constant< bool, fast >) |
Calculates the multiindices associated with the given entity. | |
template<bool fast = false> | |
void | bind (const typename Traits::Element &e, std::integral_constant< bool, fast > fast_=std::integral_constant< bool, fast >{}) |
bind local function space to entity | |
template<class FE > | |
void | bindFiniteElement (FE &&fe) |
Binds a finite element to the local space If the finite element is lvalue, the caller (i.e. FEM) must guarantee the lifetime of the object since we only keep a view on it. On the other hand, if it is rvalue, we store it locally but we require the object to be fully movable. More... | |
void | unbindFiniteElement () noexcept |
Release view of the bound finite element. | |
Traits::IndexContainer::size_type | size () const |
number of degrees of freedom contained in this lfs node | |
Traits::IndexContainer::size_type | maxSize () const |
get maximum possible size (which is maxLocalSize from grid function space) | |
Traits::IndexContainer::size_type | localVectorSize () const |
get size of an appropriate local vector object More... | |
Traits::IndexContainer::size_type | localIndex (typename Traits::IndexContainer::size_type index) const |
map index in this local function space to root local function space | |
const Traits::DOFIndex & | dofIndex (typename Traits::IndexContainer::size_type index) const |
Maps given index in this local function space to its corresponding global MultiIndex. More... | |
void | debug () const |
print debug information about this local function space | |
const GFS & | gridFunctionSpace () const |
Returns the GridFunctionSpace underlying this LocalFunctionSpace. | |
Static Public Attributes | |
static const bool | isLeaf = true |
Mark this class as a leaf in a dune-typetree. | |
static const bool | isPower = false |
Mark this class as a non power in the dune-typetree. | |
static const bool | isComposite = false |
Mark this class as a non composite in the dune-typetree. | |
Protected Member Functions | |
template<typename NodeType , bool fast = false> | |
void | bind (NodeType &node, const typename Traits::Element &e, std::integral_constant< bool, fast >=std::integral_constant< bool, fast >{}) |
bind local function space to entity More... | |
Detailed Description
class Dune::PDELab::LeafLocalFunctionSpaceNode< GFS, DOFIndex >
single component local function space
Member Function Documentation
◆ bind()
|
protectedinherited |
bind local function space to entity
This is a generic implementation of the bind function. It is parametrized with the NodeType, which the type of the derived LocalFunctionSpaceNode. Handing the NodeType as a parammeter avoid the need for the CRTP construct, but all derived classes have to add a method bind, which forward to this method.
- Parameters
-
node reference to the derived node, the address must be the same as this e entity to bind to
◆ bindFiniteElement()
|
inline |
Binds a finite element to the local space If the finite element is lvalue, the caller (i.e. FEM) must guarantee the lifetime of the object since we only keep a view on it. On the other hand, if it is rvalue, we store it locally but we require the object to be fully movable.
- Template Parameters
-
FE rvalue or lvalue type of the finite element
- Parameters
-
fe the finite element to be bound
◆ dofIndex()
|
inlineinherited |
Maps given index in this local function space to its corresponding global MultiIndex.
- Parameters
-
index The local index value from the range 0,...,size()-1
- Returns
- A const reference to the associated, globally unique MultiIndex. Note that the returned object may (and must) be copied if it needs to be stored beyond the time of the next call to bind() on this LocalFunctionSpace (e.g. when the MultiIndex is used as a DOF identifier in a constraints container).
◆ localVectorSize()
|
inlineinherited |
get size of an appropriate local vector object
this is the number of dofs of the complete local function space tree, i.e. the size() of the root node. The local vector objects must always have this size and the localIndex method maps into the range [0,localVectorSize()[
The documentation for this class was generated from the following file:
- dune/pdelab/gridfunctionspace/localfunctionspace.hh
