dune-pdelab  2.0.0
Classes | Public Member Functions | Static Public Attributes | Friends | List of all members
Dune::PDELab::MultiIndex< T, n > Class Template Reference

A multi-index representing a degree of freedom in a GridFunctionSpace. More...

#include <dune/pdelab/common/multiindex.hh>

Inheritance diagram for Dune::PDELab::MultiIndex< T, n >:
Inheritance graph

Classes

class  View
 

Public Member Functions

 MultiIndex ()
 
 MultiIndex (const View &view)
 
void set (typename ReservedVector< T, n >::value_type index)
 
View view () const
 
View view (std::size_t size) const
 
bool operator== (const MultiIndex &r) const
 Tests whether two MultiIndices are equal. More...
 
bool operator!= (const MultiIndex &r) const
 Tests whether two MultiIndices are not equal. More...
 

Static Public Attributes

static const std::size_t max_depth = n
 The maximum possible depth of the MultiIndex. More...
 

Friends

std::ostream & operator<< (std::ostream &s, const MultiIndex &mi)
 Writes a pretty representation of the MultiIndex to the given std::ostream. More...
 

Detailed Description

template<typename T, std::size_t n>
class Dune::PDELab::MultiIndex< T, n >

A multi-index representing a degree of freedom in a GridFunctionSpace.

A MultiIndex provides a way for identifying degrees of freedom in a (possibly nested) GridFunctionSpace without imposing any kind of ordering. For that purpose, a MultiIndex identifies a degree of freedom by recording

The length of this index tuple is limited by the template parameter n, which will usually be equal to the maximum depth of the current GridFunctionSpace tree. Moreover, there will never be two identical index tuples associated with the same grid entity.

The index tuple is oriented from left to right when traversing up the tree (i.e. towards the root node) and from right to left when drilling down from the root node towards a leaf. For non-leaf nodes, the associated index entry identifies the child GridFunctionSpace that the degree of freedom is associated with, while for leafs, it provides a way to provide multiple degrees of freedom for a single grid entity (usually, the index value for a leaf space will correspond to the LocalKey::index() value from the finite element).

Note that in general, the length of the index tuple will not be the same for all degrees of freedom in a GridFunctionSpace. Consider the following example of a Taylor-Hood element:

In this case, degrees of freedom for the velocity components will have an index tuple of length 3, while those related to pressure will only have an index tuple of length 2. For the Taylor-Hood space given above, the multiindices associated to a triangle with vertex and edge indices in the range {0,1,2} are

GeometryType entity index index tuple
Point 0 0, 0, 0
Point 1 0, 0, 0
Point 2 0, 0, 0
Line 0 0, 0, 0
Line 1 0, 0, 0
Line 2 0, 0, 0
Point 0 0, 1, 0
Point 1 0, 1, 0
Point 2 0, 1, 0
Line 0 0, 1, 0
Line 1 0, 1, 0
Line 2 0, 1, 0
Point 0 0, 1
Point 1 0, 1
Point 2 0, 1
Template Parameters
Tthe type of the index entries.
nthe maximum number of indices in the MultiIndex.

Constructor & Destructor Documentation

template<typename T, std::size_t n>
Dune::PDELab::MultiIndex< T, n >::MultiIndex ( )
inline
template<typename T, std::size_t n>
Dune::PDELab::MultiIndex< T, n >::MultiIndex ( const View view)
inline

Member Function Documentation

template<typename T, std::size_t n>
bool Dune::PDELab::MultiIndex< T, n >::operator!= ( const MultiIndex< T, n > &  r) const
inline

Tests whether two MultiIndices are not equal.

template<typename T, std::size_t n>
bool Dune::PDELab::MultiIndex< T, n >::operator== ( const MultiIndex< T, n > &  r) const
inline

Tests whether two MultiIndices are equal.

Note
Only MultiIndices of identical max_depth are comparable
template<typename T, std::size_t n>
void Dune::PDELab::MultiIndex< T, n >::set ( typename ReservedVector< T, n >::value_type  index)
inline
template<typename T, std::size_t n>
View Dune::PDELab::MultiIndex< T, n >::view ( ) const
inline
template<typename T, std::size_t n>
View Dune::PDELab::MultiIndex< T, n >::view ( std::size_t  size) const
inline

Friends And Related Function Documentation

template<typename T, std::size_t n>
std::ostream& operator<< ( std::ostream &  s,
const MultiIndex< T, n > &  mi 
)
friend

Writes a pretty representation of the MultiIndex to the given std::ostream.

Member Data Documentation

template<typename T, std::size_t n>
const std::size_t Dune::PDELab::MultiIndex< T, n >::max_depth = n
static

The maximum possible depth of the MultiIndex.


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