DUNE PDELab (git)

globaldofindex.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_PDELAB_COMMON_GLOBALDOFINDEX_HH
4#define DUNE_PDELAB_COMMON_GLOBALDOFINDEX_HH
5
6#include <dune/pdelab/common/multiindex.hh>
7
8namespace Dune {
9
10 namespace PDELab {
11
12
13 template<typename T, std::size_t tree_n, typename ID>
14 class GlobalDOFIndex
15 {
16
17 public:
18
20 static const std::size_t max_depth = tree_n;
21
22 typedef ID EntityID;
23 typedef MultiIndex<T,max_depth> TreeIndex;
24
25 typedef typename TreeIndex::size_type size_type;
26 typedef T value_type;
27
28 class View
29 {
30
31 friend class GlobalDOFIndex;
32
33 public:
34
35 static const std::size_t max_depth = tree_n;
36
37 typedef ID EntityID;
38 typedef typename MultiIndex<T,tree_n>::View TreeIndex;
39
40 const EntityID& entityID() const
41 {
42 return _entity_id;
43 }
44
45 const TreeIndex& treeIndex() const
46 {
47 return _tree_index_view;
48 }
49
50 View back_popped() const
51 {
52 return View(_entity_id,_tree_index_view.back_popped());
53 }
54
55 friend std::ostream& operator<< (std::ostream& s, const View& di)
56 {
57 s << "(";
58
59 s << di._entity_id;
60
61 s << " | "
62 << di._tree_index_view
63 << ")";
64 return s;
65 }
66
67 std::size_t size() const
68 {
69 return _tree_index_view.size();
70 };
71
72 private:
73
74 explicit View(const GlobalDOFIndex& dof_index)
75 : _entity_id(dof_index._entity_id)
76 , _tree_index_view(dof_index._tree_index.view())
77 {}
78
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))
82 {}
83
84 View(const EntityID& entity_id, const TreeIndex& tree_index)
85 : _entity_id(entity_id)
86 , _tree_index_view(tree_index)
87 {}
88
89 const EntityID& _entity_id;
90 TreeIndex _tree_index_view;
91
92 };
93
95 GlobalDOFIndex()
96 {}
97
98 GlobalDOFIndex(const EntityID& entity_id, const TreeIndex& tree_index)
99 : _entity_id(entity_id)
100 , _tree_index(tree_index)
101 {}
102
103 View view() const
104 {
105 return View(*this);
106 }
107
108 View view(std::size_t size) const
109 {
110 return View(*this,size);
111 }
112
113 void clear()
114 {
115 _entity_id = EntityID();
116 _tree_index.clear();
117 }
118
120 EntityID& entityID()
121 {
122 return _entity_id;
123 }
124
125 const EntityID& entityID() const
126 {
127 return _entity_id;
128 }
129
130 TreeIndex& treeIndex()
131 {
132 return _tree_index;
133 }
134
135 const TreeIndex& treeIndex() const
136 {
137 return _tree_index;
138 }
139
141 friend std::ostream& operator<< (std::ostream& s, const GlobalDOFIndex& di)
142 {
143 s << "(";
144
145 s << di._entity_id;
146
147 s << " | "
148 << di._tree_index
149 << ")";
150 return s;
151 }
152
154
157 bool operator== (const GlobalDOFIndex& r) const
158 {
159 return
160 _entity_id == r._entity_id && _tree_index == r._tree_index;
161 }
162
164 bool operator!= (const GlobalDOFIndex& r) const
165 {
166 return !(*this == r);
167 }
168
169
170 std::size_t size() const
171 {
172 return _tree_index.size();
173 }
174
175 private:
176
177 EntityID _entity_id;
178 TreeIndex _tree_index;
179
180 };
181
182 template<typename T, std::size_t n, typename ID>
183 inline std::size_t hash_value(const GlobalDOFIndex<T,n,ID>& di)
184 {
185 std::size_t seed = 0;
186 std::hash<ID> id_hasher;
187 hash_combine(seed,id_hasher(di.entityID()));
188 hash_range(seed,di.treeIndex().begin(),di.treeIndex().end());
189 return seed;
190 }
191
192
193
194 } // namespace PDELab
195} // namespace Dune
196
197DUNE_DEFINE_HASH(DUNE_HASH_TEMPLATE_ARGS(typename T, std::size_t n, typename ID),DUNE_HASH_TYPE(Dune::PDELab::GlobalDOFIndex<T,n,ID>))
198
199#endif // DUNE_PDELAB_COMMON_GLOBALDOFINDEX_HH
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 23, 23:29, 2024)