dune-functions 2.9.0
defaultlocalview.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_FUNCTIONS_FUNCTIONSPACEBASES_DEFAULTLOCALVIEW_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_DEFAULTLOCALVIEW_HH
5
6
7#include <tuple>
8#include <optional>
9
10#include <dune/common/concept.hh>
11#include <dune/common/hybridutilities.hh>
12#include <dune/common/reservedvector.hh>
13
17
18
19
20namespace Dune {
21namespace Functions {
22
23
24
26template<class GB>
28{
29public:
30
32 using GlobalBasis = GB;
33
35 using GridView = typename GlobalBasis::GridView;
36
38 using Element = typename GridView::template Codim<0>::Entity;
39
41 using size_type = std::size_t;
42
44 using Tree = typename GlobalBasis::PreBasis::Node;
45
46protected:
47
48 using PreBasis = typename GlobalBasis::PreBasis;
49
50 // Type used to store the multi indices of the basis vectors.
51 // In contrast to MultiIndex this always has dynamic size.
52 // It's guaranteed, that you can always cast it to MultiIndex
54 std::conditional_t<(PreBasis::minMultiIndexSize == PreBasis::maxMultiIndexSize),
56 Dune::ReservedVector<size_type, PreBasis::multiIndexBufferSize>>;
57
58public:
59
61 using MultiIndex =
62 std::conditional_t<(PreBasis::minMultiIndexSize == PreBasis::maxMultiIndexSize),
64 Dune::ReservedVector<size_type, PreBasis::multiIndexBufferSize>>;
65
66
70 tree_(globalBasis_->preBasis().makeNode())
71 {
72 static_assert(models<Concept::BasisTree<GridView>, Tree>(), "Tree type passed to DefaultLocalView does not model the BasisNode concept.");
74 }
75
81 void bind(const Element& e)
82 {
83 element_ = e;
85 indices_.resize(size());
86 globalBasis_->preBasis().indices(tree_, indices_.begin());
87 }
88
91 [[deprecated("Use the bound() method instead")]]
92 bool isBound() const {
93 return static_cast<bool>(element_);
94 }
95
98 bool bound() const
99 {
100 return static_cast<bool>(element_);
101 }
102
107 const Element& element() const
108 {
109 return *element_;
110 }
111
116 void unbind()
117 {
118 element_.reset();
119 }
120
125 const Tree& tree() const
126 {
127 return tree_;
128 }
129
133 {
134 return tree_.size();
135 }
136
144 {
145 return globalBasis_->preBasis().maxNodeSize();
146 }
147
149 const MultiIndex& index(size_type i) const
150 {
151 return indices_[i];
152 }
153
157 {
158 return *globalBasis_;
159 }
160
162 {
163 return *this;
164 }
165
166protected:
168 std::optional<Element> element_;
170 std::vector<MultiIndexStorage> indices_;
171};
172
173
174
175} // end namespace Functions
176} // end namespace Dune
177
178
179
180#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_DEFAULTLOCALVIEW_HH
Definition: polynomial.hh:10
void bindTree(Tree &tree, const Entity &entity, std::size_t offset=0)
Definition: nodes.hh:253
void initializeTree(Tree &tree, std::size_t treeIndexOffset=0)
Definition: nodes.hh:260
A statically sized MultiIndex type.
Definition: multiindex.hh:25
A dynamically sized array-like class with overflow.
Definition: overflowarray.hh:45
The restriction of a finite element basis to a single element.
Definition: defaultlocalview.hh:28
typename GlobalBasis::PreBasis PreBasis
Definition: defaultlocalview.hh:48
void unbind()
Unbind from the current element.
Definition: defaultlocalview.hh:116
bool bound() const
Return if the view is bound to a grid element.
Definition: defaultlocalview.hh:98
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: defaultlocalview.hh:35
const MultiIndex & index(size_type i) const
Maps from subtree index set [0..size-1] to a globally unique multi index in global basis.
Definition: defaultlocalview.hh:149
std::optional< Element > element_
Definition: defaultlocalview.hh:168
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: defaultlocalview.hh:38
const Tree & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: defaultlocalview.hh:125
void bind(const Element &e)
Bind the view to a grid element.
Definition: defaultlocalview.hh:81
const Element & element() const
Return the grid element that the view is bound to.
Definition: defaultlocalview.hh:107
size_type size() const
Total number of degrees of freedom on this element.
Definition: defaultlocalview.hh:132
GB GlobalBasis
The global FE basis that this is a view on.
Definition: defaultlocalview.hh:32
Tree tree_
Definition: defaultlocalview.hh:169
size_type maxSize() const
Maximum local size for any element on the GridView.
Definition: defaultlocalview.hh:143
std::size_t size_type
The type used for sizes.
Definition: defaultlocalview.hh:41
bool isBound() const
Return if the view is bound to a grid element.
Definition: defaultlocalview.hh:92
const DefaultLocalView & rootLocalView() const
Definition: defaultlocalview.hh:161
std::conditional_t<(PreBasis::minMultiIndexSize==PreBasis::maxMultiIndexSize), OverflowArray< StaticMultiIndex< size_type, PreBasis::maxMultiIndexSize >, PreBasis::multiIndexBufferSize >, Dune::ReservedVector< size_type, PreBasis::multiIndexBufferSize > > MultiIndexStorage
Definition: defaultlocalview.hh:56
std::conditional_t<(PreBasis::minMultiIndexSize==PreBasis::maxMultiIndexSize), StaticMultiIndex< size_type, PreBasis::maxMultiIndexSize >, Dune::ReservedVector< size_type, PreBasis::multiIndexBufferSize > > MultiIndex
Type used for global numbering of the basis vectors.
Definition: defaultlocalview.hh:64
std::vector< MultiIndexStorage > indices_
Definition: defaultlocalview.hh:170
typename GlobalBasis::PreBasis::Node Tree
Tree of local finite elements / local shape function sets.
Definition: defaultlocalview.hh:44
DefaultLocalView(const GlobalBasis &globalBasis)
Construct local view for a given global finite element basis.
Definition: defaultlocalview.hh:68
const GlobalBasis * globalBasis_
Definition: defaultlocalview.hh:167
const GlobalBasis & globalBasis() const
Return the global basis that we are a view on.
Definition: defaultlocalview.hh:156