DUNE-FUNCTIONS (2.7)

discreteglobalbasisfunction.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_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
4#define DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
5
6#include <memory>
7
8#include <dune/common/shared_ptr.hh>
9#include <dune/common/hybridutilities.hh>
10
11#include <dune/typetree/treecontainer.hh>
12
13#include <dune/functions/functionspacebases/hierarchicnodetorangemap.hh>
14#include <dune/functions/functionspacebases/flatvectorview.hh>
15#include <dune/functions/gridfunctions/gridviewentityset.hh>
16#include <dune/functions/gridfunctions/gridfunction.hh>
17#include <dune/functions/backends/concepts.hh>
18#include <dune/functions/backends/istlvectorbackend.hh>
19
20namespace Dune {
21namespace Functions {
22
23
24
67template<typename B, typename V,
68 typename NTRE = HierarchicNodeToRangeMap,
69 typename R = typename V::value_type>
71{
72public:
73 using Basis = B;
74 using Vector = V;
75
76 using Coefficient = std::decay_t<decltype(std::declval<Vector>()[std::declval<typename Basis::MultiIndex>()])>;
77
78 using GridView = typename Basis::GridView;
80 using Tree = typename Basis::LocalView::Tree;
81 using NodeToRangeEntry = NTRE;
82
83 using Domain = typename EntitySet::GlobalCoordinate;
84 using Range = R;
85
86 using LocalDomain = typename EntitySet::LocalCoordinate;
87 using Element = typename EntitySet::Element;
88
89 using Traits = Imp::GridFunctionTraits<Range(Domain), EntitySet, DefaultDerivativeTraits, 16>;
90
91 class LocalFunction
92 {
93 using LocalView = typename Basis::LocalView;
94 using size_type = typename Tree::size_type;
95
96 template<class Node>
97 using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
98
99 template<class Node>
100 using NodeData = typename std::vector<LocalBasisRange<Node>>;
101
102 using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData,Tree>;
103
104 public:
105
106 using GlobalFunction = DiscreteGlobalBasisFunction;
107 using Domain = LocalDomain;
108 using Range = GlobalFunction::Range;
109 using Element = GlobalFunction::Element;
110
111 LocalFunction(const DiscreteGlobalBasisFunction& globalFunction)
112 : globalFunction_(&globalFunction)
113 , localView_(globalFunction.basis().localView())
114 , bound_(false)
115 {}
116
117 LocalFunction(const LocalFunction& other)
118 : globalFunction_(other.globalFunction_)
119 , localView_(other.localView_)
120 , bound_(other.bound_)
121 {}
122
123 LocalFunction operator=(const LocalFunction& other)
124 {
125 globalFunction_ = other.globalFunction_;
126 localView_ = other.localView_;
127 bound_ = other.bound_;
128 }
129
136 void bind(const Element& element)
137 {
138 localView_.bind(element);
139 bound_ = true;
140 }
141
142 void unbind()
143 {
144 localView_.unbind();
145 bound_ = false;
146 }
147
151 bool bound() const {
152 return bound_;
153 }
154
164 Range operator()(const Domain& x) const
165 {
166 auto y = Range(0);
167
168 TypeTree::forEachLeafNode(localView_.tree(), [&](auto&& node, auto&& treePath) {
169 const auto& dofs = globalFunction_->dofs();
170 const auto& nodeToRangeEntry = globalFunction_->nodeToRangeEntry();
171 const auto& fe = node.finiteElement();
172 const auto& localBasis = fe.localBasis();
173 auto& shapeFunctionValues = evaluationBuffer_[treePath];
174
175 localBasis.evaluateFunction(x, shapeFunctionValues);
176
177 // Get range entry associated to this node
178 auto re = flatVectorView(nodeToRangeEntry(node, treePath, y));
179
180 for (size_type i = 0; i < localBasis.size(); ++i)
181 {
182 const auto& multiIndex = localView_.index(node.localIndex(i));
183
184 // Get coefficient associated to i-th shape function
185 auto c = flatVectorView(dofs[multiIndex]);
186
187 // Get value of i-th shape function
188 auto v = flatVectorView(shapeFunctionValues[i]);
189
190 // Notice that the range entry re, the coefficient c, and the shape functions
191 // value v may all be scalar, vector, matrix, or general container valued.
192 // The matching of their entries is done via the multistage procedure described
193 // in the class documentation of DiscreteGlobalBasisFunction.
194 auto&& dimC = c.size();
195 auto dimV = v.size();
196 assert(dimC*dimV == re.size());
197 for(size_type j=0; j<dimC; ++j)
198 {
199 auto&& c_j = c[j];
200 for(size_type k=0; k<dimV; ++k)
201 re[j*dimV + k] += c_j*v[k];
202 }
203 }
204 });
205
206 return y;
207 }
208
209 const Element& localContext() const
210 {
211 return localView_.element();
212 }
213
214 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction& t)
215 {
216 DUNE_THROW(NotImplemented,"not implemented");
217 }
218
219 private:
220
221 const DiscreteGlobalBasisFunction* globalFunction_;
222 LocalView localView_;
223 mutable PerNodeEvaluationBuffer evaluationBuffer_;
224 bool bound_ = false;
225 };
226
227 template<class B_T, class V_T, class NTRE_T>
228 DiscreteGlobalBasisFunction(B_T && basis, V_T && coefficients, NTRE_T&& nodeToRangeEntry) :
229 entitySet_(basis.gridView()),
230 basis_(wrap_or_move(std::forward<B_T>(basis))),
231 coefficients_(wrap_or_move(std::forward<V_T>(coefficients))),
232 nodeToRangeEntry_(wrap_or_move(std::forward<NTRE_T>(nodeToRangeEntry)))
233 {}
234
235 DiscreteGlobalBasisFunction(std::shared_ptr<const Basis> basis, std::shared_ptr<const V> coefficients, std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry) :
236 entitySet_(basis->gridView()),
237 basis_(basis),
238 coefficients_(coefficients),
239 nodeToRangeEntry_(nodeToRangeEntry)
240 {}
241
242 const Basis& basis() const
243 {
244 return *basis_;
245 }
246
247 const V& dofs() const
248 {
249 return *coefficients_;
250 }
251
252 const NodeToRangeEntry& nodeToRangeEntry() const
253 {
254 return *nodeToRangeEntry_;
255 }
256
257 // TODO: Implement this using hierarchic search
258 Range operator() (const Domain& x) const
259 {
260 DUNE_THROW(NotImplemented,"not implemented");
261 }
262
263 friend typename Traits::DerivativeInterface derivative(const DiscreteGlobalBasisFunction& t)
264 {
265 DUNE_THROW(NotImplemented,"not implemented");
266 }
267
283 friend LocalFunction localFunction(const DiscreteGlobalBasisFunction& t)
284 {
285 return LocalFunction(t);
286 }
287
291 const EntitySet& entitySet() const
292 {
293 return entitySet_;
294 }
295
296private:
297
298 EntitySet entitySet_;
299 std::shared_ptr<const Basis> basis_;
300 std::shared_ptr<const V> coefficients_;
301 std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry_;
302};
303
304
305
316template<typename... TT>
318
319
320
321template<typename R, typename B, typename V>
322auto makeDiscreteGlobalBasisFunction(B&& basis, V&& vector)
323{
324 using Basis = std::decay_t<B>;
325 using NTREM = HierarchicNodeToRangeMap;
326
327 // Small helper functions to wrap vectors using istlVectorBackend
328 // if they do not already satisfy the VectorBackend interface.
329 auto toConstVectorBackend = [&](auto&& v) -> decltype(auto) {
330 return Dune::Hybrid::ifElse(models<Concept::ConstVectorBackend<Basis>, decltype(v)>(),
331 [&](auto id) -> decltype(auto) {
332 return std::forward<decltype(v)>(v);
333 }, [&](auto id) -> decltype(auto) {
334 return istlVectorBackend(id(v));
335 });
336 };
337
338 using Vector = std::decay_t<decltype(toConstVectorBackend(std::forward<V>(vector)))>;
340 std::forward<B>(basis),
341 toConstVectorBackend(std::forward<V>(vector)),
343}
344
345
346
347} // namespace Functions
348} // namespace Dune
349
350#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
A grid function induced by a global basis and a coefficient vector.
Definition: discreteglobalbasisfunction.hh:71
const EntitySet & entitySet() const
Get associated EntitySet.
Definition: discreteglobalbasisfunction.hh:291
GridView::template Codim< codim >::Entity Element
Type of Elements contained in this EntitySet.
Definition: gridviewentityset.hh:32
Element::Geometry::LocalCoordinate LocalCoordinate
Type of local coordinates with respect to the Element.
Definition: gridviewentityset.hh:35
friend LocalFunction localFunction(const DiscreteGlobalBasisFunction &t)
Construct local function from a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:283
void localFunction(DiscreteGlobalBasisFunction< TT... > &&t)=delete
Construction of local functions from a temporary DiscreteGlobalBasisFunction (forbidden)
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:39
auto istlVectorBackend(Vector &v)
Return a vector backend wrapping non-const ISTL like containers.
Definition: istlvectorbackend.hh:345
Definition: polynomial.hh:10
Default implementation for derivative traits.
Definition: defaultderivativetraits.hh:37
A simple node to range map using the nested tree indices.
Definition: hierarchicnodetorangemap.hh:30
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)