DUNE PDELab (2.7)

dunefunctionslocalfunctionspace.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_PDELAB_GRIDFUNCTIONSPACE_DUNEFUNCTIONSLOCALFUNCTIONSPACE_HH
4#define DUNE_PDELAB_GRIDFUNCTIONSPACE_DUNEFUNCTIONSLOCALFUNCTIONSPACE_HH
5
6#include<vector>
7
10
11#include <dune/geometry/referenceelements.hh>
12
13#include <dune/localfunctions/common/interfaceswitch.hh>
14#include <dune/localfunctions/common/localkey.hh>
15
16#include <dune/typetree/typetree.hh>
17
18#include <dune/pdelab/gridfunctionspace/tags.hh>
20
21namespace Dune {
22 namespace PDELab {
23
27
28 namespace Experimental {
29
30 template<typename LFS>
31 struct LeafLFSMixin
32 : public TypeTree::LeafNode
33 {
34
35 const auto& finiteElement() const
36 {
37 return static_cast<const LFS*>(this)->tree().finiteElement();
38 }
39
40 template<typename Tree>
41 struct Traits
42 {
43 using FiniteElement = typename Tree::FiniteElement;
44 using FiniteElementType = FiniteElement;
45 };
46 };
47
48 template<typename GFS, typename TreePath = TypeTree::HybridTreePath<>>
49 class LocalFunctionSpace
50 : public LeafLFSMixin<LocalFunctionSpace<GFS,TreePath>>
51 {
52
53 public:
54
55 using Basis = typename GFS::Basis;
56 using LocalView = typename Basis::LocalView;
57 using Tree = TypeTree::ChildForTreePath<typename LocalView::Tree,TreePath>;
58 using DOFIndex = typename GFS::Ordering::Traits::DOFIndex;
59
60 template<typename LFS, typename C, typename Tag, bool fast>
61 friend class Dune::PDELab::LFSIndexCacheBase;
62
63 struct Traits
64 : public LeafLFSMixin<LocalFunctionSpace<GFS,TreePath>>::template Traits<Tree>
65 {
66
67 using GridFunctionSpace = GFS;
68 using GridView = typename GFS::Traits::GridView;
69 using SizeType = std::size_t;
70 using DOFIndex = typename GFS::Ordering::Traits::DOFIndex;
71 using ConstraintsType = typename GFS::Traits::ConstraintsType;
72
73 };
74
75 using size_type = std::size_t;
76
77 LocalFunctionSpace(std::shared_ptr<const GFS> gfs, TreePath tree_path = TreePath(), size_type offset = 0)
78 : _gfs(gfs)
79 , _local_view(gfs->basis().localView())
80 , _tree_path(tree_path)
81 , _tree(TypeTree::child(_local_view.tree(),tree_path))
82 {}
83
84 size_type subSpaceDepth() const
85 {
86 return 0;
87 }
88
90 size_type size () const
91 {
92 return _local_view.size();
93 }
94
95 size_type maxSize () const
96 {
97 // _dof_indices is always as large as the max local size of the root GFS
98 return _local_view.maxSize();
99 }
100
102 size_type localIndex (size_type index) const
103 {
104 return _tree.localIndex(index);
105 }
106
107 // index: local dof index for the given element
108 DOFIndex dofIndex(size_type index) const
109 {
110 auto refElement = Dune::ReferenceElements<double,Basis::GridView::dimension>::general(_local_view.element().type());
111
112 auto localKey = _local_view.tree().finiteElement().localCoefficients().localKey(index);
113
114 const auto& indexSet = _gfs->basis().gridView().indexSet();
115
116 // get geometry type of subentity
117 auto gt = refElement.type(localKey.subEntity(), localKey.codim());
118
119 // evaluate consecutive index of subentity
120 auto indexOnEntity = indexSet.subIndex(_local_view.element(),
121 localKey.subEntity(),
122 localKey.codim());
123
124
125 DOFIndex result;
126 GFS::Ordering::Traits::DOFIndexAccessor::store(result,gt,indexOnEntity,localKey.index());
127 return result;
128 }
129
130 // index: local dof index for the given element
131 auto containerIndex(size_type index) const
132 {
133 MultiIndex<std::size_t,1> result;
134 result.set({_local_view.index(_tree.localIndex(index))});
135 return result;
136 }
137
139 const GFS& gridFunctionSpace() const
140 {
141 return *_gfs;
142 }
143
144 void bind(const typename GFS::Traits::EntitySet::template Codim<0>::Entity& e)
145 {
146 _local_view.bind(e);
147 }
148
149 const typename Traits::ConstraintsType& constraints() const
150 {
151 return _gfs->constraints();
152 }
153
154 const Tree& tree() const
155 {
156 return _tree;
157 }
158
159 private:
160
161 typename GFS::Ordering::Traits::ContainerIndex containerIndex(const DOFIndex& i) const
162 {
163 return _gfs->ordering().containerIndex(i);
164 }
165
166 std::shared_ptr<const GFS> _gfs;
167 LocalView _local_view;
168 TreePath _tree_path;
169 const Tree& _tree;
170
171 };
172
173 // forward declare GridFunctionSpace
174 template<typename DFBasis, typename V, typename CE=NoConstraints>
175 class GridFunctionSpace;
176
177
178 } // namespace Experimental
179
180
181 template<typename DFBasis, typename V, typename CE, typename TAG>
182 class LocalFunctionSpace<Experimental::GridFunctionSpace<DFBasis,V,CE>,TAG>
183 : public Experimental::LocalFunctionSpace<Experimental::GridFunctionSpace<DFBasis,V,CE>>
184 {
185
186 using GFS = Experimental::GridFunctionSpace<DFBasis,V,CE>;
187
188 public:
189
190 LocalFunctionSpace(std::shared_ptr<const GFS> gfs)
191 : Experimental::LocalFunctionSpace<GFS>(gfs)
192 {}
193
194 LocalFunctionSpace(const GFS& gfs)
195 : Experimental::LocalFunctionSpace<GFS>(stackobject_to_shared_ptr(gfs))
196 {}
197
198 };
199
200 template<typename DFBasis, typename V, typename CE>
201 class LocalFunctionSpace<Experimental::GridFunctionSpace<DFBasis,V,CE>,AnySpaceTag>
202 : public Experimental::LocalFunctionSpace<Experimental::GridFunctionSpace<DFBasis,V,CE>>
203 {
204
205 using GFS = Experimental::GridFunctionSpace<DFBasis,V,CE>;
206
207 public:
208
209 LocalFunctionSpace(std::shared_ptr<const GFS> gfs)
210 : Experimental::LocalFunctionSpace<GFS>(gfs)
211 {}
212
213 LocalFunctionSpace(const GFS& gfs)
214 : Experimental::LocalFunctionSpace<GFS>(stackobject_to_shared_ptr(gfs))
215 {}
216
217 };
218
220 } // namespace PDELab
221} // namespace Dune
222
223#endif // DUNE_PDELAB_GRIDFUNCTIONSPACE_DUNEFUNCTIONSLOCALFUNCTIONSPACE_HH
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:156
void constraints(const GFS &gfs, CG &cg, const bool verbose=false)
construct constraints
Definition: constraints.hh:749
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:179
Dune namespace.
Definition: alignedallocator.hh:14
shared_ptr< T > stackobject_to_shared_ptr(T &t)
Create a shared_ptr for a stack-allocated object.
Definition: shared_ptr.hh:75
This file implements the class shared_ptr (a reference counting pointer), for those systems that don'...
Standard Dune debug streams.
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:196
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)