3#ifndef DUNE_PDELAB_ORDERING_INTERLEAVEDORDERING_HH
4#define DUNE_PDELAB_ORDERING_INTERLEAVEDORDERING_HH
9#include <dune/typetree/compositenode.hh>
10#include <dune/typetree/powernode.hh>
14#include <dune/pdelab/gridfunctionspace/tags.hh>
15#include <dune/pdelab/ordering/utility.hh>
16#include <dune/pdelab/ordering/orderingbase.hh>
24 namespace interleaved_ordering {
27 template<
typename DI,
typename CI,
typename Node>
29 :
public OrderingBase<DI,CI>
32 typedef OrderingBase<DI,CI> BaseT;
36 typedef typename OrderingBase<DI,CI>::Traits Traits;
40 static const bool consume_tree_index =
true;
48 Base(Node& node,
bool container_blocked,
const OrderingTag& ordering_tag,
typename BaseT::GFSData* gfs_data)
49 : BaseT(node,container_blocked,ordering_tag.offsets(),gfs_data,nullptr)
54 if (node.CHILDREN != ordering_tag.
offsets().size() - 1)
56 "Invalid block structure for InterleavedOrdering: "
57 << node.CHILDREN <<
" children, but "
58 << (ordering_tag.
offsets().size() - 1) <<
" block sizes.");
61 template<
typename ItIn,
typename ItOut>
62 void map_lfs_indices(
const ItIn begin,
const ItIn end, ItOut out)
const
64 typedef typename Traits::SizeType size_type;
65 if (this->_container_blocked)
67 for (ItIn in = begin; in != end; ++in, ++out)
69 size_type child_index = in->treeIndex().back();
70 size_type child_block_offset = this->_child_block_merge_offsets[child_index];
71 size_type child_block_size = this->_child_block_merge_offsets[child_index + 1] - child_block_offset;
72 size_type index = out->back();
73 size_type block_index = index / child_block_size;
74 size_type offset = index % child_block_size;
75 out->back() = child_block_offset + offset;
76 out->push_back(block_index);
81 for (ItIn in = begin; in != end; ++in, ++out)
83 size_type child_index = in->treeIndex().back();
84 size_type child_block_offset = this->_child_block_merge_offsets[child_index];
85 size_type child_block_size = this->_child_block_merge_offsets[child_index + 1] - child_block_offset;
86 size_type block_size = this->_child_block_merge_offsets.back();
87 size_type index = out->back();
88 size_type block_index = index / child_block_size;
89 size_type offset = index % child_block_size;
90 out->back() = block_index * block_size + child_block_offset + offset;
95 template<
typename CIOutIterator,
typename DIOutIterator = DummyDOFIndexIterator>
96 typename Traits::SizeType
97 extract_entity_indices(
const typename Traits::DOFIndex::EntityIndex& ei,
98 typename Traits::SizeType child_index,
99 CIOutIterator ci_out,
const CIOutIterator ci_end)
const
101 typedef typename Traits::SizeType size_type;
102 if (this->_container_blocked)
104 for (; ci_out != ci_end; ++ci_out)
106 size_type child_block_offset = this->_child_block_merge_offsets[child_index];
107 size_type child_block_size = this->_child_block_merge_offsets[child_index + 1] - child_block_offset;
108 size_type index = ci_out->back();
109 size_type block_index = index / child_block_size;
110 size_type offset =index % child_block_size;
111 ci_out->back() = child_block_offset + offset;
112 ci_out->push_back(block_index);
117 for (; ci_out != ci_end; ++ci_out)
119 size_type child_block_offset = this->_child_block_merge_offsets[child_index];
120 size_type child_block_size = this->_child_block_merge_offsets[child_index + 1] - child_block_offset;
121 size_type block_size = this->_child_block_merge_offsets.back();
122 size_type index = ci_out->back();
123 size_type block_index = index / child_block_size;
124 size_type offset =index % child_block_size;
125 ci_out->back() = block_index * block_size + child_block_offset + offset;
138 template<
typename DI,
typename CI,
typename Child, std::
size_t k>
139 class PowerInterleavedOrdering
140 :
public TypeTree::PowerNode<Child, k>
141 ,
public interleaved_ordering::Base<DI,
143 PowerInterleavedOrdering<DI,CI,Child,k>
146 typedef TypeTree::PowerNode<Child, k> Node;
148 typedef interleaved_ordering::Base<DI,
150 PowerInterleavedOrdering<DI,CI,Child,k>
164 PowerInterleavedOrdering(
bool container_blocked,
const InterleavedOrderingTag& ordering_tag,
const typename Node::NodeStorage& children,
typename Base::GFSData* gfs_data)
166 , Base(*this,container_blocked,ordering_tag,gfs_data)
171 for (std::size_t i = 0; i < k; ++i)
173 this->
child(i).update();
178 std::string name()
const {
return "PowerInterleavedOrdering"; }
182 template<
typename GFS,
typename Transformation>
183 struct power_gfs_to_interleaved_ordering_descriptor
186 static const bool recursive =
true;
188 template<
typename TC>
192 typedef PowerInterleavedOrdering<
193 typename Transformation::DOFIndex,
194 typename Transformation::ContainerIndex,
196 TypeTree::StaticDegree<GFS>::value
199 typedef std::shared_ptr<type> storage_type;
203 template<
typename TC>
204 static typename result<TC>::type transform(
const GFS& gfs,
const Transformation& t,
const std::array<std::shared_ptr<TC>,TypeTree::StaticDegree<GFS>::value>& children)
206 return typename result<TC>::type(gfs.backend().blocked(gfs),gfs.orderingTag(),children,
const_cast<GFS*
>(&gfs));
209 template<
typename TC>
210 static typename result<TC>::storage_type transform_storage(std::shared_ptr<const GFS> gfs,
const Transformation& t,
const std::array<std::shared_ptr<TC>,TypeTree::StaticDegree<GFS>::value>& children)
212 return std::make_shared<typename result<TC>::type>(gfs->backend().blocked(*gfs),gfs->orderingTag(),children,
const_cast<GFS*
>(gfs.get()));
217 template<
typename GFS,
typename Transformation>
218 power_gfs_to_interleaved_ordering_descriptor<GFS,Transformation>
219 register_power_gfs_to_ordering_descriptor(GFS*,Transformation*,InterleavedOrderingTag*);
223 template<
typename DI,
typename CI,
typename... Children>
224 class CompositeInterleavedOrdering :
225 public TypeTree::CompositeNode<Children...>,
226 public interleaved_ordering::Base<DI,
228 CompositeInterleavedOrdering<
235 typedef TypeTree::CompositeNode<Children...> Node;
237 typedef interleaved_ordering::Base<
240 CompositeInterleavedOrdering<
257 CompositeInterleavedOrdering(
bool backend_blocked,
const InterleavedOrderingTag& ordering_tag,
typename Base::GFSData* gfs_data, std::shared_ptr<Children>... children)
259 , Base(*this,backend_blocked,ordering_tag,gfs_data)
262 std::string name()
const {
return "CompositeInterleavedOrdering"; }
271 template<
typename GFS,
typename Transformation>
272 struct composite_gfs_to_interleaved_ordering_descriptor
275 static const bool recursive =
true;
277 template<
typename... TC>
281 typedef CompositeInterleavedOrdering<
282 typename Transformation::DOFIndex,
283 typename Transformation::ContainerIndex,
287 typedef std::shared_ptr<type> storage_type;
291 template<
typename... TC>
292 static typename result<TC...>::type transform(
const GFS& gfs,
const Transformation& t, std::shared_ptr<TC>... children)
294 return typename result<TC...>::type(gfs.backend().blocked(gfs),gfs.orderingTag(),
const_cast<GFS*
>(&gfs),children...);
297 template<
typename... TC>
298 static typename result<TC...>::storage_type transform_storage(std::shared_ptr<const GFS> gfs,
const Transformation& t, std::shared_ptr<TC>... children)
300 return std::make_shared<
typename result<TC...>::type>(gfs->backend().blocked(*gfs),gfs.orderingTag(),
const_cast<GFS*
>(gfs.get()),children...);
305 template<
typename GFS,
typename Transformation>
306 composite_gfs_to_interleaved_ordering_descriptor<GFS,Transformation>
307 register_composite_gfs_to_ordering_descriptor(GFS*,Transformation*,InterleavedOrderingTag*);
Error related to the logical structure of an Ordering.
Definition: exceptions.hh:46
Interface for merging index spaces.
Definition: interleavedordering.hh:30
Base(Node &node, bool container_blocked, const OrderingTag &ordering_tag, typename BaseT::GFSData *gfs_data)
Construct ordering object.
Definition: interleavedordering.hh:48
std::array< ChildStorageType, k > NodeStorage
The type used for storing the children.
Definition: powernode.hh:84
PDELab-specific exceptions.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void applyToTree(Tree &&tree, Visitor &&visitor)
Apply visitor to TypeTree.
Definition: traversal.hh:213
Dune namespace.
Definition: alignedallocator.hh:14
Indicate interleaved ordering of the unknowns of non-leaf grid function spaces according to a given b...
Definition: tags.hh:79
const std::vector< std::size_t > & offsets() const
Returns a list of offsets for the child blocks.
Definition: tags.hh:115