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#include <dune/common/shared_ptr.hh>
15
16#include <dune/grid/utility/hierarchicsearch.hh>
17
18#include <dune/common/typetree/traversal.hh>
19#include <dune/common/typetree/treecontainer.hh>
20
21#include <dune/functions/functionspacebases/hierarchicnodetorangemap.hh>
22#include <dune/functions/functionspacebases/flatvectorview.hh>
23#include <dune/functions/gridfunctions/gridviewentityset.hh>
24#include <dune/functions/gridfunctions/gridfunction.hh>
25#include <dune/functions/backends/concepts.hh>
26#include <dune/functions/backends/istlvectorbackend.hh>
27
28namespace Dune {
29namespace Functions {
30
31
32namespace ImplDoc {
33
34template<typename B, typename V, typename NTRE>
35class DiscreteGlobalBasisFunctionBase
36{
37public:
38 using Basis = B;
39 using Vector = V;
40
41 // In order to make the cache work for proxy-references
42 // we have to use AutonomousValue<T> instead of std::decay_t<T>
43 using Coefficient = Dune::AutonomousValue<decltype(std::declval<Vector>()[std::declval<typename Basis::MultiIndex>()])>;
44
45 using GridView = typename Basis::GridView;
46 using EntitySet = GridViewEntitySet<GridView, 0>;
47 using Tree = typename Basis::LocalView::Tree;
48 using NodeToRangeEntry = NTRE;
49
50 using Domain = typename EntitySet::GlobalCoordinate;
51
52 using LocalDomain = typename EntitySet::LocalCoordinate;
53 using Element = typename EntitySet::Element;
54
55protected:
56
57 // This collects all data that is shared by all related
58 // global and local functions. This way we don't need to
59 // keep track of it individually.
60 struct Data
61 {
62 EntitySet entitySet;
63 std::shared_ptr<const Basis> basis;
64 std::shared_ptr<const Vector> coefficients;
65 std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry;
66 };
67
68public:
69 class LocalFunctionBase
70 {
71 using LocalView = typename Basis::LocalView;
72 using size_type = typename Tree::size_type;
73
74 public:
75 using Domain = LocalDomain;
76 using Element = typename EntitySet::Element;
77
78 protected:
79 LocalFunctionBase(const std::shared_ptr<const Data>& data)
80 : data_(data)
81 , localView_(data_->basis->localView())
82 {
83 localDoFs_.reserve(localView_.maxSize());
84 }
85
92 LocalFunctionBase(const LocalFunctionBase& other)
93 : data_(other.data_)
94 , localView_(other.localView_)
95 {
96 localDoFs_.reserve(localView_.maxSize());
97 if (bound())
98 localDoFs_ = other.localDoFs_;
99 }
100
108 LocalFunctionBase& operator=(const LocalFunctionBase& other)
109 {
110 data_ = other.data_;
111 localView_ = other.localView_;
112 if (bound())
113 localDoFs_ = other.localDoFs_;
114 return *this;
115 }
116
117 public:
124 void bind(const Element& element)
125 {
126 localView_.bind(element);
127 // Use cache of full local view size. For a subspace basis,
128 // this may be larger than the number of local DOFs in the
129 // tree. In this case only cache entries associated to local
130 // DOFs in the subspace are filled. Cache entries associated
131 // to local DOFs which are not contained in the subspace will
132 // not be touched.
133 //
134 // Alternatively one could use a cache that exactly fits
135 // the size of the tree. However, this would require to
136 // subtract an offset from localIndex(i) on each cache
137 // access in operator().
138 localDoFs_.resize(localView_.size());
139 const auto& dofs = *data_->coefficients;
140 for (size_type i = 0; i < localView_.tree().size(); ++i)
141 {
142 // For a subspace basis the index-within-tree i
143 // is not the same as the localIndex within the
144 // full local view.
145 size_t localIndex = localView_.tree().localIndex(i);
146 localDoFs_[localIndex] = dofs[localView_.index(localIndex)];
147 }
148 }
149
151 void unbind()
152 {
153 localView_.unbind();
154 }
155
157 bool bound() const
158 {
159 return localView_.bound();
160 }
161
163 const Element& localContext() const
164 {
165 return localView_.element();
166 }
167
168 protected:
169
170 template<class To, class From>
171 void assignWith(To& to, const From& from) const
172 {
173 auto from_flat = flatVectorView(from);
174 auto to_flat = flatVectorView(to);
175 assert(from_flat.size() == to_flat.size());
176 for (size_type i = 0; i < to_flat.size(); ++i)
177 to_flat[i] = from_flat[i];
178 }
179
180 template<class Node, class TreePath, class Range>
181 decltype(auto) nodeToRangeEntry(const Node& node, const TreePath& treePath, Range& y) const
182 {
183 return (*data_->nodeToRangeEntry)(node, treePath, y);
184 }
185
186 std::shared_ptr<const Data> data_;
187 LocalView localView_;
188 std::vector<Coefficient> localDoFs_;
189 };
190
191protected:
192 DiscreteGlobalBasisFunctionBase(const std::shared_ptr<const Data>& data)
193 : data_(data)
194 {
195 /* Nothing. */
196 }
197
198public:
199
201 const Basis& basis() const
202 {
203 return *data_->basis;
204 }
205
207 const Vector& dofs() const
208 {
209 return *data_->coefficients;
210 }
211
213 const NodeToRangeEntry& nodeToRangeEntry() const
214 {
215 return *data_->nodeToRangeEntry;
216 }
217
219 const EntitySet& entitySet() const
220 {
221 return data_->entitySet;
222 }
223
224protected:
225 std::shared_ptr<const Data> data_;
226};
227
228} // namespace ImplDoc
229
230
231
232template<typename DGBF>
233class DiscreteGlobalBasisFunctionDerivative;
234
278template<typename B, typename V,
279 typename NTRE = HierarchicNodeToRangeMap,
280 typename R = typename V::value_type>
282 : public ImplDoc::DiscreteGlobalBasisFunctionBase<B, V, NTRE>
283{
284 using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<B, V, NTRE>;
285 using Data = typename Base::Data;
286
287public:
288 using Basis = typename Base::Basis;
289 using Vector = typename Base::Vector;
290
291 using Domain = typename Base::Domain;
292 using Range = R;
293
294 using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
295
296private:
297
298 template<class Node>
299 using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
300 template<class Node>
301 using NodeData = typename std::vector<LocalBasisRange<Node>>;
302 using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
303
304public:
305 class LocalFunction
306 : public Base::LocalFunctionBase
307 {
308 using LocalBase = typename Base::LocalFunctionBase;
309 using size_type = typename Base::Tree::size_type;
310 using LocalBase::nodeToRangeEntry;
311
312 public:
313
314 using GlobalFunction = DiscreteGlobalBasisFunction;
315 using Domain = typename LocalBase::Domain;
316 using Range = GlobalFunction::Range;
317 using Element = typename LocalBase::Element;
318
320 LocalFunction(const DiscreteGlobalBasisFunction& globalFunction)
321 : LocalBase(globalFunction.data_)
322 , evaluationBuffer_(this->localView_.tree())
323 {
324 /* Nothing. */
325 }
326
336 Range operator()(const Domain& x) const
337 {
338 Range y;
339 istlVectorBackend(y) = 0;
340
341 TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
342 if (node.empty())
343 return;
344 const auto& fe = node.finiteElement();
345 const auto& localBasis = fe.localBasis();
346 auto& shapeFunctionValues = evaluationBuffer_[treePath];
347
348 localBasis.evaluateFunction(x, shapeFunctionValues);
349
350 // Compute linear combinations of basis function jacobian.
351 // Non-scalar coefficients of dimension coeffDim are handled by
352 // processing the coeffDim linear combinations independently
353 // and storing them as entries of an array.
354 using Value = LocalBasisRange< std::decay_t<decltype(node)> >;
355 static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
356 auto values = std::array<Value, coeffDim>{};
357 istlVectorBackend(values) = 0;
358 for (size_type i = 0; i < localBasis.size(); ++i)
359 {
360 auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
361 for (std::size_t j = 0; j < coeffDim; ++j)
362 values[j].axpy(c[j], shapeFunctionValues[i]);
363 }
364
365 // Assign computed values to node entry of range.
366 // Types are matched using the lexicographic ordering provided by flatVectorView.
367 LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), values);
368 });
369
370 return y;
371 }
372
375 {
377 if (lf.bound())
378 dlf.bind(lf.localContext());
379 return dlf;
380 }
381
382 private:
383 mutable PerNodeEvaluationBuffer evaluationBuffer_;
384 };
385
387 template<class B_T, class V_T, class NTRE_T>
388 DiscreteGlobalBasisFunction(B_T && basis, V_T && coefficients, NTRE_T&& nodeToRangeEntry)
389 : 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))}))
390 {}
391
393 DiscreteGlobalBasisFunction(std::shared_ptr<const Basis> basis, std::shared_ptr<const V> coefficients, std::shared_ptr<const typename Base::NodeToRangeEntry> nodeToRangeEntry)
394 : Base(std::make_shared<Data>(Data{{basis->gridView()}, basis, coefficients, nodeToRangeEntry}))
395 {}
396
402 Range operator() (const Domain& x) const
403 {
404 HierarchicSearch search(this->data_->basis->gridView().grid(), this->data_->basis->gridView().indexSet());
405
406 const auto e = search.findEntity(x);
407 auto localThis = localFunction(*this);
408 localThis.bind(e);
409 return localThis(e.geometry().local(x));
410 }
411
414 {
416 }
417
426 friend LocalFunction localFunction(const DiscreteGlobalBasisFunction& t)
427 {
428 return LocalFunction(t);
429 }
430};
431
432
457template<typename R, typename B, typename V>
458auto makeDiscreteGlobalBasisFunction(B&& basis, V&& vector)
459{
460 using Basis = std::decay_t<B>;
461 using NTREM = HierarchicNodeToRangeMap;
462
463 // Small helper functions to wrap vectors using istlVectorBackend
464 // if they do not already satisfy the VectorBackend interface.
465 auto toConstVectorBackend = [&](auto&& v) -> decltype(auto) {
466 if constexpr (models<Concept::ConstVectorBackend<Basis>, decltype(v)>()) {
467 return std::forward<decltype(v)>(v);
468 } else {
469 return istlVectorBackend(v);
470 }
471 };
472
473 using Vector = std::decay_t<decltype(toConstVectorBackend(std::forward<V>(vector)))>;
475 std::forward<B>(basis),
476 toConstVectorBackend(std::forward<V>(vector)),
478}
479
480
495template<typename DGBF>
497 : public ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>
498{
499 using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>;
500 using Data = typename Base::Data;
501
502public:
503 using DiscreteGlobalBasisFunction = DGBF;
504
505 using Basis = typename Base::Basis;
506 using Vector = typename Base::Vector;
507
508 using Domain = typename Base::Domain;
510
511 using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
512
513private:
514
515 template<class Node>
516 using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::JacobianType;
517 template<class Node>
518 using NodeData = typename std::vector< LocalBasisRange<Node> >;
519 using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
520
521public:
522
531 : public Base::LocalFunctionBase
532 {
533 using LocalBase = typename Base::LocalFunctionBase;
534 using size_type = typename Base::Tree::size_type;
535 using LocalBase::nodeToRangeEntry;
536
537 public:
539 using Domain = typename LocalBase::Domain;
540 using Range = GlobalFunction::Range;
541 using Element = typename LocalBase::Element;
542
544 LocalFunction(const GlobalFunction& globalFunction)
545 : LocalBase(globalFunction.data_)
546 , evaluationBuffer_(this->localView_.tree())
547 {
548 /* Nothing. */
549 }
550
557 void bind(const Element& element)
558 {
559 LocalBase::bind(element);
560 geometry_.emplace(element.geometry());
561 }
562
564 void unbind()
565 {
566 geometry_.reset();
567 LocalBase::unbind();
568 }
569
583 Range operator()(const Domain& x) const
584 {
585 Range y;
586 istlVectorBackend(y) = 0;
587
588 const auto& jacobianInverse = geometry_->jacobianInverse(x);
589
590 TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
591 if (node.empty())
592 return;
593 const auto& fe = node.finiteElement();
594 const auto& localBasis = fe.localBasis();
595 auto& shapeFunctionJacobians = evaluationBuffer_[treePath];
596
597 localBasis.evaluateJacobian(x, shapeFunctionJacobians);
598
599 // Compute linear combinations of basis function jacobian.
600 // Non-scalar coefficients of dimension coeffDim are handled by
601 // processing the coeffDim linear combinations independently
602 // and storing them as entries of an array.
603 using RefJacobian = LocalBasisRange< std::decay_t<decltype(node)> >;
604 static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
605 auto refJacobians = std::array<RefJacobian, coeffDim>{};
606 istlVectorBackend(refJacobians) = 0;
607 for (size_type i = 0; i < localBasis.size(); ++i)
608 {
609 auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
610 for (std::size_t j = 0; j < coeffDim; ++j)
611 refJacobians[j].axpy(c[j], shapeFunctionJacobians[i]);
612 }
613
614 // Transform Jacobians form local to global coordinates.
615 using Jacobian = decltype(refJacobians[0] * jacobianInverse);
616 auto jacobians = std::array<Jacobian, coeffDim>{};
617 std::transform(
618 refJacobians.begin(), refJacobians.end(), jacobians.begin(),
619 [&](const auto& refJacobian) { return refJacobian * jacobianInverse; });
620
621 // Assign computed Jacobians to node entry of range.
622 // Types are matched using the lexicographic ordering provided by flatVectorView.
623 LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), jacobians);
624 });
625
626 return y;
627 }
628
630 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction&)
631 {
632 DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
633 }
634
635 private:
636 mutable PerNodeEvaluationBuffer evaluationBuffer_;
637 std::optional<typename Element::Geometry> geometry_;
638 };
639
646 DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr<const Data>& data)
647 : Base(data)
648 {
649 /* Nothing. */
650 }
651
657 Range operator()(const Domain& x) const
658 {
659 HierarchicSearch search(this->data_->basis->gridView().grid(), this->data_->basis->gridView().indexSet());
660
661 const auto e = search.findEntity(x);
662 auto localThis = localFunction(*this);
663 localThis.bind(e);
664 return localThis(e.geometry().local(x));
665 }
666
667 friend typename Traits::DerivativeInterface derivative(const DiscreteGlobalBasisFunctionDerivative& f)
668 {
669 DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
670 }
671
674 {
675 return LocalFunction(f);
676 }
677};
678
679
680} // namespace Functions
681} // namespace Dune
682
683#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
local function evaluating the derivative in reference coordinates
Definition: discreteglobalbasisfunction.hh:532
Range operator()(const Domain &x) const
Evaluate this local-function in coordinates x in the bound element.
Definition: discreteglobalbasisfunction.hh:583
friend Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction &)
Not implemented.
Definition: discreteglobalbasisfunction.hh:630
void unbind()
Unbind the local-function.
Definition: discreteglobalbasisfunction.hh:564
LocalFunction(const GlobalFunction &globalFunction)
Create a local function from the associated grid function.
Definition: discreteglobalbasisfunction.hh:544
void bind(const Element &element)
Bind LocalFunction to grid element.
Definition: discreteglobalbasisfunction.hh:557
Derivative of a DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:498
Range operator()(const Domain &x) const
Evaluate the discrete grid-function derivative in global coordinates.
Definition: discreteglobalbasisfunction.hh:657
friend LocalFunction localFunction(const DiscreteGlobalBasisFunctionDerivative &f)
Construct local function from a DiscreteGlobalBasisFunctionDerivative
Definition: discreteglobalbasisfunction.hh:673
DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr< const Data > &data)
create object from DiscreateGlobalBasisFunction data
Definition: discreteglobalbasisfunction.hh:646
A grid function induced by a global basis and a coefficient vector.
Definition: discreteglobalbasisfunction.hh:283
friend DiscreteGlobalBasisFunctionDerivative< DiscreteGlobalBasisFunction > derivative(const DiscreteGlobalBasisFunction &f)
Derivative of the DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:413
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:388
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:393
friend LocalFunction localFunction(const DiscreteGlobalBasisFunction &t)
Construct local function from a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:426
Range operator()(const Domain &x) const
Evaluate at a point given in world coordinates.
Definition: discreteglobalbasisfunction.hh:402
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
std::vector< Child > Vector
Descriptor for vectors with all children of the same type and dynamic size.
Definition: containerdescriptors.hh:112
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:458
auto istlVectorBackend(Vector &v)
Return a vector backend wrapping non-const ISTL like containers.
Definition: istlvectorbackend.hh:350
Definition: monomialset.hh:19
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 & Uni Heidelberg  |  generated with Hugo v0.111.3 (Jan 9, 23:34, 2026)