DUNE-FUNCTIONS (unstable)

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
4// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
5// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
6
7#ifndef DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
8#define DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
9
10#include <memory>
11#include <optional>
12
13#include <dune/common/typetraits.hh>
14
15#include <dune/grid/utility/hierarchicsearch.hh>
16
17#include <dune/typetree/treecontainer.hh>
18
19#include <dune/functions/functionspacebases/hierarchicnodetorangemap.hh>
20#include <dune/functions/functionspacebases/flatvectorview.hh>
21#include <dune/functions/gridfunctions/gridviewentityset.hh>
22#include <dune/functions/gridfunctions/gridfunction.hh>
23#include <dune/functions/backends/concepts.hh>
24#include <dune/functions/backends/istlvectorbackend.hh>
25
26namespace Dune {
27namespace Functions {
28
29
30namespace ImplDoc {
31
32template<typename B, typename V, typename NTRE>
33class DiscreteGlobalBasisFunctionBase
34{
35public:
36 using Basis = B;
37 using Vector = V;
38
39 // In order to make the cache work for proxy-references
40 // we have to use AutonomousValue<T> instead of std::decay_t<T>
41 using Coefficient = Dune::AutonomousValue<decltype(std::declval<Vector>()[std::declval<typename Basis::MultiIndex>()])>;
42
43 using GridView = typename Basis::GridView;
44 using EntitySet = GridViewEntitySet<GridView, 0>;
45 using Tree = typename Basis::LocalView::Tree;
46 using NodeToRangeEntry = NTRE;
47
48 using Domain = typename EntitySet::GlobalCoordinate;
49
50 using LocalDomain = typename EntitySet::LocalCoordinate;
51 using Element = typename EntitySet::Element;
52
53protected:
54
55 // This collects all data that is shared by all related
56 // global and local functions. This way we don't need to
57 // keep track of it individually.
58 struct Data
59 {
60 EntitySet entitySet;
61 std::shared_ptr<const Basis> basis;
62 std::shared_ptr<const Vector> coefficients;
63 std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry;
64 };
65
66public:
67 class LocalFunctionBase
68 {
69 using LocalView = typename Basis::LocalView;
70 using size_type = typename Tree::size_type;
71
72 public:
73 using Domain = LocalDomain;
74 using Element = typename EntitySet::Element;
75
76 protected:
77 LocalFunctionBase(const std::shared_ptr<const Data>& data)
78 : data_(data)
79 , localView_(data_->basis->localView())
80 {
81 localDoFs_.reserve(localView_.maxSize());
82 }
83
90 LocalFunctionBase(const LocalFunctionBase& other)
91 : data_(other.data_)
92 , localView_(other.localView_)
93 {
94 localDoFs_.reserve(localView_.maxSize());
95 if (bound())
96 localDoFs_ = other.localDoFs_;
97 }
98
106 LocalFunctionBase& operator=(const LocalFunctionBase& other)
107 {
108 data_ = other.data_;
109 localView_ = other.localView_;
110 if (bound())
111 localDoFs_ = other.localDoFs_;
112 return *this;
113 }
114
115 public:
122 void bind(const Element& element)
123 {
124 localView_.bind(element);
125 // Use cache of full local view size. For a subspace basis,
126 // this may be larger than the number of local DOFs in the
127 // tree. In this case only cache entries associated to local
128 // DOFs in the subspace are filled. Cache entries associated
129 // to local DOFs which are not contained in the subspace will
130 // not be touched.
131 //
132 // Alternatively one could use a cache that exactly fits
133 // the size of the tree. However, this would require to
134 // subtract an offset from localIndex(i) on each cache
135 // access in operator().
136 localDoFs_.resize(localView_.size());
137 const auto& dofs = *data_->coefficients;
138 for (size_type i = 0; i < localView_.tree().size(); ++i)
139 {
140 // For a subspace basis the index-within-tree i
141 // is not the same as the localIndex within the
142 // full local view.
143 size_t localIndex = localView_.tree().localIndex(i);
144 localDoFs_[localIndex] = dofs[localView_.index(localIndex)];
145 }
146 }
147
149 void unbind()
150 {
151 localView_.unbind();
152 }
153
155 bool bound() const
156 {
157 return localView_.bound();
158 }
159
161 const Element& localContext() const
162 {
163 return localView_.element();
164 }
165
166 protected:
167
168 template<class To, class From>
169 void assignWith(To& to, const From& from) const
170 {
171 auto from_flat = flatVectorView(from);
172 auto to_flat = flatVectorView(to);
173 assert(from_flat.size() == to_flat.size());
174 for (size_type i = 0; i < to_flat.size(); ++i)
175 to_flat[i] = from_flat[i];
176 }
177
178 template<class Node, class TreePath, class Range>
179 decltype(auto) nodeToRangeEntry(const Node& node, const TreePath& treePath, Range& y) const
180 {
181 return (*data_->nodeToRangeEntry)(node, treePath, y);
182 }
183
184 std::shared_ptr<const Data> data_;
185 LocalView localView_;
186 std::vector<Coefficient> localDoFs_;
187 };
188
189protected:
190 DiscreteGlobalBasisFunctionBase(const std::shared_ptr<const Data>& data)
191 : data_(data)
192 {
193 /* Nothing. */
194 }
195
196public:
197
199 const Basis& basis() const
200 {
201 return *data_->basis;
202 }
203
205 const Vector& dofs() const
206 {
207 return *data_->coefficients;
208 }
209
211 const NodeToRangeEntry& nodeToRangeEntry() const
212 {
213 return *data_->nodeToRangeEntry;
214 }
215
217 const EntitySet& entitySet() const
218 {
219 return data_->entitySet;
220 }
221
222protected:
223 std::shared_ptr<const Data> data_;
224};
225
226} // namespace ImplDoc
227
228
229
230template<typename DGBF>
231class DiscreteGlobalBasisFunctionDerivative;
232
276template<typename B, typename V,
277 typename NTRE = HierarchicNodeToRangeMap,
278 typename R = typename V::value_type>
280 : public ImplDoc::DiscreteGlobalBasisFunctionBase<B, V, NTRE>
281{
282 using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<B, V, NTRE>;
283 using Data = typename Base::Data;
284
285public:
286 using Basis = typename Base::Basis;
287 using Vector = typename Base::Vector;
288
289 using Domain = typename Base::Domain;
290 using Range = R;
291
292 using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
293
294private:
295
296 template<class Node>
297 using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
298 template<class Node>
299 using NodeData = typename std::vector<LocalBasisRange<Node>>;
300 using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
301
302public:
303 class LocalFunction
304 : public Base::LocalFunctionBase
305 {
306 using LocalBase = typename Base::LocalFunctionBase;
307 using size_type = typename Base::Tree::size_type;
308 using LocalBase::nodeToRangeEntry;
309
310 public:
311
312 using GlobalFunction = DiscreteGlobalBasisFunction;
313 using Domain = typename LocalBase::Domain;
314 using Range = GlobalFunction::Range;
315 using Element = typename LocalBase::Element;
316
318 LocalFunction(const DiscreteGlobalBasisFunction& globalFunction)
319 : LocalBase(globalFunction.data_)
320 , evaluationBuffer_(this->localView_.tree())
321 {
322 /* Nothing. */
323 }
324
334 Range operator()(const Domain& x) const
335 {
336 Range y;
337 istlVectorBackend(y) = 0;
338
339 TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
340 const auto& fe = node.finiteElement();
341 const auto& localBasis = fe.localBasis();
342 auto& shapeFunctionValues = evaluationBuffer_[treePath];
343
344 localBasis.evaluateFunction(x, shapeFunctionValues);
345
346 // Compute linear combinations of basis function jacobian.
347 // Non-scalar coefficients of dimension coeffDim are handled by
348 // processing the coeffDim linear combinations independently
349 // and storing them as entries of an array.
350 using Value = LocalBasisRange< std::decay_t<decltype(node)> >;
351 static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
352 auto values = std::array<Value, coeffDim>{};
353 istlVectorBackend(values) = 0;
354 for (size_type i = 0; i < localBasis.size(); ++i)
355 {
356 auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
357 for (std::size_t j = 0; j < coeffDim; ++j)
358 values[j].axpy(c[j], shapeFunctionValues[i]);
359 }
360
361 // Assign computed values to node entry of range.
362 // Types are matched using the lexicographic ordering provided by flatVectorView.
363 LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), values);
364 });
365
366 return y;
367 }
368
371 {
373 if (lf.bound())
374 dlf.bind(lf.localContext());
375 return dlf;
376 }
377
378 private:
379 mutable PerNodeEvaluationBuffer evaluationBuffer_;
380 };
381
383 template<class B_T, class V_T, class NTRE_T>
384 DiscreteGlobalBasisFunction(B_T && basis, V_T && coefficients, NTRE_T&& nodeToRangeEntry)
385 : Base(std::make_shared<Data>(Data{{basis.gridView()}, wrap_or_move(std::forward<B_T>(basis)), wrap_or_move(std::forward<V_T>(coefficients)), wrap_or_move(std::forward<NTRE_T>(nodeToRangeEntry))}))
386 {}
387
389 DiscreteGlobalBasisFunction(std::shared_ptr<const Basis> basis, std::shared_ptr<const V> coefficients, std::shared_ptr<const typename Base::NodeToRangeEntry> nodeToRangeEntry)
390 : Base(std::make_shared<Data>(Data{{basis->gridView()}, basis, coefficients, nodeToRangeEntry}))
391 {}
392
398 Range operator() (const Domain& x) const
399 {
400 HierarchicSearch search(this->data_->basis->gridView().grid(), this->data_->basis->gridView().indexSet());
401
402 const auto e = search.findEntity(x);
403 auto localThis = localFunction(*this);
404 localThis.bind(e);
405 return localThis(e.geometry().local(x));
406 }
407
410 {
412 }
413
422 friend LocalFunction localFunction(const DiscreteGlobalBasisFunction& t)
423 {
424 return LocalFunction(t);
425 }
426};
427
428
453template<typename R, typename B, typename V>
454auto makeDiscreteGlobalBasisFunction(B&& basis, V&& vector)
455{
456 using Basis = std::decay_t<B>;
457 using NTREM = HierarchicNodeToRangeMap;
458
459 // Small helper functions to wrap vectors using istlVectorBackend
460 // if they do not already satisfy the VectorBackend interface.
461 auto toConstVectorBackend = [&](auto&& v) -> decltype(auto) {
462 if constexpr (models<Concept::ConstVectorBackend<Basis>, decltype(v)>()) {
463 return std::forward<decltype(v)>(v);
464 } else {
465 return istlVectorBackend(v);
466 }
467 };
468
469 using Vector = std::decay_t<decltype(toConstVectorBackend(std::forward<V>(vector)))>;
471 std::forward<B>(basis),
472 toConstVectorBackend(std::forward<V>(vector)),
474}
475
476
491template<typename DGBF>
493 : public ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>
494{
495 using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>;
496 using Data = typename Base::Data;
497
498public:
499 using DiscreteGlobalBasisFunction = DGBF;
500
501 using Basis = typename Base::Basis;
502 using Vector = typename Base::Vector;
503
504 using Domain = typename Base::Domain;
506
507 using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
508
509private:
510
511 template<class Node>
512 using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::JacobianType;
513 template<class Node>
514 using NodeData = typename std::vector< LocalBasisRange<Node> >;
515 using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
516
517public:
518
527 : public Base::LocalFunctionBase
528 {
529 using LocalBase = typename Base::LocalFunctionBase;
530 using size_type = typename Base::Tree::size_type;
531 using LocalBase::nodeToRangeEntry;
532
533 public:
535 using Domain = typename LocalBase::Domain;
536 using Range = GlobalFunction::Range;
537 using Element = typename LocalBase::Element;
538
540 LocalFunction(const GlobalFunction& globalFunction)
541 : LocalBase(globalFunction.data_)
542 , evaluationBuffer_(this->localView_.tree())
543 {
544 /* Nothing. */
545 }
546
553 void bind(const Element& element)
554 {
555 LocalBase::bind(element);
556 geometry_.emplace(element.geometry());
557 }
558
560 void unbind()
561 {
562 geometry_.reset();
563 LocalBase::unbind();
564 }
565
579 Range operator()(const Domain& x) const
580 {
581 Range y;
582 istlVectorBackend(y) = 0;
583
584 const auto& jacobianInverse = geometry_->jacobianInverse(x);
585
586 TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
587 const auto& fe = node.finiteElement();
588 const auto& localBasis = fe.localBasis();
589 auto& shapeFunctionJacobians = evaluationBuffer_[treePath];
590
591 localBasis.evaluateJacobian(x, shapeFunctionJacobians);
592
593 // Compute linear combinations of basis function jacobian.
594 // Non-scalar coefficients of dimension coeffDim are handled by
595 // processing the coeffDim linear combinations independently
596 // and storing them as entries of an array.
597 using RefJacobian = LocalBasisRange< std::decay_t<decltype(node)> >;
598 static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
599 auto refJacobians = std::array<RefJacobian, coeffDim>{};
600 istlVectorBackend(refJacobians) = 0;
601 for (size_type i = 0; i < localBasis.size(); ++i)
602 {
603 auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
604 for (std::size_t j = 0; j < coeffDim; ++j)
605 refJacobians[j].axpy(c[j], shapeFunctionJacobians[i]);
606 }
607
608 // Transform Jacobians form local to global coordinates.
609 using Jacobian = decltype(refJacobians[0] * jacobianInverse);
610 auto jacobians = std::array<Jacobian, coeffDim>{};
611 std::transform(
612 refJacobians.begin(), refJacobians.end(), jacobians.begin(),
613 [&](const auto& refJacobian) { return refJacobian * jacobianInverse; });
614
615 // Assign computed Jacobians to node entry of range.
616 // Types are matched using the lexicographic ordering provided by flatVectorView.
617 LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), jacobians);
618 });
619
620 return y;
621 }
622
624 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction&)
625 {
626 DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
627 }
628
629 private:
630 mutable PerNodeEvaluationBuffer evaluationBuffer_;
631 std::optional<typename Element::Geometry> geometry_;
632 };
633
640 DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr<const Data>& data)
641 : Base(data)
642 {
643 /* Nothing. */
644 }
645
651 Range operator()(const Domain& x) const
652 {
653 HierarchicSearch search(this->data_->basis->gridView().grid(), this->data_->basis->gridView().indexSet());
654
655 const auto e = search.findEntity(x);
656 auto localThis = localFunction(*this);
657 localThis.bind(e);
658 return localThis(e.geometry().local(x));
659 }
660
661 friend typename Traits::DerivativeInterface derivative(const DiscreteGlobalBasisFunctionDerivative& f)
662 {
663 DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
664 }
665
668 {
669 return LocalFunction(f);
670 }
671};
672
673
674} // namespace Functions
675} // namespace Dune
676
677#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
local function evaluating the derivative in reference coordinates
Definition: discreteglobalbasisfunction.hh:528
Range operator()(const Domain &x) const
Evaluate this local-function in coordinates x in the bound element.
Definition: discreteglobalbasisfunction.hh:579
friend Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction &)
Not implemented.
Definition: discreteglobalbasisfunction.hh:624
void unbind()
Unbind the local-function.
Definition: discreteglobalbasisfunction.hh:560
LocalFunction(const GlobalFunction &globalFunction)
Create a local function from the associated grid function.
Definition: discreteglobalbasisfunction.hh:540
void bind(const Element &element)
Bind LocalFunction to grid element.
Definition: discreteglobalbasisfunction.hh:553
Derivative of a DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:494
Range operator()(const Domain &x) const
Evaluate the discrete grid-function derivative in global coordinates.
Definition: discreteglobalbasisfunction.hh:651
friend LocalFunction localFunction(const DiscreteGlobalBasisFunctionDerivative &f)
Construct local function from a DiscreteGlobalBasisFunctionDerivative
Definition: discreteglobalbasisfunction.hh:667
DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr< const Data > &data)
create object from DiscreateGlobalBasisFunction data
Definition: discreteglobalbasisfunction.hh:640
A grid function induced by a global basis and a coefficient vector.
Definition: discreteglobalbasisfunction.hh:281
friend DiscreteGlobalBasisFunctionDerivative< DiscreteGlobalBasisFunction > derivative(const DiscreteGlobalBasisFunction &f)
Derivative of the DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:409
DiscreteGlobalBasisFunction(B_T &&basis, V_T &&coefficients, NTRE_T &&nodeToRangeEntry)
Create a grid-function, by wrapping the arguments in std::shared_ptr.
Definition: discreteglobalbasisfunction.hh:384
DiscreteGlobalBasisFunction(std::shared_ptr< const Basis > basis, std::shared_ptr< const V > coefficients, std::shared_ptr< const typename Base::NodeToRangeEntry > nodeToRangeEntry)
Create a grid-function, by moving the arguments in std::shared_ptr.
Definition: discreteglobalbasisfunction.hh:389
friend LocalFunction localFunction(const DiscreteGlobalBasisFunction &t)
Construct local function from a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:422
Range operator()(const Domain &x) const
Evaluate at a point given in world coordinates.
Definition: discreteglobalbasisfunction.hh:398
GridView::template Codim< codim >::Entity Element
Type of Elements contained in this EntitySet.
Definition: gridviewentityset.hh:36
Element::Geometry::LocalCoordinate LocalCoordinate
Type of local coordinates with respect to the Element.
Definition: gridviewentityset.hh:39
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:43
auto makeDiscreteGlobalBasisFunction(B &&basis, V &&vector)
Generate a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:454
auto istlVectorBackend(Vector &v)
Return a vector backend wrapping non-const ISTL like containers.
Definition: istlvectorbackend.hh:350
Definition: polynomial.hh:17
Default implementation for derivative traits.
Definition: defaultderivativetraits.hh:41
A simple node to range map using the nested tree indices.
Definition: hierarchicnodetorangemap.hh:34
Helper class to deduce the signature of a callable.
Definition: signature.hh:60
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Aug 13, 22:30, 2024)