DUNE PDELab (git)

subspacelocalview.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_SUBSPACELOCALVIEW_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_SUBSPACELOCALVIEW_HH
5
6
7#include <tuple>
8
10
11#include <dune/typetree/childextraction.hh>
12
13#include <dune/functions/functionspacebases/concepts.hh>
14
15
16
17namespace Dune {
18namespace Functions {
19
20
21
22template<class RB, class PP>
23class SubspaceBasis;
24
25
26
28template<class RLV, class PP>
30{
31 using PrefixPath = PP;
32
33public:
34
35 using RootLocalView = RLV;
36
38 using GlobalBasis = SubspaceBasis<typename RootLocalView::GlobalBasis, PrefixPath>;
39
41 using GridView = typename GlobalBasis::GridView;
42
44 using Element = typename GridView::template Codim<0>::Entity;
45
47 using size_type = std::size_t;
48
50 using RootTree = typename RootLocalView::Tree;
51
54
56 using MultiIndex = typename RootLocalView::MultiIndex;
57
59 SubspaceLocalView(const GlobalBasis& globalBasis, const PrefixPath& /*prefixPath*/) :
60 globalBasis_(&globalBasis),
61 rootLocalView_(globalBasis.rootBasis().localView())
62 {
63// static_assert(models<Concept::BasisTree<GridView>, Tree>(), "Tree type passed to SubspaceLocalView does not model the BasisNode concept.");
64 }
65
71 void bind(const Element& e)
72 {
73 rootLocalView_.bind(e);
74 }
75
80 const Element& element() const
81 {
82 return rootLocalView_.element();
83 }
84
89 void unbind()
90 {
91 rootLocalView_.unbind();
92 }
93
96 bool bound() const
97 {
98 return rootLocalView_.bound();
99 }
100
105 const Tree& tree() const
106 {
107 return TypeTree::child(rootLocalView_.tree(), globalBasis_->prefixPath());
108 }
109
113 {
114 return rootLocalView_.size();
115 }
116
124 {
125 return rootLocalView_.maxSize();
126 }
127
130 {
131 return rootLocalView_.index(i);
132 }
133
137 {
138 return *globalBasis_;
139 }
140
141 const RootLocalView& rootLocalView() const
142 {
143 return rootLocalView_;
144 }
145
146protected:
147 const GlobalBasis* globalBasis_;
148 RootLocalView rootLocalView_;
149};
150
151
152
153} // end namespace Functions
154} // end namespace Dune
155
156
157
158#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_SUBSPACELOCALVIEW_HH
The restriction of a finite element basis to a single element.
Definition: subspacelocalview.hh:30
typename TypeTree::ChildForTreePath< RootTree, PrefixPath > Tree
Tree of local finite elements / local shape function sets.
Definition: subspacelocalview.hh:53
void unbind()
Unbind from the current element.
Definition: subspacelocalview.hh:89
const Element & element() const
Return the grid element that the view is bound to.
Definition: subspacelocalview.hh:80
bool bound() const
Return if the view is bound to a grid element.
Definition: subspacelocalview.hh:96
typename RootLocalView::MultiIndex MultiIndex
Type used for global numbering of the basis vectors.
Definition: subspacelocalview.hh:56
SubspaceBasis< typename RootLocalView::GlobalBasis, PrefixPath > GlobalBasis
The global FE basis that this is a view on.
Definition: subspacelocalview.hh:38
size_type size() const
Total number of degrees of freedom on this element.
Definition: subspacelocalview.hh:112
void bind(const Element &e)
Bind the view to a grid element.
Definition: subspacelocalview.hh:71
typename GlobalBasis::GridView GridView
The grid view the global FE basis lives on.
Definition: subspacelocalview.hh:41
SubspaceLocalView(const GlobalBasis &globalBasis, const PrefixPath &)
Construct local view for a given global finite element basis.
Definition: subspacelocalview.hh:59
std::size_t size_type
The type used for sizes.
Definition: subspacelocalview.hh:47
typename GridView::template Codim< 0 >::Entity Element
Type of the grid element we are bound to.
Definition: subspacelocalview.hh:44
const GlobalBasis & globalBasis() const
Return the global basis that we are a view on.
Definition: subspacelocalview.hh:136
size_type maxSize() const
Maximum local size for any element on the GridView.
Definition: subspacelocalview.hh:123
const Tree & tree() const
Return the local ansatz tree associated to the bound entity.
Definition: subspacelocalview.hh:105
MultiIndex index(size_type i) const
Maps from subtree index set [0..size-1] to a globally unique multi index in global basis.
Definition: subspacelocalview.hh:129
typename RootLocalView::Tree RootTree
Tree of local finite elements / local shape function sets.
Definition: subspacelocalview.hh:50
Infrastructure for concepts.
ImplementationDefined child(Node &&node, Indices... indices)
Extracts the child of a node given by a sequence of compile-time and run-time indices.
Definition: childextraction.hh:128
typename impl::_ChildForTreePath< Node, TreePath >::type ChildForTreePath
Template alias for the type of a child node given by a TreePath or a HybridTreePath type.
Definition: childextraction.hh:252
Dune namespace.
Definition: alignedallocator.hh:13
Static tag representing a codimension.
Definition: dimension.hh:24
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)