DUNE PDELab (2.7)

defaultlocalview.hh
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
10#include <dune/common/std/optional.hh>
11#include <dune/common/hybridutilities.hh>
12
13#include <dune/functions/functionspacebases/concepts.hh>
14
15
16
17namespace Dune {
18namespace Functions {
19
20
21
23template<class GB>
25{
26 // Node index set provided by PreBasis
27 using NodeIndexSet = typename GB::PreBasis::IndexSet;
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
47 using MultiIndex = typename NodeIndexSet::MultiIndex;
48
49private:
50
51 template<typename NodeIndexSet_>
52 using hasIndices = decltype(std::declval<NodeIndexSet_>().indices(std::declval<std::vector<typename NodeIndexSet_::MultiIndex>>().begin()));
53
54public:
55
58 globalBasis_(&globalBasis),
59 tree_(globalBasis_->preBasis().makeNode()),
60 nodeIndexSet_(globalBasis_->preBasis().makeIndexSet())
61 {
62 static_assert(models<Concept::BasisTree<GridView>, Tree>(), "Tree type passed to DefaultLocalView does not model the BasisNode concept.");
63 initializeTree(tree_);
64 }
65
71 void bind(const Element& e)
72 {
73 element_ = e;
74 bindTree(tree_, *element_);
75 nodeIndexSet_.bind(tree_);
76 indices_.resize(size());
79 [&](auto id) {
80 id(nodeIndexSet_).indices(indices_.begin());
81 },
82 [&](auto id) {
83 for (size_type i = 0 ; i < this->size() ; ++i)
84 indices_[i] = id(nodeIndexSet_).index(i);
85 });
86 }
87
90 bool isBound() const {
91 return static_cast<bool>(element_);
92 }
93
98 const Element& element() const
99 {
100 return *element_;
101 }
102
107 void unbind()
108 {
109 nodeIndexSet_.unbind();
110 element_.reset();
111 }
112
117 const Tree& tree() const
118 {
119 return tree_;
120 }
121
125 {
126 return tree_.size();
127 }
128
136 {
137 return globalBasis_->preBasis().maxNodeSize();
138 }
139
142 {
143 return indices_[i];
144 }
145
149 {
150 return *globalBasis_;
151 }
152
153 const DefaultLocalView& rootLocalView() const
154 {
155 return *this;
156 }
157
158protected:
159 const GlobalBasis* globalBasis_;
160 Std::optional<Element> element_;
161 Tree tree_;
162 NodeIndexSet nodeIndexSet_;
163 std::vector<MultiIndex> indices_;
164};
165
166
167
168} // end namespace Functions
169} // end namespace Dune
170
171
172
173#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_DEFAULTLOCALVIEW_HH
The restriction of a finite element basis to a single element.
Definition: defaultlocalview.hh:25
void unbind()
Unbind from the current element.
Definition: defaultlocalview.hh:107
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: defaultlocalview.hh:35
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:117
void bind(const Element &e)
Bind the view to a grid element.
Definition: defaultlocalview.hh:71
const Element & element() const
Return the grid element that the view is bound to.
Definition: defaultlocalview.hh:98
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:141
size_type size() const
Total number of degrees of freedom on this element.
Definition: defaultlocalview.hh:124
GB GlobalBasis
The global FE basis that this is a view on.
Definition: defaultlocalview.hh:32
size_type maxSize() const
Maximum local size for any element on the GridView.
Definition: defaultlocalview.hh:135
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:90
typename NodeIndexSet::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: defaultlocalview.hh:47
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:57
const GlobalBasis & globalBasis() const
Return the global basis that we are a view on.
Definition: defaultlocalview.hh:148
Infrastructure for concepts.
typename detected_or< nonesuch, Op, Args... >::value_t is_detected
Detects whether Op<Args...> is valid.
Definition: type_traits.hh:340
decltype(auto) ifElse(const Condition &condition, IfFunc &&ifFunc, ElseFunc &&elseFunc)
A conditional expression.
Definition: hybridutilities.hh:355
Dune namespace.
Definition: alignedallocator.hh:14
Static tag representing a codimension.
Definition: dimension.hh:22
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)