4#ifndef DUNE_PDELAB_ORDERING_ORDERINGBASE_HH
5#define DUNE_PDELAB_ORDERING_ORDERINGBASE_HH
8#include <dune/pdelab/ordering/utility.hh>
9#include <dune/pdelab/gridfunctionspace/gridfunctionspacebase.hh>
19 template<
typename DI,
typename CI>
23 template<
typename size_type>
24 friend struct ::Dune::PDELab::impl::update_ordering_data;
28 typedef OrderingTraits<DI,CI,MultiIndexOrder::Inner2Outer> Traits;
32 typedef Dune::PDELab::impl::GridFunctionSpaceOrderingData<typename Traits::SizeType> GFSData;
36 typedef HierarchicContainerAllocationTag ContainerAllocationTag;
38 typedef DefaultLFSCacheTag CacheTag;
40 static const bool has_dynamic_ordering_children =
true;
42 static const bool consume_tree_index =
true;
44 typename Traits::ContainerIndex mapIndex(
const typename Traits::DOFIndex& di)
const
46 typename Traits::ContainerIndex ci;
47 mapIndex(di.view(),ci);
51 void mapIndex(
typename Traits::DOFIndexView di,
typename Traits::ContainerIndex& ci)
const
54 _delegate->map_index_dynamic(di,ci);
61 typename Traits::SizeType
size()
const
66 typename Traits::SizeType blockCount()
const
71 typename Traits::SizeType
size(
const typename Traits::SizeType child_index)
const
73 return _child_size_offsets[child_index + 1] - _child_size_offsets[child_index];
76 typename Traits::SizeType sizeOffset(
const typename Traits::SizeType child_index)
const
78 return _child_size_offsets[child_index];
81 typename Traits::SizeType blockOffset(
const typename Traits::SizeType child_index)
const
84 return _child_block_offsets[child_index];
87 typename Traits::SizeType maxLocalSize()
const
89 return _max_local_size;
99 std::fill(_child_size_offsets.begin(),_child_size_offsets.end(),0);
100 std::fill(_child_block_offsets.begin(),_child_block_offsets.end(),0);
102 _codim_fixed_size.set();
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)
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();
116 if (_container_blocked)
119 _block_count = _child_count;
123 if (_child_block_offsets.back() % _child_block_merge_offsets.back() != 0)
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()
132 _block_count = _child_block_offsets.back() / _child_block_merge_offsets.back();
135 _block_count = _child_block_offsets.back();
136 _size = _child_size_offsets.back();
139 template<
typename Node>
140 OrderingBase(Node& node,
141 bool container_blocked,
143 VirtualOrderingBase<DI,CI>* delegate =
nullptr)
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)
154 , _delegate(delegate)
155 , _gfs_data(gfs_data)
160 template<
typename Node>
161 OrderingBase(Node& node,
162 bool container_blocked,
163 const std::vector<std::size_t>& merge_offsets,
165 VirtualOrderingBase<DI,CI>* delegate =
nullptr)
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)
177 , _delegate(delegate)
178 , _gfs_data(gfs_data)
184 bool containerBlocked()
const
186 return _container_blocked;
189 std::size_t childOrderingCount()
const
194 OrderingBase& childOrdering(
typename Traits::SizeType i)
196 return *_children[i];
199 const OrderingBase& childOrdering(
typename Traits::SizeType i)
const
201 return *_children[i];
204 bool contains(
typename Traits::SizeType codim)
const
206 return _codim_used.test(codim);
214 bool fixedSize(
typename Traits::SizeType codim)
const
216 return _codim_fixed_size.test(codim);
227 void setDelegate(
const VirtualOrderingBase<DI,CI>* delegate)
229 _delegate = delegate;
232 void _mapIndex(
typename Traits::DOFIndexView di,
typename Traits::ContainerIndex& ci)
const
234 typedef typename Traits::SizeType size_type;
235 size_type child_index = di.treeIndex().back();
236 _children[child_index]->mapIndex(di.back_popped(),ci);
239 if (_container_blocked)
240 ci.push_back(child_index);
242 ci.back() += blockOffset(child_index);
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)
253 ci.back() = block_offset;
254 ci.push_back(block_index);
258 size_type block_size = _child_block_merge_offsets.back();
259 ci.back() = block_index * block_size + block_offset;
266 bool update_gfs_data_size(
typename Traits::SizeType&
size,
typename Traits::SizeType& block_count)
const
269 block_count = _block_count;
276 const bool _container_blocked;
279 const std::size_t _child_count;
280 std::vector<OrderingBase*> _children;
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;
288 std::size_t _max_local_size;
290 std::size_t _block_count;
292 const VirtualOrderingBase<DI,CI>* _delegate;
PDELab-specific exceptions.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
std::size_t degree(const Node &node)
Returns the degree of node as run time information.
Definition: nodeinterface.hh:79
void applyToTree(Tree &&tree, Visitor &&visitor)
Apply visitor to TypeTree.
Definition: traversal.hh:239
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
constexpr std::bool_constant<((II==value)||...)> contains(std::integer_sequence< T, II... >, std::integral_constant< T, value >)
Checks whether or not a given sequence contains a value.
Definition: integersequence.hh:137
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:264