dune-pdelab  2.4.1
datahandleprovider.hh
Go to the documentation of this file.
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_DATAHANDLEPROVIDER_HH
4 #define DUNE_PDELAB_DATAHANDLEPROVIDER_HH
5 
6 #include <vector>
7 #include <stack>
8 
9 #include <dune/common/typetraits.hh>
10 #include <dune/common/reservedvector.hh>
11 #include <dune/common/std/constexpr.hh>
12 #include <dune/typetree/visitor.hh>
13 
15 
16 namespace Dune {
17  namespace PDELab {
18 
19  namespace {
20 
21  template<typename EntityIndex>
22  struct get_size_for_entity
23  : public TypeTree::TreeVisitor
24  , public TypeTree::DynamicTraversal
25  {
26 
27  template<typename Ordering, typename TreePath>
28  void leaf(const Ordering& ordering, TreePath tp)
29  {
30  _size += ordering.size(_entity_index);
31  }
32 
33  get_size_for_entity(const EntityIndex& entity_index)
34  : _size(0)
35  , _entity_index(entity_index)
36  {}
37 
38  std::size_t size() const
39  {
40  return _size;
41  }
42 
43  private:
44 
45  std::size_t _size;
46  const EntityIndex& _entity_index;
47 
48  };
49 
50 
51  template<typename EntityIndex, typename OffsetIterator>
52  struct get_leaf_offsets_for_entity
53  : public TypeTree::TreeVisitor
54  , public TypeTree::DynamicTraversal
55  {
56 
57  template<typename Ordering, typename TreePath>
58  void leaf(const Ordering& ordering, TreePath tp)
59  {
60  *(++_oit) = ordering.size(_entity_index);
61  }
62 
63  get_leaf_offsets_for_entity(const EntityIndex& entity_index, OffsetIterator oit)
64  : _oit(oit)
65  , _entity_index(entity_index)
66  {}
67 
69  OffsetIterator offsetIterator() const
70  {
71  return _oit;
72  }
73 
74  private:
75 
76  OffsetIterator _oit;
77  const EntityIndex& _entity_index;
78 
79  };
80 
81 
82  template<typename DOFIndex, typename ContainerIndex, std::size_t tree_depth, bool map_dof_indices = false>
83  struct indices_for_entity
84  : public TypeTree::TreeVisitor
85  , public TypeTree::DynamicTraversal
86  {
87 
88  typedef std::size_t size_type;
89  typedef typename DOFIndex::EntityIndex EntityIndex;
90  typedef typename std::vector<ContainerIndex>::iterator CIIterator;
91  typedef typename std::conditional<
92  map_dof_indices,
93  typename std::vector<DOFIndex>::iterator,
94  DummyDOFIndexIterator
95  >::type DIIterator;
96 
97 
98  template<typename Ordering, typename Child, typename TreePath, typename ChildIndex>
99  void beforeChild(const Ordering& ordering, const Child& child, TreePath tp, ChildIndex childIndex)
100  {
101  _stack.push(std::make_pair(_ci_it,_di_it));
102  }
103 
104  template<typename Ordering, typename TreePath>
105  void leaf(const Ordering& ordering, TreePath tp)
106  {
107  size_type size = ordering.extract_entity_indices(_entity_index,
108  tp.back(),
109  _ci_it,
110  _ci_end,
111  _di_it);
112 
113  _ci_end += size;
114  _ci_it = _ci_end;
115  _di_end += size;
116  _di_it = _di_end;
117  }
118 
119  template<typename Ordering, typename Child, typename TreePath, typename ChildIndex>
120  void afterChild(const Ordering& ordering, const Child& child, TreePath tp, ChildIndex childIndex)
121  {
122  // pop
123  ordering.extract_entity_indices(_entity_index,
124  childIndex,
125  _stack.top().first,
126  _ci_end);
127 
128  if (Ordering::consume_tree_index)
129  for (DIIterator it = _stack.top().second;
130  it != _di_end;
131  ++it)
132  it->treeIndex().push_back(childIndex);
133 
134  _stack.pop();
135  }
136 
137 
138  indices_for_entity(const EntityIndex& entity_index,
139  CIIterator ci_begin,
140  DIIterator di_begin = DIIterator())
141  : _entity_index(entity_index)
142  , _ci_it(ci_begin)
143  , _ci_end(ci_begin)
144  , _di_it(di_begin)
145  , _di_end(di_begin)
146  {}
147 
148 
149  // Exposed for multidomain support
150  CIIterator ci_end() const
151  {
152  return _ci_end;
153  }
154 
155  // Exposed for multidomain support
156  DIIterator di_end() const
157  {
158  return _di_end;
159  }
160 
161  private:
162 
163  const EntityIndex& _entity_index;
164  CIIterator _ci_it;
165  CIIterator _ci_end;
166  DIIterator _di_it;
167  DIIterator _di_end;
168 
169  std::stack<
170  std::pair<
171  CIIterator,
172  DIIterator
173  >,
174  ReservedVector<
175  std::pair<
176  CIIterator,
177  DIIterator
178  >,
179  tree_depth
180  >
181  > _stack;
182  };
183 
184  } // anonymous namespace
185 
186 
187  template<typename GFS>
189  {
190 
191  public:
192 
193  typedef std::size_t size_type;
194 
195  //------------------------------
196  // generic data handle interface
197  //------------------------------
198 
200  bool dataHandleContains (int codim) const
201  {
202  return gfs().ordering().contains(codim);
203  }
204 
206  bool dataHandleFixedSize (int codim) const
207  {
208  return gfs().ordering().fixedSize(codim);
209  }
210 
212 
226  DUNE_CONSTEXPR bool sendLeafSizes() const
227  {
228  return false;
229  }
230 
235  template<typename Entity>
236  size_type dataHandleSize (const Entity& e) const
237  {
238  typedef typename GFS::Ordering Ordering;
239 
240  typedef typename Ordering::Traits::DOFIndex::EntityIndex EntityIndex;
241  EntityIndex ei;
242 
243  Ordering::Traits::DOFIndexAccessor::GeometryIndex::store(
244  ei,
245  e.type(),
246  gfs().gridView().indexSet().index(e)
247  );
248 
249  get_size_for_entity<EntityIndex> get_size(ei);
250  TypeTree::applyToTree(gfs().ordering(),get_size);
251 
252  return get_size.size();
253  }
254 
255  template<typename V, typename EntityIndex>
256  void setup_dof_indices(V& v, size_type n, const EntityIndex& ei, std::integral_constant<bool,true>) const
257  {
258  v.resize(n);
259  for (typename V::iterator it = v.begin(),
260  endit = v.end();
261  it != endit;
262  ++it)
263  {
264  it->treeIndex().clear();
265  it->entityIndex() = ei;
266  }
267  }
268 
269  template<typename V, typename EntityIndex>
270  void setup_dof_indices(V& v, size_type n, const EntityIndex& ei, std::integral_constant<bool,false>) const
271  {}
272 
273  template<typename V>
274  typename V::iterator dof_indices_begin(V& v, std::integral_constant<bool,true>) const
275  {
276  return v.begin();
277  }
278 
279  template<typename V>
280  DummyDOFIndexIterator dof_indices_begin(V& v, std::integral_constant<bool,false>) const
281  {
282  return DummyDOFIndexIterator();
283  }
284 
286  template<typename Entity, typename ContainerIndex, typename DOFIndex, typename OffsetIterator, bool map_dof_indices>
287  void dataHandleIndices (const Entity& e,
288  std::vector<ContainerIndex>& container_indices,
289  std::vector<DOFIndex>& dof_indices,
290  OffsetIterator oit,
291  std::integral_constant<bool,map_dof_indices> map_dof_indices_value
292  ) const
293  {
294  typedef typename GFS::Ordering Ordering;
295 
297  "dataHandleContainerIndices() called with invalid ContainerIndex type.");
298 
299  typedef typename Ordering::Traits::DOFIndex::EntityIndex EntityIndex;
300  EntityIndex ei;
301 
302  Ordering::Traits::DOFIndexAccessor::GeometryIndex::store(
303  ei,
304  e.type(),
305  gfs().entitySet().indexSet().index(e)
306  );
307 
308  get_leaf_offsets_for_entity<EntityIndex,OffsetIterator> get_offsets(ei,oit);
309  TypeTree::applyToTree(gfs().ordering(),get_offsets);
310  OffsetIterator end_oit = oit + (TypeTree::TreeInfo<Ordering>::leafCount + 1);
311 
312  // convert sizes to offsets - last entry contains total size
313  std::partial_sum(oit,end_oit,oit);
314  size_type size = *(oit + TypeTree::TreeInfo<Ordering>::leafCount);
315 
316  container_indices.resize(size);
317  // Clear index state
318  for (typename std::vector<ContainerIndex>::iterator it = container_indices.begin(),
319  endit = container_indices.end();
320  it != endit;
321  ++it)
322  it->clear();
323 
324  setup_dof_indices(dof_indices,size,ei,map_dof_indices_value);
325 
326  indices_for_entity<
327  DOFIndex,
328  ContainerIndex,
329  TypeTree::TreeInfo<Ordering>::depth,
330  map_dof_indices
331  > extract_indices(ei,container_indices.begin(),dof_indices_begin(dof_indices,map_dof_indices_value));
332  TypeTree::applyToTree(gfs().ordering(),extract_indices);
333 
334  }
335 
336  protected:
337 
338  const GFS& gfs() const
339  {
340  return static_cast<const GFS&>(*this);
341  }
342 
343  };
344 
345  } // namespace PDELab
346 } // namespace Dune
347 
348 #endif // DUNE_PDELAB_DATAHANDLEPROVIDER
void setup_dof_indices(V &v, size_type n, const EntityIndex &ei, std::integral_constant< bool, false >) const
Definition: datahandleprovider.hh:270
bool dataHandleContains(int codim) const
returns true if data for this codim should be communicated
Definition: datahandleprovider.hh:200
std::array< T, entity_capacity > EntityIndex
Definition: dofindex.hh:157
static const unsigned int value
Definition: gridfunctionspace/tags.hh:177
Dummy iterator type over DOF indices.
Definition: ordering/utility.hh:288
size_type dataHandleSize(const Entity &e) const
Definition: datahandleprovider.hh:236
Definition: adaptivity.hh:27
const GFS & gfs() const
Definition: datahandleprovider.hh:338
Definition: datahandleprovider.hh:188
const Entity & e
Definition: localfunctionspace.hh:111
V::iterator dof_indices_begin(V &v, std::integral_constant< bool, true >) const
Definition: datahandleprovider.hh:274
std::size_t size_type
Definition: datahandleprovider.hh:193
DummyDOFIndexIterator dof_indices_begin(V &v, std::integral_constant< bool, false >) const
Definition: datahandleprovider.hh:280
DUNE_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:226
A multi-index representing a degree of freedom in a GridFunctionSpace.
Definition: dofindex.hh:146
void setup_dof_indices(V &v, size_type n, const EntityIndex &ei, std::integral_constant< bool, true >) const
Definition: datahandleprovider.hh:256
bool dataHandleFixedSize(int codim) const
returns true if size per entity of given dim and codim is a constant
Definition: datahandleprovider.hh:206
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:287