DUNE PDELab (2.8)

orderingbase.hh
1// -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=8 sw=2 sts=2:
3
4#ifndef DUNE_PDELAB_ORDERING_ORDERINGBASE_HH
5#define DUNE_PDELAB_ORDERING_ORDERINGBASE_HH
6
8#include <dune/pdelab/ordering/utility.hh>
9#include <dune/pdelab/gridfunctionspace/gridfunctionspacebase.hh>
10
11#include <vector>
12
13namespace Dune {
14 namespace PDELab {
15
18
19 template<typename DI, typename CI>
20 class OrderingBase
21 {
22
23 template<typename size_type>
24 friend struct ::Dune::PDELab::impl::update_ordering_data;
25
26 public:
27
28 typedef OrderingTraits<DI,CI,MultiIndexOrder::Inner2Outer> Traits;
29
30 protected:
31
32 typedef Dune::PDELab::impl::GridFunctionSpaceOrderingData<typename Traits::SizeType> GFSData;
33
34 public:
35
36 typedef HierarchicContainerAllocationTag ContainerAllocationTag;
37
38 typedef DefaultLFSCacheTag CacheTag;
39
40 static const bool has_dynamic_ordering_children = true;
41
42 static const bool consume_tree_index = true;
43
44 typename Traits::ContainerIndex mapIndex(const typename Traits::DOFIndex& di) const
45 {
46 typename Traits::ContainerIndex ci;
47 mapIndex(di.view(),ci);
48 return ci;
49 }
50
51 void mapIndex(typename Traits::DOFIndexView di, typename Traits::ContainerIndex& ci) const
52 {
53 if (_delegate)
54 _delegate->map_index_dynamic(di,ci);
55 else
56 {
57 _mapIndex(di,ci);
58 }
59 }
60
61 typename Traits::SizeType size() const
62 {
63 return _size;
64 }
65
66 typename Traits::SizeType blockCount() const
67 {
68 return _block_count;
69 }
70
71 typename Traits::SizeType size(const typename Traits::SizeType child_index) const
72 {
73 return _child_size_offsets[child_index + 1] - _child_size_offsets[child_index];
74 }
75
76 typename Traits::SizeType sizeOffset(const typename Traits::SizeType child_index) const
77 {
78 return _child_size_offsets[child_index];
79 }
80
81 typename Traits::SizeType blockOffset(const typename Traits::SizeType child_index) const
82 {
83 assert(_merge_mode == MergeMode::lexicographic);
84 return _child_block_offsets[child_index];
85 }
86
87 typename Traits::SizeType maxLocalSize() const
88 {
89 return _max_local_size;
90 }
91
92 MergeMode::type mergeMode() const
93 {
94 return _merge_mode;
95 }
96
97 void update()
98 {
99 std::fill(_child_size_offsets.begin(),_child_size_offsets.end(),0);
100 std::fill(_child_block_offsets.begin(),_child_block_offsets.end(),0);
101 _codim_used.reset();
102 _codim_fixed_size.set();
103 _max_local_size = 0;
104 _block_count = 0;
105 typename Traits::SizeType block_carry = 0;
106 typename Traits::SizeType size_carry = 0;
107 for (typename Traits::SizeType i = 0; i < _child_count; ++i)
108 {
109 _child_block_offsets[i+1] = (block_carry += _children[i]->blockCount());
110 _child_size_offsets[i+1] = (size_carry += _children[i]->size());
111 _codim_used |= _children[i]->_codim_used;
112 _codim_fixed_size &= _children[i]->_codim_fixed_size;
113 _block_count += _children[i]->blockCount();
114 _max_local_size += _children[i]->maxLocalSize();
115 }
116 if (_container_blocked)
117 if (_merge_mode == MergeMode::lexicographic)
118 {
119 _block_count = _child_count;
120 }
121 else
122 {
123 if (_child_block_offsets.back() % _child_block_merge_offsets.back() != 0)
124 DUNE_THROW(OrderingStructureError,
125 "Invalid ordering structure: "
126 << "total number of blocks ("
127 << _child_block_offsets.back()
128 << ") is not a multiple of the interleaved block size ("
129 << _child_block_merge_offsets.back()
130 << ")."
131 );
132 _block_count = _child_block_offsets.back() / _child_block_merge_offsets.back();
133 }
134 else
135 _block_count = _child_block_offsets.back();
136 _size = _child_size_offsets.back();
137 }
138
139 template<typename Node>
140 OrderingBase(Node& node,
141 bool container_blocked,
142 GFSData* gfs_data,
143 VirtualOrderingBase<DI,CI>* delegate = nullptr)
144 : _fixed_size(false)
145 , _container_blocked(container_blocked)
146 , _merge_mode(MergeMode::lexicographic)
147 , _child_count(Node::has_dynamic_ordering_children ? TypeTree::degree(node) : 0)
148 , _children(_child_count,nullptr)
149 , _child_size_offsets((_child_count > 0 ? _child_count + 1 : 0),0)
150 , _child_block_offsets((_child_count > 0 ? _child_count + 1 : 0),0)
151 , _max_local_size(0)
152 , _size(0)
153 , _block_count(0)
154 , _delegate(delegate)
155 , _gfs_data(gfs_data)
156 {
157 TypeTree::applyToTree(node,extract_child_bases<OrderingBase>(_children));
158 }
159
160 template<typename Node>
161 OrderingBase(Node& node,
162 bool container_blocked,
163 const std::vector<std::size_t>& merge_offsets,
164 GFSData* gfs_data,
165 VirtualOrderingBase<DI,CI>* delegate = nullptr)
166 : _fixed_size(false)
167 , _container_blocked(container_blocked)
168 , _merge_mode(MergeMode::interleaved)
169 , _child_count(Node::has_dynamic_ordering_children ? TypeTree::degree(node) : 0)
170 , _children(_child_count,nullptr)
171 , _child_size_offsets((_child_count > 0 ? _child_count + 1 : 0),0)
172 , _child_block_offsets((_child_count > 0 ? _child_count + 1 : 0),0)
173 , _child_block_merge_offsets(merge_offsets)
174 , _max_local_size(0)
175 , _size(0)
176 , _block_count(0)
177 , _delegate(delegate)
178 , _gfs_data(gfs_data)
179 {
180 TypeTree::applyToTree(node,extract_child_bases<OrderingBase>(_children));
181 }
182
183
184 bool containerBlocked() const
185 {
186 return _container_blocked;
187 }
188
189 std::size_t childOrderingCount() const
190 {
191 return _child_count;
192 }
193
194 OrderingBase& childOrdering(typename Traits::SizeType i)
195 {
196 return *_children[i];
197 }
198
199 const OrderingBase& childOrdering(typename Traits::SizeType i) const
200 {
201 return *_children[i];
202 }
203
204 bool contains(typename Traits::SizeType codim) const
205 {
206 return _codim_used.test(codim);
207 }
208
209 bool fixedSize() const
210 {
211 return _fixed_size;
212 }
213
214 bool fixedSize(typename Traits::SizeType codim) const
215 {
216 return _codim_fixed_size.test(codim);
217 }
218
219 protected:
220
222
227 void setDelegate(const VirtualOrderingBase<DI,CI>* delegate)
228 {
229 _delegate = delegate;
230 }
231
232 void _mapIndex(typename Traits::DOFIndexView di, typename Traits::ContainerIndex& ci) const
233 {
234 typedef typename Traits::SizeType size_type;
235 size_type child_index = di.treeIndex().back();
236 _children[child_index]->mapIndex(di.back_popped(),ci);
237 if (_merge_mode == MergeMode::lexicographic)
238 {
239 if (_container_blocked)
240 ci.push_back(child_index);
241 else
242 ci.back() += blockOffset(child_index);
243 }
244 else
245 {
246 size_type child_block_offset = _child_block_merge_offsets[child_index];
247 size_type child_block_size = _child_block_merge_offsets[child_index + 1] - child_block_offset;
248 size_type block_index = ci.back() / child_block_size;
249 size_type offset = ci.back() % child_block_size;
250 size_type block_offset = child_block_offset + offset;
251 if (_container_blocked)
252 {
253 ci.back() = block_offset;
254 ci.push_back(block_index);
255 }
256 else
257 {
258 size_type block_size = _child_block_merge_offsets.back();
259 ci.back() = block_index * block_size + block_offset;
260 }
261 }
262 }
263
264 private:
265
266 bool update_gfs_data_size(typename Traits::SizeType& size, typename Traits::SizeType& block_count) const
267 {
268 size = _size;
269 block_count = _block_count;
270 return true;
271 }
272
273 public:
274
275 bool _fixed_size;
276 const bool _container_blocked;
277 const MergeMode::type _merge_mode;
278
279 const std::size_t _child_count;
280 std::vector<OrderingBase*> _children;
281
282 std::vector<typename Traits::SizeType> _child_size_offsets;
283 std::vector<typename Traits::SizeType> _child_block_offsets;
284 std::vector<typename Traits::SizeType> _child_block_merge_offsets;
285 typename Traits::CodimFlag _codim_used;
286 typename Traits::CodimFlag _codim_fixed_size;
287
288 std::size_t _max_local_size;
289 std::size_t _size;
290 std::size_t _block_count;
291
292 const VirtualOrderingBase<DI,CI>* _delegate;
293 GFSData* _gfs_data;
294
295 };
296
298
299 } // namespace PDELab
300} // namespace Dune
301
302#endif // DUNE_PDELAB_ORDERING_ORDERINGBASE_HH
PDELab-specific exceptions.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
std::size_t degree(const Node &node)
Returns the degree of node as run time information.
Definition: nodeinterface.hh:76
void applyToTree(Tree &&tree, Visitor &&visitor)
Apply visitor to TypeTree.
Definition: traversal.hh:237
Dune namespace.
Definition: alignedallocator.hh:11
type
Definition: utility.hh:24
@ lexicographic
Lexicographically ordered ([i1,i2],[j1,j2] -> [i1,i2,j1,j2]).
Definition: utility.hh:25
std::size_t fixedSize
The number of data items per index if it is fixed, 0 otherwise.
Definition: variablesizecommunicator.hh:272
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)