AMDiS  0.1
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/functions/NodeCache.hpp>
11 #include <amdis/functions/Nodes.hpp>
12 
13 // NOTE: This is a variant of dune-functions DefaultLocalView
14 
15 namespace AMDiS
16 {
18  template <class GB>
19  class LocalView
20  {
21  using PrefixPath = Dune::TypeTree::HybridTreePath<>;
22 
23  // Node index set provided by PreBasis
24  using NodeIndexSet = NodeIndexSet_t<typename GB::PreBasis, PrefixPath>;
25 
26  public:
28  using GlobalBasis = GB;
29 
31  using GridView = typename GlobalBasis::GridView;
32 
34  using Element = typename GridView::template Codim<0>::Entity;
35 
37  using size_type = std::size_t;
38 
40  using Tree = Node_t<typename GlobalBasis::PreBasis, PrefixPath>;
41 
44 
47 
48  private:
49  template <class NIS>
50  using hasIndices = decltype(std::declval<NIS>().indices(std::declval<std::vector<typename NIS::MultiIndex>>().begin()));
51 
52  public:
55  : globalBasis_(&globalBasis)
56  , tree_(makeNode(globalBasis_->preBasis(), PrefixPath{}))
57  , treeCache_(makeNodeCache(tree_))
58  , nodeIndexSet_(makeNodeIndexSet(globalBasis_->preBasis(), PrefixPath{}))
59  {
60  static_assert(Dune::models<Dune::Functions::Concept::BasisTree<GridView>, Tree>());
61  initializeTree(tree_);
62  }
63 
65 
70  void bind (Element const& element)
71  {
72  element_ = element;
73  bindTree(tree_, *element_);
74  nodeIndexSet_.bind(tree_);
75  indices_.resize(size());
76 
77  if constexpr (Dune::Std::is_detected<hasIndices,NodeIndexSet>{})
78  nodeIndexSet_.indices(indices_.begin());
79  else
80  for (size_type i = 0; i < size(); ++i)
81  indices_[i] = nodeIndexSet_.index(i);
82  }
83 
85  bool isBound () const
86  {
87  return bool(element_);
88  }
89 
91  Element const& element () const
92  {
93  return *element_;
94  }
95 
97 
100  void unbind ()
101  {
102  nodeIndexSet_.unbind();
103  element_.reset();
104  }
105 
107  Tree const& tree () const
108  {
109  return tree_;
110  }
111 
113  TreeCache const& treeCache () const
114  {
115  return treeCache_;
116  }
117 
119  size_type size () const
120  {
121  return tree_.size();
122  }
123 
125 
129  {
130  return globalBasis_->preBasis().maxNodeSize();
131  }
132 
136  {
137  return indices_[i];
138  }
139 
141  GlobalBasis const& globalBasis () const
142  {
143  return *globalBasis_;
144  }
145 
147  LocalView const& rootLocalView () const
148  {
149  return *this;
150  }
151 
152  protected:
153  GlobalBasis const* globalBasis_;
154  std::optional<Element> element_;
155  Tree tree_;
156  TreeCache treeCache_;
157  NodeIndexSet nodeIndexSet_;
158  std::vector<MultiIndex> indices_;
159  };
160 
161 } // end namespace AMDiS
size_type size() const
Total number of degrees of freedom on this element.
Definition: LocalView.hpp:119
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:128
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: LocalView.hpp:31
auto makeNodeCache(Node const &node)
Construct a new local-basis cache from a basis-node.
Definition: NodeCache.hpp:35
LocalView const & rootLocalView() const
Return this local-view.
Definition: LocalView.hpp:147
bool isBound() const
Return if the view is bound to a grid element.
Definition: LocalView.hpp:85
Contains all classes needed for solving linear and non linear equation systems.
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:135
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: LocalView.hpp:34
The restriction of a finite element basis to a single element.
Definition: LocalView.hpp:19
typename Impl::NodeCacheFactory< Node >::type NodeCache_t
Defines the type of a node cache associated to a given Node.
Definition: NodeCache.hpp:31
Node_t< typename GlobalBasis::PreBasis, PrefixPath > Tree
Tree of local finite elements / local shape function sets.
Definition: LocalView.hpp:40
TreeCache const & treeCache() const
Cached version of the local ansatz tree.
Definition: LocalView.hpp:113
Tree const & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: LocalView.hpp:107
GlobalBasis const & globalBasis() const
Return the global basis that we are a view on.
Definition: LocalView.hpp:141
Element const & element() const
Return the grid element that the view is bound to.
Definition: LocalView.hpp:91
void bind(Element const &element)
Bind the view to a grid element.
Definition: LocalView.hpp:70
typename PreBasis::GridView GridView
The grid view that the FE space is defined on.
Definition: GlobalBasis.hpp:66
void unbind()
Unbind from the current element.
Definition: LocalView.hpp:100
typename NodeIndexSet::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: LocalView.hpp:46
LocalView(GlobalBasis const &globalBasis)
Construct local view for a given global finite element basis.
Definition: LocalView.hpp:54
GB GlobalBasis
The global FE basis that this is a view on.
Definition: LocalView.hpp:28
NodeCache_t< Tree > TreeCache
Cached basis-tree.
Definition: LocalView.hpp:43
std::size_t size_type
The type used for sizes.
Definition: LocalView.hpp:37