AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
InterpolationDataTransfer.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#include <dune/common/typetree/treecontainer.hh>
18
19#include <dune/grid/common/geometry.hh>
20#include <dune/grid/common/mcmgmapper.hh>
21#include <dune/grid/common/rangegenerators.hh>
22
23#include <amdis/Output.hpp>
24#include <amdis/common/ConcurrentCache.hpp>
25#include <amdis/functions/EntitySet.hpp>
26#include <amdis/operations/Assigner.hpp>
27#include <amdis/typetree/Traversal.hpp>
28
29namespace AMDiS {
30namespace Impl {
31
32// Hash function for cache container
33struct CoordHasher
34{
35 template <class LocalCoord>
36 std::size_t operator()(LocalCoord const& coord) const
37 {
38 std::size_t seed = 0;
39 for (std::size_t i = 0; i < coord.size(); ++i)
40 Dune::hash_combine(seed, coord[i]);
41 return seed;
42 }
43};
44
45} // end namespace Impl
46
47
48template <class B, class C>
50preAdapt(B const& basis, C const& coeff, bool mightCoarsen)
51{
52 GridView gv = basis.gridView();
53 LocalView lv = basis.localView();
54 auto const& grid = gv.grid();
55 auto const& idSet = grid.localIdSet();
56
57 nodeDataTransfer_.resize(lv.tree());
58 Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
59 nodeDataTransfer_[tp].preAdaptInit(lv, coeff, node);
60 });
61
62 // Make persistent DoF container
63 persistentContainer_.clear(); // Redundant if postAdapt was correctly called last cycle
64 for (const auto& e : elements(gv))
65 {
66 auto it = persistentContainer_.emplace(idSet.id(e),
67 Dune::TypeTree::makeTreeContainer(lv.tree(), [](auto&& node) { return NodeElementData<TYPEOF(node)>{}; }));
68
69 lv.bind(e);
70 auto& treeContainer = it.first->second;
71 Traversal::forEachLeafNode(lv.tree(), [&](auto const& /*node*/, auto const& tp) {
72 nodeDataTransfer_[tp].cacheLocal(treeContainer[tp]);
73 });
74 }
75
76 if (!mightCoarsen)
77 return;
78
79 // Interpolate from possibly vanishing elements
80 auto maxLevel = grid.maxLevel();
81 using std::sqrt;
82 typename Grid::ctype const checkInsideTolerance = sqrt(std::numeric_limits<typename Grid::ctype>::epsilon());
83 using Seed = typename Grid::template Codim<0>::EntitySeed;
84 std::vector<Seed> seeds(maxLevel+1);
85
86 for (auto const& e : entitySet(basis))
87 {
88 auto father = e;
89 while (father.mightVanish() && father.hasFather())
90 {
91 father = father.father();
92 auto it = persistentContainer_.emplace(idSet.id(father),
93 Dune::TypeTree::makeTreeContainer(lv.tree(), [](auto&& node) { return NodeElementData<TYPEOF(node)>{}; }));
94 if (!it.second)
95 continue;
96
97 auto& treeContainer = it.first->second;
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 {
103 seeds[hIt->level()] = hIt->seed(); // Save element in hierarchy to access geometryInFather at each step later
104
105 if (!hIt->isLeaf())
106 continue;
107
108 auto const& child = *hIt;
109 auto search = persistentContainer_.find(idSet.id(child));
110 assert(search != persistentContainer_.end());
111 auto const& childContainer = search->second;
112 lv.bind(child);
113
114 using BoolCoordPair = std::pair<bool, LocalCoordinate>;
115 using CacheImp = std::unordered_map<LocalCoordinate, BoolCoordPair, Impl::CoordHasher>;
116 using ChildCache = ConcurrentCache<LocalCoordinate, BoolCoordPair, ConsecutivePolicy, CacheImp>;
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 y = x;
122 for (int i = father.level()+1; i <= child.level(); ++i)
123 {
124 auto currentElement = grid.entity(seeds[i]);
125 auto const& geoInFather = currentElement.geometryInFather();
126 y = geoInFather.local(y);
127 // TODO(FM): Using an implementation detail as workaround for insufficient
128 // tolerance, see https://gitlab.dune-project.org/core/dune-grid/issues/84
129 bool isInside = Dune::Geo::Impl::checkInside(Dune::referenceElement(geoInFather).type().id(), Geometry::mydimension, y, checkInsideTolerance);
130 if (!isInside)
131 return BoolCoordPair(false, std::move(y));
132 }
133 return BoolCoordPair(true, std::move(y));
134 };
135 // Cache result of xInChild for subsequent calls from other basis nodes
136 // TODO(FM): Disable for single-node basis
137 ChildCache childCache;
138 auto xInChildCached = [&](LocalCoordinate const& x) -> BoolCoordPair {
139 return childCache.get(x, [&](LocalCoordinate const& x) { return xInChild(x); });
140 };
141
142 restrictLocalCompleted = true;
143 Traversal::forEachLeafNode(lv.tree(), [&](auto const& /*node*/, auto const& tp) {
144 restrictLocalCompleted &=
145 nodeDataTransfer_[tp].restrictLocal(father, treeContainer[tp], xInChildCached,
146 childContainer[tp], init);
147 });
148 init = false;
149 }
150 // test if restrictLocal was completed on all nodes
151 assert(restrictLocalCompleted);
152
153 } // end while (father.mightVanish)
154 } // end for (elements)
155}
156
157
158template <class B, class C>
159void InterpolationDataTransfer<B,C>::adapt(B const& basis, C& coeff)
160{
161 coeff.resize(basis);
162
163 // No data was saved before adapting the grid, make
164 // sure to call DataTransfer::preAdapt before calling adapt() on the grid
165 if (persistentContainer_.empty())
166 return;
167
168 GridView gv = basis.gridView();
169 LocalView lv = basis.localView();
170 auto const& idSet = gv.grid().localIdSet();
171 Traversal::forEachLeafNode(lv.tree(), [&](auto const& node, auto const& tp) {
172 nodeDataTransfer_[tp].adaptInit(lv, coeff, node);
173 });
174
175 using Mapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
176 Mapper mapper{gv, Dune::mcmgElementLayout()};
177
178 std::vector<bool> finished(mapper.size(), false);
179 for (const auto& e : entitySet(basis))
180 {
181 auto index = mapper.index(e);
182 if (finished[index])
183 continue;
184
185 auto it = persistentContainer_.find(idSet.id(e));
186
187 // Data already exists and no interpolation is required
188 if (it != persistentContainer_.end()) {
189 lv.bind(e);
190 auto const& treeContainer = it->second;
191 Traversal::forEachLeafNode(lv.tree(), [&](auto const& /*node*/, auto const& tp) {
192 nodeDataTransfer_[tp].copyLocal(treeContainer[tp]);
193 });
194 finished[index] = true;
195 continue;
196 }
197
198 // Data needs to be interpolated
199 auto father = e;
200 while (father.hasFather() && father.isNew())
201 father = father.father();
202
203 auto maxLevel = gv.grid().maxLevel();
204 bool init = true; // init flag for first call on new father element
205
206 auto father_it = persistentContainer_.find(idSet.id(father));
207 assert(father_it != persistentContainer_.end());
208 auto const& treeContainer = father_it->second;
209
210 auto hItEnd = father.hend(maxLevel);
211 for (auto hIt = father.hbegin(maxLevel); hIt != hItEnd; ++hIt) {
212 if (!hIt->isLeaf())
213 continue;
214
215 auto const& child = *hIt;
216 lv.bind(child);
217
218 // coordinate transform from child to father element
219 auto xInFather = [&](LocalCoordinate const& x) -> LocalCoordinate
220 {
221 auto y = x;
222 auto currentElement = child;
223 while (currentElement.level() != father.level())
224 {
225 y = currentElement.geometryInFather().global(y);
226 currentElement = currentElement.father();
227 }
228 return y;
229 };
230
231 Traversal::forEachLeafNode(lv.tree(), [&](auto const& /*node*/, auto const& tp) {
232 nodeDataTransfer_[tp].prolongLocal(father, treeContainer[tp], xInFather, init);
233 });
234
235 finished[mapper.index(child)] = true;
236 init = false;
237 }
238 } // end for (elements)
239
240 coeff.finish();
241}
242
243
244template <class B, class C>
246{
247 persistentContainer_.clear();
248}
249
250
254template <class Node, class Container, class Basis>
256{
257 using T = typename Container::value_type;
258 using LocalView = typename Basis::LocalView;
259 using Element = typename Node::Element;
260
261 using LocalBasis = typename Node::FiniteElement::Traits::LocalBasisType;
262 using LBRangeType = typename LocalBasis::Traits::RangeType;
263 using LocalInterpolation = typename Node::FiniteElement::Traits::LocalBasisType;
264 using LIDomainType = typename LocalInterpolation::Traits::DomainType;
265 using LIRangeType = typename LocalInterpolation::Traits::RangeType;
266
267public:
268 using NodeElementData = std::vector<T>;
269
270public:
271 NodeDataTransfer() = default;
272
274 void preAdaptInit(LocalView const& lv, Container const& coeff, Node const& node)
275 {
276 lv_ = &lv;
277 node_ = &node;
278 fatherNode_ = std::make_shared<Node>(node);
279 constCoeff_ = &coeff;
280 }
281
283
287 // [[expects: preAdaptInit to be called before]]
288 void cacheLocal(NodeElementData& dofs) const
289 {
290 constCoeff_->gather(*lv_, *node_, dofs);
291 }
292
307 // [[expects: preAdaptInit to be called before]]
308 template <class Trafo>
309 bool restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
310 NodeElementData const& childDOFs, bool init);
311
312
314 void adaptInit(LocalView const& lv, Container& coeff, Node const& node)
315 {
316 lv_ = &lv;
317 node_ = &node;
318 fatherNode_ = std::make_shared<Node>(node);
319 coeff_ = &coeff;
320 }
321
323 // [[expects: adaptInit to be called before]]
324 void copyLocal(NodeElementData const& dofs) const
325 {
326 coeff_->scatter(*lv_, *node_, dofs, Assigner::assign{});
327 }
328
340 // [[expects: adaptInit to be called before]]
341 template <class Trafo>
342 void prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
343 Trafo const& trafo, bool init);
344
345private:
346 LocalView const* lv_ = nullptr;
347 Node const* node_ = nullptr;
348 std::shared_ptr<Node> fatherNode_;
349 Container const* constCoeff_ = nullptr;
350 Container* coeff_ = nullptr;
351 std::vector<bool> finishedDOFs_;
352 NodeElementData fatherDOFsTemp_;
353};
354
355
356template <class N, class C, class B>
357 template <class Trafo>
359restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
360 NodeElementData const& childDOFs, bool init)
361{
362 auto& fatherNode = *fatherNode_;
363 std::size_t currentDOF = 0;
364 if (init)
365 {
366 // TODO(FM): This is UB, see https://gitlab.com/amdis/amdis/-/issues/16 (case 2)
367 bindTree(fatherNode, father);
368 }
369 auto const& childNode = *node_;
370 auto const& childFE = childNode.finiteElement();
371 auto const& fatherFE = fatherNode.finiteElement();
372
373 if (init) {
374 finishedDOFs_.assign(fatherFE.size(), false);
375 fatherDOFsTemp_.assign(fatherFE.size(), 0);
376 }
377
378 auto evalLeaf = [&](LIDomainType const& x) -> LIRangeType {
379 if (!finishedDOFs_[currentDOF])
380 {
381 auto const& insideLocal = trafo(x);
382 bool isInside = insideLocal.first;
383 if (isInside)
384 {
385 auto const& local = insideLocal.second;
386 thread_local std::vector<LBRangeType> shapeValues;
387 childFE.localBasis().evaluateFunction(local, shapeValues);
388
389 assert(childDOFs.size() == shapeValues.size());
390
391 LIRangeType y(0);
392 for (std::size_t i = 0; i < shapeValues.size(); ++i)
393 y += shapeValues[i] * childDOFs[i];
394
395 fatherDOFsTemp_[currentDOF] = T(y);
396 finishedDOFs_[currentDOF++] = true;
397 return y;
398 }
399 }
400 return fatherDOFsTemp_[currentDOF++];
401 };
402
403 fatherFE.localInterpolation().interpolate(evalLeaf, fatherDOFs);
404
405 // Return true if all father DOFs have been evaluated
406 return std::accumulate(finishedDOFs_.begin(), finishedDOFs_.end(), true,
407 std::logical_and<bool>());
408}
409
410
411template <class N, class C, class B>
412 template <class Trafo>
414prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
415 Trafo const& trafo, bool init)
416{
417 auto& fatherNode = *fatherNode_;
418 if (init)
419 {
420 // TODO(FM): This is UB, see https://gitlab.com/amdis/amdis/-/issues/16 (case 1)
421 bindTree(fatherNode, father);
422 }
423 auto const& childNode = *node_;
424
425 // evaluate father in child coordinate x
426 auto evalFather = [&](LIDomainType const& x) -> LIRangeType
427 {
428 thread_local std::vector<LBRangeType> shapeValues;
429 fatherNode.finiteElement().localBasis().evaluateFunction(trafo(x), shapeValues);
430 assert(shapeValues.size() == fatherDOFs.size());
431
432 LIRangeType y(0);
433 for (std::size_t i = 0; i < shapeValues.size(); ++i)
434 y += shapeValues[i] * fatherDOFs[i];
435
436 return y;
437 };
438
439 auto const& childFE = childNode.finiteElement();
440 thread_local std::vector<T> childDOFs;
441 childFE.localInterpolation().interpolate(evalFather, childDOFs);
442
443 coeff_->scatter(*lv_, childNode, childDOFs, Assigner::assign{});
444}
445
446} // end namespace AMDiS
Definition InterpolationDataTransfer.hpp:37
void preAdapt(Basis const &basis, Container const &coeff, bool mightCoarsen)
Definition InterpolationDataTransfer.inc.hpp:50
The restriction of a finite element basis to a single element.
Definition LocalView.hpp:16
Definition InterpolationDataTransfer.inc.hpp:256
void copyLocal(NodeElementData const &dofs) const
Copy already existing data to element bound to node_.
Definition InterpolationDataTransfer.inc.hpp:324
void cacheLocal(NodeElementData &dofs) const
Cache data on the element bound to node_.
Definition InterpolationDataTransfer.inc.hpp:288
void adaptInit(LocalView const &lv, Container &coeff, Node const &node)
To be called once before copyLocal/prolongLocal are called within the adapt step.
Definition InterpolationDataTransfer.inc.hpp:314
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 InterpolationDataTransfer.inc.hpp:274
Definition Assigner.hpp:8