3#ifndef DUNE_PDELAB_COMMON_GLOBALDOFINDEX_HH
4#define DUNE_PDELAB_COMMON_GLOBALDOFINDEX_HH
6#include <dune/pdelab/common/multiindex.hh>
13 template<
typename T, std::
size_t tree_n,
typename ID>
20 static const std::size_t max_depth = tree_n;
23 typedef MultiIndex<T,max_depth> TreeIndex;
31 friend class GlobalDOFIndex;
35 static const std::size_t max_depth = tree_n;
38 typedef typename MultiIndex<T,tree_n>::View TreeIndex;
40 const EntityID& entityID()
const
45 const TreeIndex& treeIndex()
const
47 return _tree_index_view;
50 View back_popped()
const
52 return View(_entity_id,_tree_index_view.back_popped());
55 friend std::ostream& operator<< (std::ostream& s,
const View& di)
62 << di._tree_index_view
67 std::size_t
size()
const
69 return _tree_index_view.size();
74 explicit View(
const GlobalDOFIndex& dof_index)
75 : _entity_id(dof_index._entity_id)
76 , _tree_index_view(dof_index._tree_index.view())
79 View(
const GlobalDOFIndex& dof_index, std::size_t
size)
80 : _entity_id(dof_index._entity_id)
81 , _tree_index_view(dof_index._tree_index.view(
size))
84 View(
const EntityID& entity_id,
const TreeIndex& tree_index)
85 : _entity_id(entity_id)
86 , _tree_index_view(tree_index)
89 const EntityID& _entity_id;
90 TreeIndex _tree_index_view;
98 GlobalDOFIndex(
const EntityID& entity_id,
const TreeIndex& tree_index)
99 : _entity_id(entity_id)
100 , _tree_index(tree_index)
108 View view(std::size_t
size)
const
110 return View(*
this,
size);
115 _entity_id = EntityID();
125 const EntityID& entityID()
const
130 TreeIndex& treeIndex()
135 const TreeIndex& treeIndex()
const
141 friend std::ostream& operator<< (std::ostream& s,
const GlobalDOFIndex& di)
157 bool operator== (
const GlobalDOFIndex& r)
const
160 _entity_id == r._entity_id && _tree_index == r._tree_index;
164 bool operator!= (
const GlobalDOFIndex& r)
const
166 return !(*
this == r);
170 std::size_t
size()
const
172 return _tree_index.size();
178 TreeIndex _tree_index;
182 template<
typename T, std::
size_t n,
typename ID>
183 inline std::size_t hash_value(
const GlobalDOFIndex<T,n,ID>& di)
185 std::size_t seed = 0;
186 std::hash<ID> id_hasher;
188 hash_range(seed,di.treeIndex().begin(),di.treeIndex().end());
storage_type::size_type size_type
An unsigned integral type.
Definition: reservedvector.hh:65
EnableIfInterOperable< T1, T2, bool >::type operator==(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for equality.
Definition: iteratorfacades.hh:238
EnableIfInterOperable< T1, T2, bool >::type operator!=(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for inequality.
Definition: iteratorfacades.hh:260
#define DUNE_DEFINE_HASH(template_args, type)
Defines the required struct specialization to make type hashable via Dune::hash.
Definition: hash.hh:100
#define DUNE_HASH_TYPE(...)
Wrapper macro for the type to be hashed in DUNE_DEFINE_HASH.
Definition: hash.hh:117
#define DUNE_HASH_TEMPLATE_ARGS(...)
Wrapper macro for the template arguments in DUNE_DEFINE_HASH.
Definition: hash.hh:109
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
std::size_t hash_range(It first, It last)
Hashes all elements in the range [first,last) and returns the combined hash.
Definition: hash.hh:322
void hash_combine(std::size_t &seed, const T &arg)
Calculates the hash value of arg and combines it in-place with seed.
Definition: hash.hh:307