dune-pdelab  2.4-dev
orderinginterface.hh
Go to the documentation of this file.
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_ORDERINGINTERFACE_HH
5 #define DUNE_PDELAB_ORDERING_ORDERINGINTERFACE_HH
6 
7 #include <cstddef>
8 
9 #include <dune/common/documentation.hh>
10 
11 #include <dune/typetree/nodeinterface.hh>
12 
13 namespace Dune {
14  namespace PDELab {
15 
18 
19  //===============================================================
20  // Utilities for the power and composite gfs
21  // ===============================================================
22 
25  public TypeTree::NodeInterface
26  {
27  typedef ImplementationDefined SizeType;
28 
30 
39  (const TypeTree::NodeInterface::NodeStorage &children);
40 
42 
51 
53 
59  void update();
60 
62  bool blocked() const;
63 
66 
71  bool fixedSize() const;
72 
74 
82  SizeType subMap(SizeType child, SizeType indexInChild) const;
83 
85  SizeType size() const;
86 
89 
94  SizeType maxLocalSize() const;
95 
98 
110  template<class Entity>
111  SizeType entitySize(const Entity &e) const;
113 
119  template<class Element>
120  SizeType entitySize(const Element &e, std::size_t codim,
121  std::size_t subentity) const;
123  template<class Intersection>
124  SizeType intersectionSize(const Intersection &i) const;
125 
128 
141  template<class Entity>
142  SizeType entityOffset(const Entity &e) const;
145 
152  template<class Element>
153  SizeType entityOffset(const Element &e, std::size_t codim,
154  std::size_t subentity) const;
156  template<class Intersection>
157  SizeType intersectionOffset(const Intersection &i) const;
158  };
159 
161  } // namespace PDELab
162 } // namespace Dune
163 
164 #endif // DUNE_PDELAB_ORDERING_ORDERINGINTERFACE_HH
GridFunctionSpaceOrderingInterface()
Construct ordering object.
SizeType entitySize(const Entity &e) const
number of indices attached to a given entity (of arbitrary codimension)
SizeType maxLocalSize() const
maximum number of dofs attached to any given element and all of its subentities and intersections ...
const E & e
Definition: interpolate.hh:172
Interface for merging index spaces.
Definition: orderinginterface.hh:24
SizeType size() const
number of indices in this ordering
bool blocked() const
whether dofs are blocked per entity/intersection
SizeType subMap(SizeType child, SizeType indexInChild) const
map a global dof index from a child
SizeType entityOffset(const Entity &e) const
offset of the block of dofs attached to a given entity (of arbitrary codimension) ...
void update()
update internal data structures
Definition: adaptivity.hh:27
SizeType intersectionOffset(const Intersection &i) const
offset of the block of dofs attached to a given intersection
bool fixedSize() const
whether all entites of the same geometry type/all intersections have the same number of dofs ...
ImplementationDefined SizeType
Definition: orderinginterface.hh:27
SizeType intersectionSize(const Intersection &i) const
number of indices attached to a given intersection