3 #ifndef DUNE_PDELAB_COMMON_DOFINDEX_HH
4 #define DUNE_PDELAB_COMMON_DOFINDEX_HH
146 template<
typename T, std::
size_t tree_n, std::
size_t entity_n = 1>
179 return _entity_index_view;
184 return _tree_index_view;
189 return View(_entity_index_view,_tree_index_view.back_popped());
196 for (
typename std::remove_reference<EntityIndex>::type::const_iterator it = di._entity_index_view.begin(); it != di._entity_index_view.end(); ++it)
197 s << std::setw(4) << *it;
200 << di._tree_index_view
207 return _tree_index_view.size();
213 : _entity_index_view(dof_index._entity_index)
214 , _tree_index_view(dof_index._tree_index.
view())
218 : _entity_index_view(dof_index._entity_index)
219 , _tree_index_view(dof_index._tree_index.
view(
size))
223 : _entity_index_view(entity_index)
224 , _tree_index_view(tree_index)
237 : _entity_index(
view._entity_index_view)
238 , _tree_index(
view._tree_index_view)
253 std::fill(_entity_index.begin(),_entity_index.end(),0);
260 return _entity_index;
266 return _entity_index;
286 for (
typename EntityIndex::const_iterator it = di._entity_index.begin(); it != di._entity_index.end(); ++it)
287 s << std::setw(4) << *it;
302 std::equal(_entity_index.begin(),_entity_index.end(),r._entity_index.begin()) &&
303 _tree_index == r._tree_index;
309 return !(*
this == r);
313 bool operator< (
const DOFIndex& r)
const
316 return _c.size() < _r.size();
317 return std::lexicographical_compare(_c.begin(),_c.end(),r._c.begin(),r._c.end());
323 return _tree_index.size();
333 template<
typename T, std::
size_t n1, std::
size_t n2>
336 std::size_t seed = 0;
const std::string s
Definition: function.hh:843
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
std::size_t hash_value(const DOFIndex< T, n1, n2 > &di)
Definition: dofindex.hh:334
A multi-index representing a degree of freedom in a GridFunctionSpace.
Definition: dofindex.hh:148
const EntityIndex & entityIndex() const
Returns the index of the grid entity associated with the DOF.
Definition: dofindex.hh:264
const TreeIndex & treeIndex() const
Returns the tuple of entity-local indices associated with the DOF.
Definition: dofindex.hh:276
EntityIndex & entityIndex()
Returns the index of the grid entity associated with the DOF.
Definition: dofindex.hh:258
std::size_t size() const
Definition: dofindex.hh:321
View view() const
Definition: dofindex.hh:241
bool operator!=(const DOFIndex &r) const
Tests whether two DOFIndices are not equal.
Definition: dofindex.hh:307
static const std::size_t max_depth
The maximum possible length of the tuple of entity-local indices.
Definition: dofindex.hh:156
TreeIndex & treeIndex()
Returns the tuple of entity-local indices associated with the DOF.
Definition: dofindex.hh:270
DOFIndex(const View &view)
Definition: dofindex.hh:236
friend std::ostream & operator<<(std::ostream &s, const DOFIndex &di)
Writes a pretty representation of the DOFIndex to the given std::ostream.
Definition: dofindex.hh:282
std::array< T, entity_capacity > EntityIndex
Definition: dofindex.hh:158
View view(std::size_t size) const
Definition: dofindex.hh:246
MultiIndex< T, max_depth > TreeIndex
Definition: dofindex.hh:159
static const std::size_t entity_capacity
The maximum possible length of a tuple which represents the entity index.
Definition: dofindex.hh:153
T value_type
Definition: dofindex.hh:162
void clear()
Definition: dofindex.hh:251
TreeIndex::size_type size_type
Definition: dofindex.hh:161
DOFIndex()
Default constructor.
Definition: dofindex.hh:233
bool operator==(const DOFIndex &r) const
Tests whether two DOFIndices are equal.
Definition: dofindex.hh:299
Definition: dofindex.hh:165
static const std::size_t max_depth
Definition: dofindex.hh:171
static const std::size_t entity_capacity
Definition: dofindex.hh:172
std::size_t size() const
Definition: dofindex.hh:205
EntityIndex & entityIndex() const
Definition: dofindex.hh:177
friend std::ostream & operator<<(std::ostream &s, const View &di)
Definition: dofindex.hh:192
const TreeIndex & treeIndex() const
Definition: dofindex.hh:182
View back_popped() const
Definition: dofindex.hh:187
MultiIndex< T, tree_n >::View TreeIndex
Definition: dofindex.hh:175
const std::array< T, entity_n > & EntityIndex
Definition: dofindex.hh:174