AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
DiscreteLocalFunction.inc.hpp
1 #pragma once
2 
3 #include <amdis/Output.hpp>
4 #include <amdis/common/DerivativeTraits.hpp>
5 #include <amdis/common/FieldMatVec.hpp>
6 #include <amdis/common/Math.hpp>
7 #include <amdis/common/RecursiveForEach.hpp>
8 #include <amdis/functions/HierarchicNodeToRangeMap.hpp>
9 #include <amdis/functions/NodeIndices.hpp>
10 #include <amdis/functions/Order.hpp>
11 #include <amdis/typetree/FiniteElementType.hpp>
12 #include <amdis/utility/LocalToGlobalAdapter.hpp>
13 
14 #include <dune/common/ftraits.hh>
15 
16 namespace AMDiS {
17 namespace Impl {
18 
19 // specialization for containers with gather method
20 template <class Coeff, class LocalView, class LocalCoeff,
21  class = decltype(std::declval<Coeff>().gather(std::declval<LocalView>(), std::declval<LocalCoeff&>()))>
22 void gather(Coeff const& coeff, LocalView const& localView, LocalCoeff& localCoeff, Dune::PriorityTag<2>)
23 {
24  coeff.gather(localView, localCoeff);
25 }
26 
27 // implementation for containers with multi-index access
28 template <class Coeff, class LocalView, class LocalCoeff,
29  class = decltype(std::declval<Coeff>()[std::declval<typename LocalView::MultiIndex>()])>
30 void gather(Coeff const& coeff, LocalView const& localView, LocalCoeff& localCoeff, Dune::PriorityTag<1>)
31 {
32  localCoeff.resize(localView.size());
33  auto it = localCoeff.begin();
34  for (auto const& idx : nodeIndices(localView))
35  *it++ = coeff[idx];
36 }
37 
38 // fallback implementation
39 template <class Coeff, class LocalView, class LocalCoeff,
40  class = decltype(std::declval<Coeff>()[std::integral_constant<std::size_t,0>{}])>
41 void gather(Coeff const& coeff, LocalView const& localView, LocalCoeff& localCoeff, Dune::PriorityTag<0>)
42 {
43  auto backend = Dune::Functions::istlVectorBackend(coeff);
44  localCoeff.resize(localView.size());
45  auto it = localCoeff.begin();
46  for (auto const& idx : nodeIndices(localView))
47  *it++ = backend[idx];
48 }
49 
50 } // end namespace Impl
51 
52 
53 template <class C, class GB, class TP, class R>
54  template <class Type>
56 {
57  using DomainType = typename DiscreteFunction::Domain;
58  using RangeType = typename DiscreteFunction::Range;
59 
60  template <class T>
61  using DerivativeRange = typename DerivativeTraits<RangeType(DomainType), T>::Range;
62 
63 public:
64  using Domain = typename EntitySet::LocalCoordinate;
65  using Range = DerivativeRange<Type>;
66 
67  enum { hasDerivative = std::is_same<Type, tag::value>::value };
68 
69 private:
70  using LocalView = typename GlobalBasis::LocalView;
71  using Element = typename EntitySet::Element;
72  using Geometry = typename Element::Geometry;
74 
75 public:
77  template <class LV>
78  LocalFunction(std::shared_ptr<LV> localView, TP const& treePath, C const& coefficients, Type type)
79  : localView_(std::move(localView))
80  , treePath_(treePath)
81  , coefficients_(coefficients)
82  , type_(type)
83  {}
84 
86  void bind(Element const& element)
87  {
88  localView_->bind(element);
89  geometry_.emplace(element.geometry());
90  Impl::gather(coefficients_, *localView_, localCoefficients_, Dune::PriorityTag<4>{});
91  bound_ = true;
92  }
93 
95  void unbind()
96  {
97  localView_->unbind();
98  geometry_.reset();
99  bound_ = false;
100  }
101 
103  Range operator() (const Domain& local) const
104  {
105  assert(bound_);
106  if constexpr (std::is_same_v<Type, tag::value>)
107  return evaluateFunction(local);
108  else if constexpr (std::is_same_v<Type, tag::gradient>)
109  return evaluateJacobian(local);
110  else if constexpr (std::is_same_v<Type, tag::partial>)
111  return evaluatePartialDerivative(local);
112  else if constexpr (std::is_same_v<Type, tag::divergence>) {
113  static constexpr bool isVectorNode
114  = SubTree::isPower && SubTree::degree() == GridView::dimensionworld;
115  static_assert(isVectorNode);
116  if (isVectorNode)
117  return evaluateDivergence(local);
118  }
119 
120  return Range{};
121  }
122 
124  template <class DerType>
126  {
127  static_assert(std::is_same_v<Type, tag::value>);
128  return {localView_, treePath_, coefficients_, type};
129  }
130 
132  auto order() const
133  -> decltype(AMDiS::order(std::declval<SubTree>()))
134  {
135  assert(bound_);
136  return AMDiS::order(subTree());
137  }
138 
140  Element const& localContext() const
141  {
142  assert(bound_);
143  return localView_->element();
144  }
145 
146 private:
147 
148  template <class LocalCoordinate>
149  auto evaluateFunction (const LocalCoordinate& local) const;
150 
151  template <class LocalCoordinate>
152  auto evaluateJacobian (const LocalCoordinate& local) const;
153 
154  template <class LocalCoordinate>
155  auto evaluateDivergence (const LocalCoordinate& local) const;
156 
157  template <class LocalCoordinate>
158  auto evaluatePartialDerivative (const LocalCoordinate& local) const;
159 
160  // get the geometry of the bound element, that is constructed in the bind() method
161  Geometry const& geometry() const
162  {
163  assert(bound_);
164  return *geometry_;
165  }
166 
167  // return the sub-tree of the basis-tree associated to the tree-path
168  SubTree const& subTree() const
169  {
170  return Dune::TypeTree::child(localView_->tree(), treePath_);
171  }
172 
173  // return the sub-tree of the basis-tree-cache associated to the tree-path
174  // NOTE: this is only return in case the local-view does provide a treeCache() function
175  template <class LV,
176  class = decltype(std::declval<LV>().treeCache())>
177  auto const& subTreeCache(Dune::PriorityTag<1>) const
178  {
179  return Dune::TypeTree::child(localView_->treeCache(), treePath_);
180  }
181 
182  // construct a new tree-cache that stores a reference to the sub-tree
183  template <class>
184  auto subTreeCache(Dune::PriorityTag<0>) const
185  {
186  return makeNodeCache(subTree());
187  }
188 
189  decltype(auto) subTreeCache() const
190  {
191  return subTreeCache<LocalView>(Dune::PriorityTag<2>());
192  }
193 
194 private:
195  std::shared_ptr<LocalView> localView_;
196  TP treePath_;
197  Coefficients const& coefficients_;
198  Type type_;
199  NodeToRangeMap nodeToRangeMap_;
200 
201  std::optional<Geometry> geometry_;
202  std::vector<Coefficient> localCoefficients_;
203  bool bound_ = false;
204 };
205 
206 
207 template <class C, class GB, class TP, class R>
208  template <class Type>
209  template <class LocalCoordinate>
211 ::evaluateFunction(LocalCoordinate const& local) const
212 {
213  Range y;
214  Recursive::forEach(y, [](auto& v) { v = 0; });
215 
216  Traversal::forEachLeafNode(this->subTreeCache(), [&](auto const& node, auto const& tp)
217  {
218  auto const& shapeFunctionValues = node.localBasisValuesAt(local);
219  std::size_t size = node.finiteElement().size();
220 
221  // Get range entry associated to this node
222  auto re = Dune::Functions::flatVectorView(nodeToRangeMap_(node, tp, y));
223 
224  for (std::size_t i = 0; i < size; ++i) {
225  // Get coefficient associated to i-th shape function
226  auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);
227 
228  // Get value of i-th shape function
229  auto v = Dune::Functions::flatVectorView(shapeFunctionValues[i]);
230 
231  std::size_t dimC = c.size();
232  std::size_t dimV = v.size();
233  assert(dimC*dimV == std::size_t(re.size()));
234  for(std::size_t j = 0; j < dimC; ++j) {
235  auto&& c_j = c[j];
236  for(std::size_t k = 0; k < dimV; ++k)
237  re[j*dimV + k] += c_j*v[k];
238  }
239  }
240  });
241 
242  return y;
243 }
244 
245 
246 template <class C, class GB, class TP, class R>
247  template <class Type>
248  template <class LocalCoordinate>
250 ::evaluateJacobian(LocalCoordinate const& local) const
251 {
252  Range dy;
253  Recursive::forEach(dy, [](auto& v) { v = 0; });
254 
255  Traversal::forEachLeafNode(this->subTreeCache(), [&](auto const& node, auto const& tp)
256  {
257  LocalToGlobalBasisAdapter localBasis(node, this->geometry());
258  auto const& gradients = localBasis.gradientsAt(local);
259 
260  // Get range entry associated to this node
261  auto re = Dune::Functions::flatVectorView(nodeToRangeMap_(node, tp, dy));
262 
263  for (std::size_t i = 0; i < localBasis.size(); ++i) {
264  // Get coefficient associated to i-th shape function
265  auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);
266 
267  // Get value of i-th transformed reference gradient
268  auto grad = Dune::Functions::flatVectorView(gradients[i]);
269 
270  std::size_t dimC = c.size();
271  std::size_t dimV = grad.size();
272  assert(dimC*dimV == std::size_t(re.size()));
273  for(std::size_t j = 0; j < dimC; ++j) {
274  auto&& c_j = c[j];
275  for(std::size_t k = 0; k < dimV; ++k)
276  re[j*dimV + k] += c_j*grad[k];
277  }
278  }
279  });
280 
281  return dy;
282 }
283 
284 
285 template <class C, class GB, class TP, class R>
286  template <class Type>
287  template <class LocalCoordinate>
289 ::evaluateDivergence(LocalCoordinate const& local) const
290 {
291  static_assert(static_size_v<Range> == 1, "Divergence implemented for scalar coefficients only.");
292 
293  Range dy;
294  Recursive::forEach(dy, [](auto& v) { v = 0; });
295 
296  auto&& node = this->subTreeCache();
297 
298  LocalToGlobalBasisAdapter localBasis(node.child(0), this->geometry());
299  auto const& gradients = localBasis.gradientsAt(local);
300 
301  auto re = Dune::Functions::flatVectorView(dy);
302  assert(int(re.size()) == 1);
303  for (std::size_t i = 0; i < localBasis.size(); ++i) {
304  auto grad = Dune::Functions::flatVectorView(gradients[i]);
305 
306  assert(int(grad.size()) == GridView::dimensionworld);
307  for (std::size_t j = 0; j < GridView::dimensionworld; ++j)
308  re[0] += localCoefficients_[node.child(j).localIndex(i)] * grad[j];
309  }
310 
311  return dy;
312 }
313 
314 
315 template <class C, class GB, class TP, class R>
316  template <class Type>
317  template <class LocalCoordinate>
319 ::evaluatePartialDerivative(LocalCoordinate const& local) const
320 {
321  std::size_t comp = this->type_.comp;
322 
323  Range dy;
324  Recursive::forEach(dy, [](auto& v) { v = 0; });
325 
326  Traversal::forEachLeafNode(this->subTreeCache(), [&](auto const& node, auto const& tp)
327  {
328  LocalToGlobalBasisAdapter localBasis(node, this->geometry());
329  auto const& partial = localBasis.partialsAt(local, comp);
330 
331  // Get range entry associated to this node
332  auto re = Dune::Functions::flatVectorView(nodeToRangeMap_(node, tp, dy));
333 
334  for (std::size_t i = 0; i < localBasis.size(); ++i) {
335  // Get coefficient associated to i-th shape function
336  auto c = Dune::Functions::flatVectorView(localCoefficients_[node.localIndex(i)]);
337 
338  // Get value of i-th transformed reference partial_derivative
339  auto d_comp = Dune::Functions::flatVectorView(partial[i]);
340 
341  std::size_t dimC = c.size();
342  std::size_t dimV = d_comp.size();
343  assert(dimC*dimV == std::size_t(re.size()));
344  for(std::size_t j = 0; j < dimC; ++j) {
345  auto&& c_j = c[j];
346  for(std::size_t k = 0; k < dimV; ++k)
347  re[j*dimV + k] += c_j*d_comp[k];
348  }
349  }
350  });
351 
352  return dy;
353 }
354 
355 
356 } // end namespace AMDiS
LocalFunction(std::shared_ptr< LV > localView, TP const &treePath, C const &coefficients, Type type)
Constructor. Stores a copy of the DiscreteFunction.
Definition: DiscreteLocalFunction.inc.hpp:78
Convert a simple (scalar) local basis into a global basis.
Definition: LocalToGlobalAdapter.hpp:64
auto makeNodeCache(Node const &node)
Construct a new local-basis cache from a basis-node.
Definition: NodeCache.hpp:35
Definition: FieldMatVec.hpp:12
AMDiS::LocalView< Self > LocalView
Type of the local view on the restriction of the basis to a single element.
Definition: GlobalBasis.hpp:70
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition: Concepts.hpp:182
A mutable view on the subspace of a DOFVector,.
Definition: DiscreteFunction.hpp:36
std::size_t size() const
Return the number of local basis functions.
Definition: LocalToGlobalAdapter.hpp:96
A simple node to range map using the nested tree indices.
Definition: HierarchicNodeToRangeMap.hpp:23
auto order(F const &f) -> decltype(&F::operator(), f.order())
polynomial order of functions
Definition: Order.hpp:11
auto const & partialsAt(typename Traits::DomainLocal const &x, std::size_t comp) const
Definition: LocalToGlobalAdapter.hpp:176
LocalFunction< DerType > makeDerivative(DerType type) const
Create a LocalFunction representing the derivative.
Definition: DiscreteLocalFunction.inc.hpp:125
void bind(Element const &element)
Bind the LocalView to the element.
Definition: DiscreteLocalFunction.inc.hpp:86
Definition: DerivativeTraits.hpp:29
void unbind()
Unbind the LocalView from the element.
Definition: DiscreteLocalFunction.inc.hpp:95
auto order() const -> decltype(AMDiS::order(std::declval< SubTree >()))
The polynomial degree of the LocalFunctions basis-functions.
Definition: DiscreteLocalFunction.inc.hpp:132
Element const & localContext() const
Return the bound element.
Definition: DiscreteLocalFunction.inc.hpp:140
auto nodeIndices(LocalView const &localView, Node const &node)
Returns a range over (flat) DOF indices on a node, given by the localView.
Definition: NodeIndices.hpp:13
Definition: DiscreteLocalFunction.inc.hpp:55
auto const & gradientsAt(typename Traits::DomainLocal const &x) const
Definition: LocalToGlobalAdapter.hpp:147