AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
LocalView.hpp
1 #pragma once
2 
3 #include <tuple>
4 #include <optional>
5 #include <vector>
6 
7 #include <dune/common/concept.hh>
8 #include <dune/functions/functionspacebases/concepts.hh>
9 
10 #include <amdis/common/OptionalNoCopy.hpp>
11 #include <amdis/common/TypeTraits.hpp>
12 #include <amdis/functions/NodeCache.hpp>
13 #include <amdis/functions/Nodes.hpp>
14 
15 // NOTE: This is a variant of dune-functions DefaultLocalView
16 
17 namespace AMDiS
18 {
20  template <class GB>
21  class LocalView
22  {
23  using PrefixPath = Dune::TypeTree::HybridTreePath<>;
24 
25  // Node index set provided by PreBasis
26  using NodeIndexSet = NodeIndexSet_t<typename GB::PreBasis, PrefixPath>;
27 
28  public:
30  using GlobalBasis = GB;
31 
33  using GridView = typename GlobalBasis::GridView;
34 
36  using Element = typename GridView::template Codim<0>::Entity;
37 
39  using size_type = std::size_t;
40 
42  using Tree = Node_t<typename GlobalBasis::PreBasis, PrefixPath>;
43 
45  using TreeCache = NodeCache_t<Tree>;
46 
49 
50  public:
53  : globalBasis_(&globalBasis)
54  , tree_(makeNode(globalBasis_->preBasis(), PrefixPath{}))
55  , nodeIndexSet_(makeNodeIndexSet(globalBasis_->preBasis(), PrefixPath{}))
56  {
57  static_assert(Dune::models<Dune::Functions::Concept::BasisTree<GridView>, Tree>());
58  initializeTree(tree_);
59  }
60 
62 
67  void bind (Element const& element)
68  {
69  element_ = element;
70  bindTree(tree_, element_);
71  indices_.resize(tree_.size());
72  bindNodeIndexSet(globalBasis_->preBasis(), nodeIndexSet_, tree_, indices_.begin(), Dune::PriorityTag<5>{});
73  bound_ = true;
74  }
75 
76  // do not bind an rvalue-reference
77  void bind (Element&& element)
78  {
79  static_assert(not std::is_rvalue_reference_v<Element&&>);
80  }
81 
83  bool isBound () const
84  {
85  return bound_;
86  }
87 
89  Element const& element () const
90  {
91  assert(isBound());
92  return element_;
93  }
94 
96 
99  void unbind ()
100  {
101  unbindNodeIndexSet(nodeIndexSet_, Dune::PriorityTag<5>{});
102  bound_ = false;
103  }
104 
106  Tree const& tree () const
107  {
108  return tree_;
109  }
110 
112  TreeCache const& treeCache () const
113  {
114  if (!treeCache_)
115  treeCache_.emplace(makeNodeCache(tree_));
116 
117  return *treeCache_;
118  }
119 
121  size_type size () const
122  {
123  assert(isBound());
124  return tree_.size();
125  }
126 
128 
132  {
133  return globalBasis_->preBasis().maxNodeSize();
134  }
135 
139  {
140  assert(isBound());
141  return indices_[i];
142  }
143 
145  GlobalBasis const& globalBasis () const
146  {
147  return *globalBasis_;
148  }
149 
151  LocalView const& rootLocalView () const
152  {
153  return *this;
154  }
155 
156  private:
157  template <class PB, class NIS, class Tree, class Iter,
158  class = decltype(std::declval<PB>().indices(std::declval<Tree>(), std::declval<Iter>()))>
159  static void bindNodeIndexSet(PB const& preBasis, NIS& /*nodeIndexSet*/, Tree const& tree, Iter it, Dune::PriorityTag<3>)
160  {
161  preBasis.indices(tree, it);
162  }
163 
164  template <class PB, class NIS, class Tree, class Iter,
165  class = decltype(std::declval<NIS>().bind(std::declval<Tree>())),
166  class = decltype(std::declval<NIS>().indices(std::declval<Iter>()))>
167  static void bindNodeIndexSet(PB const& /*preBasis*/, NIS& nodeIndexSet, Tree const& tree, Iter it, Dune::PriorityTag<2>)
168  {
169  nodeIndexSet.bind(tree);
170  nodeIndexSet.indices(it);
171  }
172 
173  template <class PB, class NIS, class Tree, class Iter,
174  class = decltype(std::declval<NIS>().bind(std::declval<Tree>())),
175  class = decltype(std::declval<NIS>().index(0u))>
176  static void bindNodeIndexSet(PB const& /*preBasis*/, NIS& nodeIndexSet, Tree const& tree, Iter it, Dune::PriorityTag<1>)
177  {
178  nodeIndexSet.bind(tree);
179  for (size_type i = 0; i < tree.size(); ++i)
180  *it++ = nodeIndexSet.index(i);
181  }
182 
183  template <class PB, class NIS, class Tree, class Iter>
184  static void bindNodeIndexSet(PB const& /*preBasis*/, NIS& /*nodeIndexSet*/, Tree const& /*tree*/, Iter /*it*/, Dune::PriorityTag<0>)
185  {
186  /* do nothing */
187  }
188 
189  template <class NIS,
190  class = decltype(std::declval<NIS>().unbind())>
191  static void unbindNodeIndexSet(NIS& nodeIndexSet, Dune::PriorityTag<1>)
192  {
193  nodeIndexSet.unbind();
194  }
195 
196  template <class NIS>
197  static void unbindNodeIndexSet(NIS& /*nodeIndexSet*/, Dune::PriorityTag<0>)
198  {
199  /* do nothing */
200  }
201 
202  protected:
203  GlobalBasis const* globalBasis_;
204  Element element_;
205  Tree tree_;
206  NodeIndexSet nodeIndexSet_;
207 
208  mutable OptionalNoCopy<TreeCache> treeCache_ = std::nullopt;
209  std::vector<MultiIndex> indices_ = {};
210  bool bound_ = false;
211  };
212 
213 } // end namespace AMDiS
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:121
constexpr bool MultiIndex
A multi-index type.
Definition: Concepts.hpp:150
size_type maxSize() const
Maximum local size for any element on the GridView.
Definition: LocalView.hpp:131
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: LocalView.hpp:33
LocalView const & rootLocalView() const
Return this local-view.
Definition: LocalView.hpp:151
typename GlobalBasis::PreBasis::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: LocalView.hpp:48
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:83
Definition: AdaptBase.hpp:6
MultiIndex index(size_type i) const
Maps from subtree index set [0..size-1] to a globally unique multi index in global basis...
Definition: LocalView.hpp:138
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: LocalView.hpp:36
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:21
Node_t< typename GlobalBasis::PreBasis, PrefixPath > Tree
Tree of local finite elements / local shape function sets.
Definition: LocalView.hpp:42
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:112
Tree const & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: LocalView.hpp:106
GlobalBasis const & globalBasis() const
Return the global basis that we are a view on.
Definition: LocalView.hpp:145
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:89
void bind(Element const &element)
Bind the view to a grid element.
Definition: LocalView.hpp:67
typename PreBasis::GridView GridView
The grid view that the FE space is defined on.
Definition: GlobalBasis.hpp:62
void unbind()
Unbind from the current element.
Definition: LocalView.hpp:99
LocalView(GlobalBasis const &globalBasis)
Construct local view for a given global finite element basis.
Definition: LocalView.hpp:52
GB GlobalBasis
The global FE basis that this is a view on.
Definition: LocalView.hpp:30
NodeCache_t< Tree > TreeCache
Cached basis-tree.
Definition: LocalView.hpp:45
std::size_t size_type
The type used for sizes.
Definition: LocalView.hpp:39