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