Dune Core Modules (2.9.0)

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 #include <optional>
8 
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 
20 namespace Dune {
21 namespace Functions {
22 
23 
24 namespace ImplDoc {
25 
26 template<typename B, typename V, typename NTRE>
27 class DiscreteGlobalBasisFunctionBase
28 {
29 public:
30  using Basis = B;
31  using Vector = V;
32 
33  // In order to make the cache work for proxy-references
34  // we have to use AutonomousValue<T> instead of std::decay_t<T>
35  using Coefficient = Dune::AutonomousValue<decltype(std::declval<Vector>()[std::declval<typename Basis::MultiIndex>()])>;
36 
37  using GridView = typename Basis::GridView;
38  using EntitySet = GridViewEntitySet<GridView, 0>;
39  using Tree = typename Basis::LocalView::Tree;
40  using NodeToRangeEntry = NTRE;
41 
42  using Domain = typename EntitySet::GlobalCoordinate;
43 
44  using LocalDomain = typename EntitySet::LocalCoordinate;
45  using Element = typename EntitySet::Element;
46 
47 protected:
48 
49  // This collects all data that is shared by all related
50  // global and local functions. This way we don't need to
51  // keep track of it individually.
52  struct Data
53  {
54  EntitySet entitySet;
55  std::shared_ptr<const Basis> basis;
56  std::shared_ptr<const Vector> coefficients;
57  std::shared_ptr<const NodeToRangeEntry> nodeToRangeEntry;
58  };
59 
60 public:
61  class LocalFunctionBase
62  {
63  using LocalView = typename Basis::LocalView;
64  using size_type = typename Tree::size_type;
65 
66  public:
67  using Domain = LocalDomain;
68  using Element = typename EntitySet::Element;
69 
70  protected:
71  LocalFunctionBase(const std::shared_ptr<const Data>& data)
72  : data_(data)
73  , localView_(data_->basis->localView())
74  {
75  localDoFs_.reserve(localView_.maxSize());
76  }
77 
84  LocalFunctionBase(const LocalFunctionBase& other)
85  : data_(other.data_)
86  , localView_(other.localView_)
87  {
88  localDoFs_.reserve(localView_.maxSize());
89  if (bound())
90  localDoFs_ = other.localDoFs_;
91  }
92 
100  LocalFunctionBase& operator=(const LocalFunctionBase& other)
101  {
102  data_ = other.data_;
103  localView_ = other.localView_;
104  if (bound())
105  localDoFs_ = other.localDoFs_;
106  return *this;
107  }
108 
109  public:
116  void bind(const Element& element)
117  {
118  localView_.bind(element);
119  // Use cache of full local view size. For a subspace basis,
120  // this may be larger than the number of local DOFs in the
121  // tree. In this case only cache entries associated to local
122  // DOFs in the subspace are filled. Cache entries associated
123  // to local DOFs which are not contained in the subspace will
124  // not be touched.
125  //
126  // Alternatively one could use a cache that exactly fits
127  // the size of the tree. However, this would require to
128  // subtract an offset from localIndex(i) on each cache
129  // access in operator().
130  localDoFs_.resize(localView_.size());
131  const auto& dofs = *data_->coefficients;
132  for (size_type i = 0; i < localView_.tree().size(); ++i)
133  {
134  // For a subspace basis the index-within-tree i
135  // is not the same as the localIndex within the
136  // full local view.
137  size_t localIndex = localView_.tree().localIndex(i);
138  localDoFs_[localIndex] = dofs[localView_.index(localIndex)];
139  }
140  }
141 
143  void unbind()
144  {
145  localView_.unbind();
146  }
147 
149  bool bound() const
150  {
151  return localView_.bound();
152  }
153 
155  const Element& localContext() const
156  {
157  return localView_.element();
158  }
159 
160  protected:
161 
162  template<class To, class From>
163  void assignWith(To& to, const From& from) const
164  {
165  auto from_flat = flatVectorView(from);
166  auto to_flat = flatVectorView(to);
167  assert(from_flat.size() == to_flat.size());
168  for (size_type i = 0; i < to_flat.size(); ++i)
169  to_flat[i] = from_flat[i];
170  }
171 
172  template<class Node, class TreePath, class Range>
173  decltype(auto) nodeToRangeEntry(const Node& node, const TreePath& treePath, Range& y) const
174  {
175  return (*data_->nodeToRangeEntry)(node, treePath, y);
176  }
177 
178  std::shared_ptr<const Data> data_;
179  LocalView localView_;
180  std::vector<Coefficient> localDoFs_;
181  };
182 
183 protected:
184  DiscreteGlobalBasisFunctionBase(const std::shared_ptr<const Data>& data)
185  : data_(data)
186  {
187  /* Nothing. */
188  }
189 
190 public:
191 
193  const Basis& basis() const
194  {
195  return *data_->basis;
196  }
197 
199  const Vector& dofs() const
200  {
201  return *data_->coefficients;
202  }
203 
205  const NodeToRangeEntry& nodeToRangeEntry() const
206  {
207  return *data_->nodeToRangeEntry;
208  }
209 
211  const EntitySet& entitySet() const
212  {
213  return data_->entitySet;
214  }
215 
216 protected:
217  std::shared_ptr<const Data> data_;
218 };
219 
220 } // namespace ImplDoc
221 
222 
223 
224 template<typename DGBF>
225 class DiscreteGlobalBasisFunctionDerivative;
226 
264 template<typename B, typename V,
265  typename NTRE = HierarchicNodeToRangeMap,
266  typename R = typename V::value_type>
268  : public ImplDoc::DiscreteGlobalBasisFunctionBase<B, V, NTRE>
269 {
270  using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<B, V, NTRE>;
271  using Data = typename Base::Data;
272 
273 public:
274  using Basis = typename Base::Basis;
275  using Vector = typename Base::Vector;
276 
277  using Domain = typename Base::Domain;
278  using Range = R;
279 
280  using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
281 
282 private:
283 
284  template<class Node>
285  using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::RangeType;
286  template<class Node>
287  using NodeData = typename std::vector<LocalBasisRange<Node>>;
288  using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
289 
290 public:
291  class LocalFunction
292  : public Base::LocalFunctionBase
293  {
294  using LocalBase = typename Base::LocalFunctionBase;
295  using size_type = typename Base::Tree::size_type;
296  using LocalBase::nodeToRangeEntry;
297 
298  public:
299 
300  using GlobalFunction = DiscreteGlobalBasisFunction;
301  using Domain = typename LocalBase::Domain;
302  using Range = GlobalFunction::Range;
303  using Element = typename LocalBase::Element;
304 
306  LocalFunction(const DiscreteGlobalBasisFunction& globalFunction)
307  : LocalBase(globalFunction.data_)
308  , evaluationBuffer_(this->localView_.tree())
309  {
310  /* Nothing. */
311  }
312 
322  Range operator()(const Domain& x) const
323  {
324  Range y;
325  istlVectorBackend(y) = 0;
326 
327  TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
328  const auto& fe = node.finiteElement();
329  const auto& localBasis = fe.localBasis();
330  auto& shapeFunctionValues = evaluationBuffer_[treePath];
331 
332  localBasis.evaluateFunction(x, shapeFunctionValues);
333 
334  // Compute linear combinations of basis function jacobian.
335  // Non-scalar coefficients of dimension coeffDim are handled by
336  // processing the coeffDim linear combinations independently
337  // and storing them as entries of an array.
338  using Value = LocalBasisRange< std::decay_t<decltype(node)> >;
339  static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
340  auto values = std::array<Value, coeffDim>{};
341  istlVectorBackend(values) = 0;
342  for (size_type i = 0; i < localBasis.size(); ++i)
343  {
344  auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
345  for (std::size_t j = 0; j < coeffDim; ++j)
346  values[j].axpy(c[j], shapeFunctionValues[i]);
347  }
348 
349  // Assign computed values to node entry of range.
350  // Types are matched using the lexicographic ordering provided by flatVectorView.
351  LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), values);
352  });
353 
354  return y;
355  }
356 
359  {
361  if (lf.bound())
362  dlf.bind(lf.localContext());
363  return dlf;
364  }
365 
366  private:
367  mutable PerNodeEvaluationBuffer evaluationBuffer_;
368  };
369 
371  template<class B_T, class V_T, class NTRE_T>
372  DiscreteGlobalBasisFunction(B_T && basis, V_T && coefficients, NTRE_T&& nodeToRangeEntry)
373  : 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))}))
374  {}
375 
377  DiscreteGlobalBasisFunction(std::shared_ptr<const Basis> basis, std::shared_ptr<const V> coefficients, std::shared_ptr<const typename Base::NodeToRangeEntry> nodeToRangeEntry)
378  : Base(std::make_shared<Data>(Data{{basis->gridView()}, basis, coefficients, nodeToRangeEntry}))
379  {}
380 
382  Range operator() (const Domain& x) const
383  {
384  // TODO: Implement this using hierarchic search
385  DUNE_THROW(NotImplemented,"not implemented");
386  }
387 
390  {
392  }
393 
402  friend LocalFunction localFunction(const DiscreteGlobalBasisFunction& t)
403  {
404  return LocalFunction(t);
405  }
406 };
407 
408 
431 template<typename R, typename B, typename V>
432 auto makeDiscreteGlobalBasisFunction(B&& basis, V&& vector)
433 {
434  using Basis = std::decay_t<B>;
435  using NTREM = HierarchicNodeToRangeMap;
436 
437  // Small helper functions to wrap vectors using istlVectorBackend
438  // if they do not already satisfy the VectorBackend interface.
439  auto toConstVectorBackend = [&](auto&& v) -> decltype(auto) {
440  if constexpr (models<Concept::ConstVectorBackend<Basis>, decltype(v)>()) {
441  return std::forward<decltype(v)>(v);
442  } else {
443  return istlVectorBackend(v);
444  }
445  };
446 
447  using Vector = std::decay_t<decltype(toConstVectorBackend(std::forward<V>(vector)))>;
449  std::forward<B>(basis),
450  toConstVectorBackend(std::forward<V>(vector)),
452 }
453 
454 
469 template<typename DGBF>
471  : public ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>
472 {
473  using Base = ImplDoc::DiscreteGlobalBasisFunctionBase<typename DGBF::Basis, typename DGBF::Vector, typename DGBF::NodeToRangeEntry>;
474  using Data = typename Base::Data;
475 
476 public:
477  using DiscreteGlobalBasisFunction = DGBF;
478 
479  using Basis = typename Base::Basis;
480  using Vector = typename Base::Vector;
481 
482  using Domain = typename Base::Domain;
484 
485  using Traits = Imp::GridFunctionTraits<Range(Domain), typename Base::EntitySet, DefaultDerivativeTraits, 16>;
486 
487 private:
488 
489  template<class Node>
490  using LocalBasisRange = typename Node::FiniteElement::Traits::LocalBasisType::Traits::JacobianType;
491  template<class Node>
492  using NodeData = typename std::vector< LocalBasisRange<Node> >;
493  using PerNodeEvaluationBuffer = typename TypeTree::TreeContainer<NodeData, typename Base::Tree>;
494 
495 public:
496 
505  : public Base::LocalFunctionBase
506  {
507  using LocalBase = typename Base::LocalFunctionBase;
508  using size_type = typename Base::Tree::size_type;
509  using LocalBase::nodeToRangeEntry;
510 
511  public:
513  using Domain = typename LocalBase::Domain;
514  using Range = GlobalFunction::Range;
515  using Element = typename LocalBase::Element;
516 
518  LocalFunction(const GlobalFunction& globalFunction)
519  : LocalBase(globalFunction.data_)
520  , evaluationBuffer_(this->localView_.tree())
521  {
522  /* Nothing. */
523  }
524 
531  void bind(const Element& element)
532  {
533  LocalBase::bind(element);
534  geometry_.emplace(element.geometry());
535  }
536 
538  void unbind()
539  {
540  geometry_.reset();
541  LocalBase::unbind();
542  }
543 
557  Range operator()(const Domain& x) const
558  {
559  Range y;
560  istlVectorBackend(y) = 0;
561 
562  const auto& jacobianInverse = geometry_->jacobianInverse(x);
563 
564  TypeTree::forEachLeafNode(this->localView_.tree(), [&](auto&& node, auto&& treePath) {
565  const auto& fe = node.finiteElement();
566  const auto& localBasis = fe.localBasis();
567  auto& shapeFunctionJacobians = evaluationBuffer_[treePath];
568 
569  localBasis.evaluateJacobian(x, shapeFunctionJacobians);
570 
571  // Compute linear combinations of basis function jacobian.
572  // Non-scalar coefficients of dimension coeffDim are handled by
573  // processing the coeffDim linear combinations independently
574  // and storing them as entries of an array.
575  using RefJacobian = LocalBasisRange< std::decay_t<decltype(node)> >;
576  static constexpr auto coeffDim = decltype(flatVectorView(this->localDoFs_[node.localIndex(0)]).size())::value;
577  auto refJacobians = std::array<RefJacobian, coeffDim>{};
578  istlVectorBackend(refJacobians) = 0;
579  for (size_type i = 0; i < localBasis.size(); ++i)
580  {
581  auto c = flatVectorView(this->localDoFs_[node.localIndex(i)]);
582  for (std::size_t j = 0; j < coeffDim; ++j)
583  refJacobians[j].axpy(c[j], shapeFunctionJacobians[i]);
584  }
585 
586  // Transform Jacobians form local to global coordinates.
587  using Jacobian = decltype(refJacobians[0] * jacobianInverse);
588  auto jacobians = std::array<Jacobian, coeffDim>{};
590  refJacobians.begin(), refJacobians.end(), jacobians.begin(),
591  [&](const auto& refJacobian) { return refJacobian * jacobianInverse; });
592 
593  // Assign computed Jacobians to node entry of range.
594  // Types are matched using the lexicographic ordering provided by flatVectorView.
595  LocalBase::assignWith(nodeToRangeEntry(node, treePath, y), jacobians);
596  });
597 
598  return y;
599  }
600 
602  friend typename Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction&)
603  {
604  DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
605  }
606 
607  private:
608  mutable PerNodeEvaluationBuffer evaluationBuffer_;
609  std::optional<typename Element::Geometry> geometry_;
610  };
611 
618  DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr<const Data>& data)
619  : Base(data)
620  {
621  /* Nothing. */
622  }
623 
625  Range operator()(const Domain& x) const
626  {
627  // TODO: Implement this using hierarchic search
628  DUNE_THROW(NotImplemented,"not implemented");
629  }
630 
631  friend typename Traits::DerivativeInterface derivative(const DiscreteGlobalBasisFunctionDerivative& f)
632  {
633  DUNE_THROW(NotImplemented, "derivative of derivative is not implemented");
634  }
635 
638  {
639  return LocalFunction(f);
640  }
641 };
642 
643 
644 } // namespace Functions
645 } // namespace Dune
646 
647 #endif // DUNE_FUNCTIONS_GRIDFUNCTIONS_DISCRETEGLOBALBASISFUNCTIONS_HH
local function evaluating the derivative in reference coordinates
Definition: discreteglobalbasisfunction.hh:506
Range operator()(const Domain &x) const
Evaluate this local-function in coordinates x in the bound element.
Definition: discreteglobalbasisfunction.hh:557
friend Traits::LocalFunctionTraits::DerivativeInterface derivative(const LocalFunction &)
Not implemented.
Definition: discreteglobalbasisfunction.hh:602
void unbind()
Unbind the local-function.
Definition: discreteglobalbasisfunction.hh:538
LocalFunction(const GlobalFunction &globalFunction)
Create a local function from the associated grid function.
Definition: discreteglobalbasisfunction.hh:518
void bind(const Element &element)
Bind LocalFunction to grid element.
Definition: discreteglobalbasisfunction.hh:531
Derivative of a DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:472
Range operator()(const Domain &x) const
Not implemented.
Definition: discreteglobalbasisfunction.hh:625
friend LocalFunction localFunction(const DiscreteGlobalBasisFunctionDerivative &f)
Construct local function from a DiscreteGlobalBasisFunctionDerivative
Definition: discreteglobalbasisfunction.hh:637
DiscreteGlobalBasisFunctionDerivative(const std::shared_ptr< const Data > &data)
create object from DiscreateGlobalBasisFunction data
Definition: discreteglobalbasisfunction.hh:618
A grid function induced by a global basis and a coefficient vector.
Definition: discreteglobalbasisfunction.hh:269
friend DiscreteGlobalBasisFunctionDerivative< DiscreteGlobalBasisFunction > derivative(const DiscreteGlobalBasisFunction &f)
Derivative of the DiscreteGlobalBasisFunction
Definition: discreteglobalbasisfunction.hh:389
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:372
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:377
friend LocalFunction localFunction(const DiscreteGlobalBasisFunction &t)
Construct local function from a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:402
Range operator()(const Domain &x) const
Not implemented.
Definition: discreteglobalbasisfunction.hh:382
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
Default exception for dummy implementations.
Definition: exceptions.hh:263
Traits for type conversions and type information.
constexpr auto models()
Check if concept is modeled by given types.
Definition: concept.hh:184
typename AutonomousValueType< T >::type AutonomousValue
Type free of internal references that T can be converted to.
Definition: typetraits.hh:558
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:39
auto makeDiscreteGlobalBasisFunction(B &&basis, V &&vector)
Generate a DiscreteGlobalBasisFunction.
Definition: discreteglobalbasisfunction.hh:432
auto istlVectorBackend(Vector &v)
Return a vector backend wrapping non-const ISTL like containers.
Definition: istlvectorbackend.hh:346
constexpr HybridTreePath< T... > treePath(const T &... t)
Constructs a new HybridTreePath from the given indices.
Definition: treepath.hh:191
void forEachLeafNode(Tree &&tree, LeafFunc &&leafFunc)
Traverse tree and visit each leaf node.
Definition: traversal.hh:304
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:304
Dune namespace.
Definition: alignedallocator.hh:13
std::vector< decltype(std::declval< Op >)(std::declval< T >))) > transform(const std::vector< T > &in, Op op)
copy a vector, performing an operation on each element
Definition: misc.hh:24
auto wrap_or_move(T &&t)
Capture R-value reference to shared_ptr.
Definition: shared_ptr.hh:96
Default implementation for derivative traits.
Definition: defaultderivativetraits.hh:37
A simple node to range map using the nested tree indices.
Definition: hierarchicnodetorangemap.hh:30
Helper class to deduce the signature of a callable.
Definition: signature.hh:56
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 25, 22:37, 2024)