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/TypeTraits.hpp>
11 #include <amdis/functions/NodeCache.hpp>
12 #include <amdis/functions/Nodes.hpp>
13 
14 // NOTE: This is a variant of dune-functions DefaultLocalView
15 
16 namespace AMDiS
17 {
19  template <class GB>
20  class LocalView
21  {
22  using PrefixPath = Dune::TypeTree::HybridTreePath<>;
23 
24  // Node index set provided by PreBasis
25  using NodeIndexSet = NodeIndexSet_t<typename GB::PreBasis, PrefixPath>;
26 
27  public:
29  using GlobalBasis = GB;
30 
32  using GridView = typename GlobalBasis::GridView;
33 
35  using Element = typename GridView::template Codim<0>::Entity;
36 
38  using size_type = std::size_t;
39 
41  using Tree = Node_t<typename GlobalBasis::PreBasis, PrefixPath>;
42 
44  using TreeCache = NodeCache_t<Tree>;
45 
48 
49  public:
52  : globalBasis_(&globalBasis)
53  , tree_(makeNode(globalBasis_->preBasis(), PrefixPath{}))
54  , treeCache_(makeNodeCache(tree_))
55  , nodeIndexSet_(makeNodeIndexSet(globalBasis_->preBasis(), PrefixPath{}))
56  {
57  static_assert(Dune::models<Dune::Functions::Concept::BasisTree<GridView>, Tree>());
58  initializeTree(tree_);
59  }
60 
62  // NOTE: needs to be implemented manually, because of the reference dependency
63  // between members tree_ and treeCache_
64  LocalView (LocalView const& other)
65  : globalBasis_(other.globalBasis_)
66  , element_(other.element_)
67  , tree_(other.tree_)
68  , treeCache_(makeNodeCache(tree_))
69  , nodeIndexSet_(other.nodeIndexSet_)
70  , indices_(other.indices_)
71  {}
72 
73  // Move constructor
74  // NOTE: needs to be implemented manually, because of the reference dependency
75  // between members tree_ and treeCache_
76  LocalView (LocalView&& other)
77  : globalBasis_(other.globalBasis_)
78  , element_(other.element_)
79  , tree_(std::move(other.tree_))
80  , treeCache_(makeNodeCache(tree_))
81  , nodeIndexSet_(std::move(other.nodeIndexSet_))
82  , indices_(std::move(other.indices_))
83  {}
84 
86 
91  void bind (Element const& element)
92  {
93  element_ = &element;
94  bindTree(tree_, element);
95  indices_.resize(tree_.size());
96  bindNodeIndexSet(globalBasis_->preBasis(), nodeIndexSet_, tree_, indices_.begin(), Dune::PriorityTag<5>{});
97  }
98 
99  // do not bind an rvalue-reference
100  void bind (Element&& element)
101  {
102  static_assert(not std::is_rvalue_reference_v<Element&&>);
103  }
104 
106  bool isBound () const
107  {
108  return element_ != nullptr;
109  }
110 
112  Element const& element () const
113  {
114  assert(isBound());
115  return *element_;
116  }
117 
119 
122  void unbind ()
123  {
124  unbindNodeIndexSet(nodeIndexSet_, Dune::PriorityTag<5>{});
125  element_ = nullptr;
126  }
127 
129  Tree const& tree () const
130  {
131  return tree_;
132  }
133 
135  TreeCache const& treeCache () const
136  {
137  return treeCache_;
138  }
139 
141  size_type size () const
142  {
143  assert(isBound());
144  return tree_.size();
145  }
146 
148 
152  {
153  return globalBasis_->preBasis().maxNodeSize();
154  }
155 
159  {
160  assert(isBound());
161  return indices_[i];
162  }
163 
165  GlobalBasis const& globalBasis () const
166  {
167  return *globalBasis_;
168  }
169 
171  LocalView const& rootLocalView () const
172  {
173  return *this;
174  }
175 
176  private:
177  template <class PB, class NIS, class Tree, class Iter,
178  class = decltype(std::declval<PB>().indices(std::declval<Tree>(), std::declval<Iter>()))>
179  static void bindNodeIndexSet(PB const& preBasis, NIS& /*nodeIndexSet*/, Tree const& tree, Iter it, Dune::PriorityTag<3>)
180  {
181  preBasis.indices(tree, it);
182  }
183 
184  template <class PB, class NIS, class Tree, class Iter,
185  class = decltype(std::declval<NIS>().bind(std::declval<Tree>())),
186  class = decltype(std::declval<NIS>().indices(std::declval<Iter>()))>
187  static void bindNodeIndexSet(PB const& /*preBasis*/, NIS& nodeIndexSet, Tree const& tree, Iter it, Dune::PriorityTag<2>)
188  {
189  nodeIndexSet.bind(tree);
190  nodeIndexSet.indices(it);
191  }
192 
193  template <class PB, class NIS, class Tree, class Iter,
194  class = decltype(std::declval<NIS>().bind(std::declval<Tree>())),
195  class = decltype(std::declval<NIS>().index(0u))>
196  static void bindNodeIndexSet(PB const& /*preBasis*/, NIS& nodeIndexSet, Tree const& tree, Iter it, Dune::PriorityTag<1>)
197  {
198  nodeIndexSet.bind(tree);
199  for (size_type i = 0; i < tree.size(); ++i)
200  *it++ = nodeIndexSet.index(i);
201  }
202 
203  template <class PB, class NIS, class Tree, class Iter>
204  static void bindNodeIndexSet(PB const& /*preBasis*/, NIS& /*nodeIndexSet*/, Tree const& /*tree*/, Iter /*it*/, Dune::PriorityTag<0>)
205  {
206  /* do nothing */
207  }
208 
209  template <class NIS,
210  class = decltype(std::declval<NIS>().unbind())>
211  static void unbindNodeIndexSet(NIS& nodeIndexSet, Dune::PriorityTag<1>)
212  {
213  nodeIndexSet.unbind();
214  }
215 
216  template <class NIS>
217  static void unbindNodeIndexSet(NIS& /*nodeIndexSet*/, Dune::PriorityTag<0>)
218  {
219  /* do nothing */
220  }
221 
222  protected:
223  GlobalBasis const* globalBasis_;
224  Element const* element_ = nullptr;
225  Tree tree_;
226  TreeCache treeCache_;
227  NodeIndexSet nodeIndexSet_;
228  std::vector<MultiIndex> indices_;
229  };
230 
231 } // end namespace AMDiS
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:141
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:151
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: LocalView.hpp:32
LocalView const & rootLocalView() const
Return this local-view.
Definition: LocalView.hpp:171
typename GlobalBasis::PreBasis::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: LocalView.hpp:47
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:106
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:158
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: LocalView.hpp:35
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:20
Node_t< typename GlobalBasis::PreBasis, PrefixPath > Tree
Tree of local finite elements / local shape function sets.
Definition: LocalView.hpp:41
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:135
Tree const & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: LocalView.hpp:129
GlobalBasis const & globalBasis() const
Return the global basis that we are a view on.
Definition: LocalView.hpp:165
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:112
LocalView(LocalView const &other)
Copy constructor.
Definition: LocalView.hpp:64
void bind(Element const &element)
Bind the view to a grid element.
Definition: LocalView.hpp:91
typename PreBasis::GridView GridView
The grid view that the FE space is defined on.
Definition: GlobalBasis.hpp:63
void unbind()
Unbind from the current element.
Definition: LocalView.hpp:122
LocalView(GlobalBasis const &globalBasis)
Construct local view for a given global finite element basis.
Definition: LocalView.hpp:51
GB GlobalBasis
The global FE basis that this is a view on.
Definition: LocalView.hpp:29
NodeCache_t< Tree > TreeCache
Cached basis-tree.
Definition: LocalView.hpp:44
std::size_t size_type
The type used for sizes.
Definition: LocalView.hpp:38