AMDiS  0.3
The Adaptive Multi-Dimensional Simulation Toolbox
DataTransfer.inc.hpp
1 #pragma once
2 
3 #include <cmath>
4 #include <functional>
5 #include <limits>
6 #include <map>
7 #include <memory>
8 #include <numeric>
9 #include <type_traits>
10 #include <unordered_map>
11 #include <utility>
12 #include <vector>
13 
14 #include <dune/common/fvector.hh>
15 #include <dune/common/hash.hh>
16 #include <dune/common/version.hh>
17 
18 #include <dune/grid/common/geometry.hh>
19 #include <dune/grid/common/rangegenerators.hh>
20 
21 #include <amdis/Output.hpp>
22 #include <amdis/common/ConcurrentCache.hpp>
23 #include <amdis/functions/FunctionFromCallable.hpp>
24 #include <amdis/operations/Assigner.hpp>
25 #include <amdis/typetree/Traversal.hpp>
26 
27 namespace AMDiS {
28 namespace Impl {
29 
30 // Hash function for cache container
31 struct CoordHasher
32 {
33  template <class LocalCoord>
34  std::size_t operator()(LocalCoord const& coord) const
35  {
36  std::size_t seed = 0;
37  for (std::size_t i = 0; i < coord.size(); ++i)
38  Dune::hash_combine(seed, coord[i]);
39  return seed;
40  }
41 };
42 
43 } // end namespace Impl
44 
45 
46 template <class C, class B>
47 DataTransfer<C,B>::DataTransfer(std::shared_ptr<B const> basis)
48  : basis_(std::move(basis))
49  , mapper_(basis_->gridView(), Dune::mcmgElementLayout())
50  , nodeDataTransfer_()
51 {}
52 
53 
54 template <class C, class B>
56 preAdapt(C const& coeff, bool mightCoarsen)
57 {
58  GridView gv = basis_->gridView();
59  LocalView lv = basis_->localView();
60  auto const& idSet = gv.grid().localIdSet();
61 
62  Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
63  nodeDataTransfer_[tp].preAdaptInit(lv, coeff, node);
64  });
65 
66  // Make persistent DoF container
67  persistentContainer_.clear(); // Redundant if postAdapt was correctly called last cycle
68  for (const auto& e : elements(gv))
69  {
70  auto it = persistentContainer_.emplace(idSet.id(e), TypeTree::treeContainer<NodeElementData,true>(lv.tree()));
71 
72  lv.bind(e);
73  auto& treeContainer = it.first->second;
74  Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
75  nodeDataTransfer_[tp].cacheLocal(treeContainer[tp]);
76  });
77  }
78 
79  if (!mightCoarsen)
80  return;
81 
82  // Interpolate from possibly vanishing elements
83  auto maxLevel = gv.grid().maxLevel();
84  using std::sqrt;
85  typename Grid::ctype const checkInsideTolerance = sqrt(std::numeric_limits<typename Grid::ctype>::epsilon());
86  for (auto const& e : elements(gv, typename C::Traits::PartitionSet{}))
87  {
88  auto father = e;
89  while (father.mightVanish() && father.hasFather())
90  {
91  father = father.father();
92  auto it = persistentContainer_.emplace(idSet.id(father), TypeTree::treeContainer<NodeElementData,true>(lv.tree()));
93  if (!it.second)
94  continue;
95 
96  auto& treeContainer = it.first->second;
97  auto geo = father.geometry();
98  bool init = true; // init flag for first call on new father element
99  bool restrictLocalCompleted = false;
100  auto hItEnd = father.hend(maxLevel);
101  for (auto hIt = father.hbegin(maxLevel); hIt != hItEnd; ++hIt) {
102  if (!hIt->isLeaf())
103  continue;
104 
105  auto const& child = *hIt;
106  auto search = persistentContainer_.find(idSet.id(child));
107  assert(search != persistentContainer_.end());
108  auto const& childContainer = search->second;
109  lv.bind(child);
110 
111  auto const& childGeo = child.geometry();
112  auto refTypeId = Dune::referenceElement(childGeo).type().id();
113 
114  using BoolCoordPair = std::pair<bool, LocalCoordinate>;
115  using CacheImp = std::unordered_map<LocalCoordinate, BoolCoordPair, Impl::CoordHasher>;
117 
118  // Transfers input father-local point x into child-local point y
119  // Returns false if x is not inside the child
120  auto xInChild = [&](LocalCoordinate const& x) -> BoolCoordPair {
121  LocalCoordinate local = childGeo.local(geo.global(x));
122  // TODO(FM): Using an implementation detail as workaround for insufficient
123  // tolerance, see https://gitlab.dune-project.org/core/dune-grid/issues/84
124  bool isInside = Dune::Geo::Impl::checkInside(refTypeId, Geometry::mydimension, local, checkInsideTolerance);
125  return BoolCoordPair(isInside, std::move(local));
126  };
127  // TODO(FM): Disable for single-node basis
128  ChildCache childCache;
129  auto xInChildCached = [&](LocalCoordinate const& x) -> BoolCoordPair {
130  return childCache.get(x, [&](LocalCoordinate const& x) { return xInChild(x); });
131  };
132 
133  restrictLocalCompleted = true;
134  Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
135  restrictLocalCompleted &=
136  nodeDataTransfer_[tp].restrictLocal(father, treeContainer[tp], xInChildCached,
137  childContainer[tp], init);
138  });
139  init = false;
140  }
141  // test if restrictLocal was completed on all nodes
142  assert(restrictLocalCompleted);
143 
144  } // end while (father.mightVanish)
145  } // end for (elements)
146 }
147 
148 
149 template <class C, class B>
151 {
152  coeff.resize();
153 
154  // No data was saved before adapting the grid, make
155  // sure to call DataTransfer::preAdapt before calling adapt() on the grid
156  if (persistentContainer_.empty())
157  return;
158 
159  GridView gv = basis_->gridView();
160  LocalView lv = basis_->localView();
161  auto const& idSet = gv.grid().localIdSet();
162  Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
163  nodeDataTransfer_[tp].adaptInit(lv, coeff, node);
164  });
165 
166  mapper_.update(gv);
167 
168  std::vector<bool> finished(mapper_.size(), false);
169  for (const auto& e : elements(gv, typename C::Traits::PartitionSet{}))
170  {
171  auto index = mapper_.index(e);
172  if (finished[index])
173  continue;
174 
175  auto it = persistentContainer_.find(idSet.id(e));
176 
177  // Data already exists and no interpolation is required
178  if (it != persistentContainer_.end()) {
179  lv.bind(e);
180  auto const& treeContainer = it->second;
181  Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
182  nodeDataTransfer_[tp].copyLocal(treeContainer[tp]);
183  });
184  finished[index] = true;
185  continue;
186  }
187 
188  // Data needs to be interpolated
189  auto father = e;
190  while (father.hasFather() && father.isNew())
191  father = father.father();
192 
193  auto maxLevel = gv.grid().maxLevel();
194  auto fatherGeo = father.geometry();
195  bool init = true; // init flag for first call on new father element
196 
197  auto father_it = persistentContainer_.find(idSet.id(father));
198  assert(father_it != persistentContainer_.end());
199  auto const& treeContainer = father_it->second;
200 
201  auto hItEnd = father.hend(maxLevel);
202  for (auto hIt = father.hbegin(maxLevel); hIt != hItEnd; ++hIt) {
203  if (!hIt->isLeaf())
204  continue;
205 
206  auto const& child = *hIt;
207  lv.bind(child);
208 
209  // coordinate transform from child to father element
210  auto xInFather = [&fatherGeo, childGeo = child.geometry()]
211  (LocalCoordinate const& x) -> LocalCoordinate
212  {
213  return fatherGeo.local(childGeo.global(x));
214  };
215 
216  Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
217  nodeDataTransfer_[tp].prolongLocal(father, treeContainer[tp], xInFather, init);
218  });
219 
220  finished[mapper_.index(child)] = true;
221  init = false;
222  }
223  } // end for (elements)
224 
225  coeff.finish();
226 }
227 
228 
229 template <class C, class B>
231 {
232  persistentContainer_.clear();
233 }
234 
235 
239 template <class Node, class Container, class Basis>
240 class NodeDataTransfer
241 {
242  using T = typename Container::value_type;
243  using LocalView = typename Basis::LocalView;
244  using Element = typename Node::Element;
245 
246  using LocalBasis = typename Node::FiniteElement::Traits::LocalBasisType;
247  using LBRangeType = typename LocalBasis::Traits::RangeType;
248  using LocalInterpolation = typename Node::FiniteElement::Traits::LocalBasisType;
249  using LIDomainType = typename LocalInterpolation::Traits::DomainType;
250  using LIRangeType = typename LocalInterpolation::Traits::RangeType;
251 
252 public:
253  using NodeElementData = std::vector<T>;
254 
255 public:
256  NodeDataTransfer() = default;
257 
259  void preAdaptInit(LocalView const& lv, Container const& coeff, Node const& node)
260  {
261  lv_ = &lv;
262  node_ = &node;
263  fatherNode_ = std::make_unique<Node>(node);
264  constCoeff_ = &coeff;
265  }
266 
268 
272  // [[expects: preAdaptInit to be called before]]
273  void cacheLocal(NodeElementData& dofs) const
274  {
275  constCoeff_->gather(*lv_, *node_, dofs);
276  }
277 
292  // [[expects: preAdaptInit to be called before]]
293  template <class Trafo>
294  bool restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
295  NodeElementData const& childDOFs, bool init);
296 
297 
299  void adaptInit(LocalView const& lv, Container& coeff, Node const& node)
300  {
301  lv_ = &lv;
302  node_ = &node;
303  fatherNode_ = std::make_unique<Node>(node);
304  coeff_ = &coeff;
305  }
306 
308  // [[expects: adaptInit to be called before]]
309  void copyLocal(NodeElementData const& dofs) const
310  {
311  coeff_->scatter(*lv_, *node_, dofs, Assigner::assign{});
312  }
313 
325  // [[expects: adaptInit to be called before]]
326  template <class Trafo>
327  void prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
328  Trafo const& trafo, bool init);
329 
330 private:
331  LocalView const* lv_ = nullptr;
332  Node const* node_ = nullptr;
333  std::unique_ptr<Node> fatherNode_;
334  Container const* constCoeff_ = nullptr;
335  Container* coeff_ = nullptr;
336  std::vector<bool> finishedDOFs_;
337  NodeElementData fatherDOFsTemp_;
338 };
339 
340 
341 template <class N, class C, class B>
342  template <class Trafo>
344 restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
345  NodeElementData const& childDOFs, bool init)
346 {
347  auto& fatherNode = *fatherNode_;
348  std::size_t currentDOF = 0;
349  if (init)
350  {
351  // TODO(FM): This is UB, replace with FE cache for father
352  bindTree(fatherNode, father);
353  }
354  auto const& childNode = *node_;
355  auto const& childFE = childNode.finiteElement();
356  auto const& fatherFE = fatherNode.finiteElement();
357 
358  if (init) {
359  finishedDOFs_.assign(fatherFE.size(), false);
360  fatherDOFsTemp_.assign(fatherFE.size(), 0);
361  }
362 
363  auto evalLeaf = [&](LIDomainType const& x) -> LIRangeType {
364  if (!finishedDOFs_[currentDOF])
365  {
366  auto const& insideLocal = trafo(x);
367  bool isInside = insideLocal.first;
368  if (isInside)
369  {
370  auto const& local = insideLocal.second;
371  thread_local std::vector<LBRangeType> shapeValues;
372  childFE.localBasis().evaluateFunction(local, shapeValues);
373 
374  assert(childDOFs.size() == shapeValues.size());
375 
376  LIRangeType y(0);
377  for (std::size_t i = 0; i < shapeValues.size(); ++i)
378  y += shapeValues[i] * childDOFs[i];
379 
380  fatherDOFsTemp_[currentDOF] = T(y);
381  finishedDOFs_[currentDOF++] = true;
382  return y;
383  }
384  }
385  return fatherDOFsTemp_[currentDOF++];
386  };
387 
388  auto evalLeafFct = functionFromCallable<LIRangeType(LIDomainType)>(evalLeaf);
389  fatherFE.localInterpolation().interpolate(evalLeafFct, fatherDOFs);
390 
391  // Return true if all father DOFs have been evaluated
392  return std::accumulate(finishedDOFs_.begin(), finishedDOFs_.end(), true,
393  std::logical_and<bool>());
394 }
395 
396 
397 template <class N, class C, class B>
398  template <class Trafo>
400 prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
401  Trafo const& trafo, bool init)
402 {
403  auto& fatherNode = *fatherNode_;
404  if (init)
405  {
406  // TODO(FM): This is UB, replace with FE cache for father
407  bindTree(fatherNode, father);
408  }
409  auto const& childNode = *node_;
410 
411  // evaluate father in child coordinate x
412  auto evalFather = [&](LIDomainType const& x) -> LIRangeType
413  {
414  thread_local std::vector<LBRangeType> shapeValues;
415  fatherNode.finiteElement().localBasis().evaluateFunction(trafo(x), shapeValues);
416  assert(shapeValues.size() == fatherDOFs.size());
417 
418  LIRangeType y(0);
419  for (std::size_t i = 0; i < shapeValues.size(); ++i)
420  y += shapeValues[i] * fatherDOFs[i];
421 
422  return y;
423  };
424 
425  auto const& childFE = childNode.finiteElement();
426  thread_local std::vector<T> childDOFs;
427  auto evalFatherFct = functionFromCallable<LIRangeType(LIDomainType)>(evalFather);
428  childFE.localInterpolation().interpolate(evalFatherFct, childDOFs);
429 
430  coeff_->scatter(*lv_, childNode, childDOFs, Assigner::assign{});
431 }
432 
433 } // end namespace AMDiS
void copyLocal(NodeElementData const &dofs) const
Copy already existing data to element bound to node_.
Definition: DataTransfer.inc.hpp:309
void adapt(Container &coeff) override
Definition: DataTransfer.inc.hpp:150
The class template ConcurrentCache describes an associative static container that allows the concurre...
Definition: ConcurrentCache.hpp:53
void prolongLocal(Element const &father, NodeElementData const &fatherDOFs, Trafo const &trafo, bool init)
Interpolate data from father onto the child element bound to node_ using the transformation trafo fro...
Definition: DataTransfer.inc.hpp:400
Definition: AdaptiveGrid.hpp:373
Definition: FieldMatVec.hpp:12
void postAdapt(Container &coeff) override
Definition: DataTransfer.inc.hpp:230
Definition: DataTransfer.hpp:74
Definition: AdaptBase.hpp:6
constexpr bool LocalView
A Dune::Functions::LocalView type.
Definition: Concepts.hpp:182
void adaptInit(LocalView const &lv, Container &coeff, Node const &node)
To be called once before copyLocal/prolongLocal are called within the adapt step. ...
Definition: DataTransfer.inc.hpp:299
void preAdapt(Container const &coeff, bool mightCoarsen) override
Definition: DataTransfer.inc.hpp:56
bool restrictLocal(Element const &father, NodeElementData &fatherDOFs, Trafo const &trafo, NodeElementData const &childDOFs, bool init)
Evaluate data on the child element bound to node_ and interpolate onto father entity using the coordi...
Definition: DataTransfer.inc.hpp:344
void cacheLocal(NodeElementData &dofs) const
Cache data on the element bound to node_.
Definition: DataTransfer.inc.hpp:273
void preAdaptInit(LocalView const &lv, Container const &coeff, Node const &node)
To be called once before cacheLocal/restrictLocal are called within the preAdapt step.
Definition: DataTransfer.inc.hpp:259
Definition: Assigner.hpp:7