DUNE PDELab (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
14
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 if (node.empty())
341 return;
342 const auto& fe = node.finiteElement();
343 const auto& localBasis = fe.localBasis();
344 auto& shapeFunctionValues = evaluationBuffer_[treePath];
345
346 localBasis.evaluateFunction(x, shapeFunctionValues);
347
348 // Compute linear combinations of basis function jacobian.
349 // Non-scalar coefficients of dimension coeffDim are handled by
350 // processing the coeffDim linear combinations independently
351 // and storing them as entries of an array.
352 using Value = LocalBasisRange< std::decay_t<decltype(node)> >;
353 static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
354 auto values = std::array<Value, coeffDim>{};
355 istlVectorBackend(values) = 0;
356 for (size_type i = 0; i < localBasis.size(); ++i)
357 {
358 auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
359 for (std::size_t j = 0; j < coeffDim; ++j)
360 values[j].axpy(c[j], shapeFunctionValues[i]);
361 }
362
363 // Assign computed values to node entry of range.
364 // Types are matched using the lexicographic ordering provided by flatVectorView.
365 LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), values);
366 });
367
368 return y;
369 }
370
373 {
375 if (lf.bound())
376 dlf.bind(lf.localContext());
377 return dlf;
378 }
379
380 private:
381 mutable PerNodeEvaluationBuffer evaluationBuffer_;
382 };
383
385 template<class B_T, class V_T, class NTRE_T>
386 DiscreteGlobalBasisFunction(B_T && basis, V_T && coefficients, NTRE_T&& nodeToRangeEntry)
387 : 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))}))
388 {}
389
391 DiscreteGlobalBasisFunction(std::shared_ptr<const Basis> basis, std::shared_ptr<const V> coefficients, std::shared_ptr<const typename Base::NodeToRangeEntry> nodeToRangeEntry)
392 : Base(std::make_shared<Data>(Data{{basis->gridView()}, basis, coefficients, nodeToRangeEntry}))
393 {}
394
400 Range operator() (const Domain& x) const
401 {
402 HierarchicSearch search(this->data_->basis->gridView().grid(), this->data_->basis->gridView().indexSet());
403
404 const auto e = search.findEntity(x);
405 auto localThis = localFunction(*this);
406 localThis.bind(e);
407 return localThis(e.geometry().local(x));
408 }
409
412 {
414 }
415
424 friend LocalFunction localFunction(const DiscreteGlobalBasisFunction& t)
425 {
426 return LocalFunction(t);
427 }
428};
429
430
455template<typename R, typename B, typename V>
456auto makeDiscreteGlobalBasisFunction(B&& basis, V&& vector)
457{
458 using Basis = std::decay_t<B>;
459 using NTREM = HierarchicNodeToRangeMap;
460
461 // Small helper functions to wrap vectors using istlVectorBackend
462 // if they do not already satisfy the VectorBackend interface.
463 auto toConstVectorBackend = [&](auto&& v) -> decltype(auto) {
464 if constexpr (models<Concept::ConstVectorBackend<Basis>, decltype(v)>()) {
465 return std::forward<decltype(v)>(v);
466 } else {
467 return istlVectorBackend(v);
468 }
469 };
470
471 using Vector = std::decay_t<decltype(toConstVectorBackend(std::forward<V>(vector)))>;
473 std::forward<B>(basis),
474 toConstVectorBackend(std::forward<V>(vector)),
476}
477
478
493template<typename DGBF>
495 : public ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>
496{
497 using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>;
498 using Data = typename Base::Data;
499
500public:
501 using DiscreteGlobalBasisFunction = DGBF;
502
503 using Basis = typename Base::Basis;
504 using Vector = typename Base::Vector;
505
506 using Domain = typename Base::Domain;
508
509 using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
510
511private:
512
513 template<class Node>
514 using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::JacobianType;
515 template<class Node>
516 using NodeData = typename std::vector< LocalBasisRange<Node> >;
517 using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
518
519public:
520
529 : public Base::LocalFunctionBase
530 {
531 using LocalBase = typename Base::LocalFunctionBase;
532 using size_type = typename Base::Tree::size_type;
533 using LocalBase::nodeToRangeEntry;
534
535 public:
537 using Domain = typename LocalBase::Domain;
538 using Range = GlobalFunction::Range;
539 using Element = typename LocalBase::Element;
540
542 LocalFunction(const GlobalFunction& globalFunction)
543 : LocalBase(globalFunction.data_)
544 , evaluationBuffer_(this->localView_.tree())
545 {
546 /* Nothing. */
547 }
548
555 void bind(const Element& element)
556 {
557 LocalBase::bind(element);
558 geometry_.emplace(element.geometry());
559 }
560
562 void unbind()
563 {
564 geometry_.reset();
565 LocalBase::unbind();
566 }
567
581 Range operator()(const Domain& x) const
582 {
583 Range y;
584 istlVectorBackend(y) = 0;
585
586 const auto& jacobianInverse = geometry_->jacobianInverse(x);
587
588 TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
589 if (node.empty())
590 return;
591 const auto& fe = node.finiteElement();
592 const auto& localBasis = fe.localBasis();
593 auto& shapeFunctionJacobians = evaluationBuffer_[treePath];
594
595 localBasis.evaluateJacobian(x, shapeFunctionJacobians);
596
597 // Compute linear combinations of basis function jacobian.
598 // Non-scalar coefficients of dimension coeffDim are handled by
599 // processing the coeffDim linear combinations independently
600 // and storing them as entries of an array.
601 using RefJacobian = LocalBasisRange< std::decay_t<decltype(node)> >;
602 static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
603 auto refJacobians = std::array<RefJacobian, coeffDim>{};
604 istlVectorBackend(refJacobians) = 0;
605 for (size_type i = 0; i < localBasis.size(); ++i)
606 {
607 auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
608 for (std::size_t j = 0; j < coeffDim; ++j)
609 refJacobians[j].axpy(c[j], shapeFunctionJacobians[i]);
610 }
611
612 // Transform Jacobians form local to global coordinates.
613 using Jacobian = decltype(refJacobians[0] * jacobianInverse);
614 auto jacobians = std::array<Jacobian, coeffDim>{};
615 std::transform(
616 refJacobians.begin(), refJacobians.end(), jacobians.begin(),
617 [&](const auto& refJacobian) { return refJacobian * jacobianInverse; });
618
619 // Assign computed Jacobians to node entry of range.
620 // Types are matched using the lexicographic ordering provided by flatVectorView.
621 LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), jacobians);
622 });
623
624 return y;
625 }
626
628 friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction&)
629 {
630 DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
631 }
632
633 private:
634 mutable PerNodeEvaluationBuffer evaluationBuffer_;
635 std::optional<typename Element::Geometry> geometry_;
636 };
637
644 DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr<const Data>& data)
645 : Base(data)
646 {
647 /* Nothing. */
648 }
649
655 Range operator()(const Domain& x) const
656 {
657 HierarchicSearch search(this->data_->basis->gridView().grid(), this->data_->basis->gridView().indexSet());
658
659 const auto e = search.findEntity(x);
660 auto localThis = localFunction(*this);
661 localThis.bind(e);
662 return localThis(e.geometry().local(x));
663 }
664
665 friend typename Traits::DerivativeInterface derivative(const DiscreteGlobalBasisFunctionDerivative& f)
666 {
667 DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
668 }
669
672 {
673 return LocalFunction(f);
674 }
675};
676
677
678} // namespace Functions
679} // namespace Dune
680
681#endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
local function evaluating the derivative in reference coordinates
Definition: discreteglobalbasisfunction.hh:530
Range operator()(const Domain &x) const
Evaluate this local-function in coordinates x in the bound element.
Definition: discreteglobalbasisfunction.hh:581
friend Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction &)
Not implemented.
Definition: discreteglobalbasisfunction.hh:628
void unbind()
Unbind the local-function.
Definition: discreteglobalbasisfunction.hh:562
LocalFunction(const GlobalFunction &globalFunction)
Create a local function from the associated grid function.
Definition: discreteglobalbasisfunction.hh:542
void bind(const Element &element)
Bind LocalFunction to grid element.
Definition: discreteglobalbasisfunction.hh:555
Derivative of a DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:496
Range operator()(const Domain &x) const
Evaluate the discrete grid-function derivative in global coordinates.
Definition: discreteglobalbasisfunction.hh:655
friend LocalFunction localFunction(const DiscreteGlobalBasisFunctionDerivative &f)
Construct local function from a DiscreteGlobalBasisFunctionDerivative
Definition: discreteglobalbasisfunction.hh:671
DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr< const Data > &data)
create object from DiscreateGlobalBasisFunction data
Definition: discreteglobalbasisfunction.hh:644
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:411
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:386
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:391
friend LocalFunction localFunction(const DiscreteGlobalBasisFunction &t)
Construct local function from a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:424
Range operator()(const Domain &x) const
Evaluate at a point given in world coordinates.
Definition: discreteglobalbasisfunction.hh:400
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
Search an IndexSet for an Entity containing a given point.
Definition: hierarchicsearch.hh:35
Entity findEntity(const FieldVector< ct, dimw > &global) const
Search the IndexSet of this HierarchicSearch for an Entity containing point global.
Definition: hierarchicsearch.hh:128
Default exception for dummy implementations.
Definition: exceptions.hh:357
Traits for type conversions and type information.
typename AutonomousValueType< T >::type AutonomousValue
Type free of internal references that T can be converted to.
Definition: typetraits.hh:589
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
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:456
auto istlVectorBackend(Vector &v)
Return a vector backend wrapping non-const ISTL like containers.
Definition: istlvectorbackend.hh:350
constexpr auto treePath(const T &... t)
Constructs a new HybridTreePath from the given indices.
Definition: treepath.hh:199
void forEachLeafNode(Tree &&tree, LeafFunc &&leafFunc)
Traverse tree and visit each leaf node.
Definition: traversal.hh:289
std::decay_t< decltype(makeTreeContainer(std::declval< const Tree & >(), std::declval< Detail::LeafToDefaultConstructibleValue< LeafToValue > >()))> TreeContainer
Alias to container type generated by makeTreeContainer for give tree type and when using LeafToValue ...
Definition: treecontainer.hh:306
Utility class for hierarchically searching for an Entity containing a given point.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
auto wrap_or_move(T &&t)
Capture R-value reference to shared_ptr.
Definition: shared_ptr.hh:96
STL namespace.
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 (Sep 4, 22:38, 2025)