3 #ifndef DUNE_PDELAB_GRIDFUNCTIONSPACE_DATAHANDLEPROVIDER_HH 4 #define DUNE_PDELAB_GRIDFUNCTIONSPACE_DATAHANDLEPROVIDER_HH 9 #include <dune/common/typetraits.hh> 10 #include <dune/common/reservedvector.hh> 11 #include <dune/typetree/visitor.hh> 20 template<
typename EntityIndex>
21 struct get_size_for_entity
22 :
public TypeTree::TreeVisitor
23 ,
public TypeTree::DynamicTraversal
26 template<
typename Ordering,
typename TreePath>
27 void leaf(
const Ordering& ordering, TreePath tp)
29 _size += ordering.size(_entity_index);
32 get_size_for_entity(
const EntityIndex& entity_index)
34 , _entity_index(entity_index)
37 std::size_t size()
const 45 const EntityIndex& _entity_index;
50 template<
typename EntityIndex,
typename OffsetIterator>
51 struct get_leaf_offsets_for_entity
52 :
public TypeTree::TreeVisitor
53 ,
public TypeTree::DynamicTraversal
56 template<
typename Ordering,
typename TreePath>
57 void leaf(
const Ordering& ordering, TreePath tp)
59 *(++_oit) = ordering.size(_entity_index);
62 get_leaf_offsets_for_entity(
const EntityIndex& entity_index, OffsetIterator oit)
64 , _entity_index(entity_index)
68 OffsetIterator offsetIterator()
const 76 const EntityIndex& _entity_index;
81 template<
typename DOFIndex,
typename ContainerIndex, std::
size_t tree_depth,
bool map_dof_indices = false>
82 struct indices_for_entity
83 :
public TypeTree::TreeVisitor
84 ,
public TypeTree::DynamicTraversal
87 typedef std::size_t size_type;
89 typedef typename std::vector<ContainerIndex>::iterator CIIterator;
90 typedef typename std::conditional<
92 typename std::vector<DOFIndex>::iterator,
97 template<
typename Ordering,
typename Child,
typename TreePath,
typename ChildIndex>
98 void beforeChild(
const Ordering& ordering,
const Child& child, TreePath tp, ChildIndex childIndex)
100 _stack.push(std::make_pair(_ci_it,_di_it));
103 template<
typename Ordering,
typename TreePath>
104 void leaf(
const Ordering& ordering, TreePath tp)
106 size_type size = ordering.extract_entity_indices(_entity_index,
118 template<
typename Ordering,
typename Child,
typename TreePath,
typename ChildIndex>
119 void afterChild(
const Ordering& ordering,
const Child& child, TreePath tp, ChildIndex childIndex)
122 ordering.extract_entity_indices(_entity_index,
127 if (Ordering::consume_tree_index)
128 for (DIIterator it = _stack.top().second;
131 it->treeIndex().push_back(childIndex);
137 indices_for_entity(
const EntityIndex& entity_index,
139 DIIterator di_begin = DIIterator())
140 : _entity_index(entity_index)
149 CIIterator ci_end()
const 155 DIIterator di_end()
const 162 const EntityIndex& _entity_index;
186 template<
typename GFS>
201 return gfs().ordering().contains(codim);
207 return gfs().ordering().fixedSize(codim);
234 template<
typename Entity>
237 typedef typename GFS::Ordering Ordering;
239 typedef typename Ordering::Traits::DOFIndex::EntityIndex EntityIndex;
242 Ordering::Traits::DOFIndexAccessor::GeometryIndex::store(
245 gfs().gridView().indexSet().index(e)
248 get_size_for_entity<EntityIndex> get_size(ei);
249 TypeTree::applyToTree(gfs().ordering(),get_size);
251 return get_size.size();
254 template<
typename V,
typename EntityIndex>
255 void setup_dof_indices(V& v, size_type n,
const EntityIndex& ei, std::integral_constant<bool,true>)
const 258 for (
typename V::iterator it = v.begin(),
263 it->treeIndex().clear();
264 it->entityIndex() = ei;
268 template<
typename V,
typename EntityIndex>
269 void setup_dof_indices(V& v, size_type n,
const EntityIndex& ei, std::integral_constant<bool,false>)
const 285 template<
typename Entity,
typename ContainerIndex,
typename DOFIndex,
typename OffsetIterator,
bool map_dof_indices>
287 std::vector<ContainerIndex>& container_indices,
288 std::vector<DOFIndex>& dof_indices,
290 std::integral_constant<bool,map_dof_indices> map_dof_indices_value
293 typedef typename GFS::Ordering Ordering;
296 "dataHandleContainerIndices() called with invalid ContainerIndex type.");
298 typedef typename Ordering::Traits::DOFIndex::EntityIndex EntityIndex;
301 Ordering::Traits::DOFIndexAccessor::GeometryIndex::store(
304 gfs().entitySet().indexSet().index(e)
307 get_leaf_offsets_for_entity<EntityIndex,OffsetIterator> get_offsets(ei,oit);
308 TypeTree::applyToTree(gfs().ordering(),get_offsets);
309 OffsetIterator end_oit = oit + (TypeTree::TreeInfo<Ordering>::leafCount + 1);
312 std::partial_sum(oit,end_oit,oit);
313 size_type size = *(oit + TypeTree::TreeInfo<Ordering>::leafCount);
315 container_indices.resize(size);
317 for (
typename std::vector<ContainerIndex>::iterator it = container_indices.begin(),
318 endit = container_indices.end();
323 setup_dof_indices(dof_indices,size,ei,map_dof_indices_value);
328 TypeTree::TreeInfo<Ordering>::depth,
330 > extract_indices(ei,container_indices.begin(),dof_indices_begin(dof_indices,map_dof_indices_value));
331 TypeTree::applyToTree(gfs().ordering(),extract_indices);
339 return static_cast<const GFS&
>(*this);
347 #endif // DUNE_PDELAB_GRIDFUNCTIONSPACE_DATAHANDLEPROVIDER_HH constexpr bool sendLeafSizes() const
Returns true if the sizes of the leaf orderings in this tree should be sent as part of the communcati...
Definition: datahandleprovider.hh:225
size_type dataHandleSize(const Entity &e) const
Definition: datahandleprovider.hh:235
void dataHandleIndices(const Entity &e, std::vector< ContainerIndex > &container_indices, std::vector< DOFIndex > &dof_indices, OffsetIterator oit, std::integral_constant< bool, map_dof_indices > map_dof_indices_value) const
return vector of global indices associated with the given entity
Definition: datahandleprovider.hh:286
DummyDOFIndexIterator dof_indices_begin(V &v, std::integral_constant< bool, false >) const
Definition: datahandleprovider.hh:279
void setup_dof_indices(V &v, size_type n, const EntityIndex &ei, std::integral_constant< bool, false >) const
Definition: datahandleprovider.hh:269
V::iterator dof_indices_begin(V &v, std::integral_constant< bool, true >) const
Definition: datahandleprovider.hh:273
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
bool dataHandleContains(int codim) const
returns true if data for this codim should be communicated
Definition: datahandleprovider.hh:199
bool dataHandleFixedSize(int codim) const
returns true if size per entity of given dim and codim is a constant
Definition: datahandleprovider.hh:205
const Entity & e
Definition: localfunctionspace.hh:111
std::size_t size_type
Definition: datahandleprovider.hh:192
Dummy iterator type over DOF indices.
Definition: ordering/utility.hh:288
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139
const GFS & gfs() const
Definition: datahandleprovider.hh:337
std::array< T, entity_capacity > EntityIndex
Definition: dofindex.hh:157
Definition: datahandleprovider.hh:187
void setup_dof_indices(V &v, size_type n, const EntityIndex &ei, std::integral_constant< bool, true >) const
Definition: datahandleprovider.hh:255
A multi-index representing a degree of freedom in a GridFunctionSpace.
Definition: dofindex.hh:146