DataTransfer.inc.hpp 13.8 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
#pragma once

#include <cmath>
#include <functional>
#include <limits>
#include <map>
#include <memory>
#include <numeric>
#include <type_traits>
#include <unordered_map>
#include <utility>
#include <vector>

#include <dune/common/fvector.hh>
#include <dune/common/hash.hh>

17
#include <dune/grid/common/geometry.hh>
18
19
20
#include <dune/grid/common/rangegenerators.hh>

#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
30
// Hash function for cache container
struct CoordHasher
31
{
32
33
  template <class LocalCoord>
  std::size_t operator()(LocalCoord const& coord) const
34
  {
35
36
37
38
39
40
    std::size_t seed = 0;
    for (std::size_t i = 0; i < coord.size(); ++i)
      Dune::hash_combine(seed, coord[i]);
    return seed;
  }
};
41

42
} // end namespace Impl
43
44


45
template <class C, class B>
46
47
DataTransfer<C,B>::
DataTransfer(std::shared_ptr<B> basis)
48
49
50
51
  : basis_(std::move(basis))
  , mapper_(basis_->gridView().grid(), Dune::mcmgElementLayout())
  , nodeDataTransfer_(makeTreeContainer<NDT>(basis_->localView().tree()))
{}
52
53


54
template <class C, class B>
55
56
void DataTransfer<C,B>::
preAdapt(C const& coeff, bool mightCoarsen)
57
58
59
60
{
  GridView gv = basis_->gridView();
  LocalView lv = basis_->localView();
  auto const& idSet = gv.grid().localIdSet();
61

62
63
64
65
66
67
68
  for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
    nodeDataTransfer_[tp].preAdaptInit(lv, coeff, node);
  });

  // Make persistent DoF container
  persistentContainer_.clear();
  for (const auto& e : elements(gv))
69
  {
70
    auto it = persistentContainer_.emplace(idSet.id(e), makeTreeContainer<NodeElementData>(lv.tree()));
71

72
73
    lv.bind(e);
    auto& treeContainer = it.first->second;
74
    for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
75
      nodeDataTransfer_[tp].cacheLocal(treeContainer[tp]);
76
    });
77
78
79
80
  }

  if (!mightCoarsen)
    return;
81

82
83
84
85
  // Interpolate from possibly vanishing elements
  auto maxLevel = gv.grid().maxLevel();
  using std::sqrt;
  typename Grid::ctype const checkInsideTolerance = sqrt(std::numeric_limits<typename Grid::ctype>::epsilon());
86
  for (auto const& e : elements(gv, typename C::Backend::Traits::PartitionSet{}))
87
88
89
  {
    auto father = e;
    while (father.mightVanish() && father.hasFather())
90
    {
91
92
93
94
      father = father.father();
      auto it = persistentContainer_.emplace(idSet.id(father), makeTreeContainer<NodeElementData>(lv.tree()));
      if (!it.second)
        continue;
95
96

      auto& treeContainer = it.first->second;
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
      auto geo = father.geometry();
      bool init = true; // init flag for first call on new father element
      bool restrictLocalCompleted = false;
      auto hItEnd = father.hend(maxLevel);
      for (auto hIt = father.hbegin(maxLevel); hIt != hItEnd; ++hIt) {
        if (!hIt->isLeaf())
          continue;

        auto const& child = *hIt;
        auto search = persistentContainer_.find(idSet.id(child));
        assert(search != persistentContainer_.end());
        auto const& childContainer = search->second;
        lv.bind(child);

        auto const& childGeo = child.geometry();
        auto refTypeId = Dune::referenceElement(childGeo).type().id();

        using BoolCoordPair = std::pair<bool, LocalCoordinate>;
        using CacheImp = std::unordered_map<LocalCoordinate, BoolCoordPair, Impl::CoordHasher>;
        using ChildCache = ConcurrentCache<LocalCoordinate, BoolCoordPair, ConsecutivePolicy, CacheImp>;

        // Transfers input father-local point x into child-local point y
        // Returns false if x is not inside the child
        auto xInChild = [&](LocalCoordinate const& x) -> BoolCoordPair {
          LocalCoordinate local = childGeo.local(geo.global(x));
          // TODO(FM): Using an implementation detail as workaround for insufficient
          //   tolerance, see https://gitlab.dune-project.org/core/dune-grid/issues/84
          bool isInside = Dune::Geo::Impl::checkInside(refTypeId, Geometry::coorddimension, local, checkInsideTolerance);
          return BoolCoordPair(isInside, std::move(local));
        };
        // TODO(FM): Disable for single-node basis
        ChildCache childCache;
        auto xInChildCached = [&](LocalCoordinate const& x) -> BoolCoordPair {
          return childCache.get(x, [&](LocalCoordinate const& x) { return xInChild(x); });
        };

        restrictLocalCompleted = true;
        for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
          restrictLocalCompleted &=
            nodeDataTransfer_[tp].restrictLocal(father, treeContainer[tp], xInChildCached,
                                                childContainer[tp], init);
        });
        init = false;
      }
      // test if restrictLocal was completed on all nodes
      assert(restrictLocalCompleted);

    } // end while (father.mightVanish)
  } // end for (elements)
}


template <class C, class B>
150
151
void DataTransfer<C,B>::
postAdapt(C& coeff, bool refined)
152
{
153
154
  coeff.init(false);

155
156
157
158
159
160
161
162
163
  GridView gv = basis_->gridView();
  LocalView lv = basis_->localView();
  auto const& idSet = gv.grid().localIdSet();
  for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
    nodeDataTransfer_[tp].postAdaptInit(lv, coeff, node);
  });

  mapper_.update();
  std::vector<bool> finished(mapper_.size(), false);
164
  for (const auto& e : elements(gv, typename C::Backend::Traits::PartitionSet{}))
165
166
167
168
169
170
171
172
173
174
175
  {
    auto index = mapper_.index(e);
    if (finished[index])
      continue;

    auto it = persistentContainer_.find(idSet.id(e));

    // Data already exists and no interpolation is required
    if (it != persistentContainer_.end()) {
      lv.bind(e);
      auto const& treeContainer = it->second;
176
      for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
177
        nodeDataTransfer_[tp].copyLocal(treeContainer[tp]);
178
      });
179
180
      finished[index] = true;
      continue;
181
182
    }

183
184
185
186
    // Data needs to be interpolated
    auto father = e;
    while (father.hasFather() && father.isNew())
      father = father.father();
187

188
189
190
    auto maxLevel = gv.grid().maxLevel();
    auto fatherGeo = father.geometry();
    bool init = true; // init flag for first call on new father element
191

192
193
194
195
196
197
198
    auto father_it = persistentContainer_.find(idSet.id(father));
    assert(father_it != persistentContainer_.end());
    auto const& treeContainer = father_it->second;

    auto hItEnd = father.hend(maxLevel);
    for (auto hIt = father.hbegin(maxLevel); hIt != hItEnd; ++hIt) {
      if (!hIt->isLeaf())
199
200
        continue;

201
202
      auto const& child = *hIt;
      lv.bind(child);
203

204
205
206
207
208
209
210
211
212
213
214
215
216
      // coordinate transform from child to father element
      auto xInFather = [&fatherGeo, childGeo = child.geometry()]
        (LocalCoordinate const& x) -> LocalCoordinate
      {
        return fatherGeo.local(childGeo.global(x));
      };

      for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
        nodeDataTransfer_[tp].prolongLocal(father, treeContainer[tp], xInFather, init);
      });

      finished[mapper_.index(child)] = true;
      init = false;
217
    }
218
  } // end for (elements)
219
220

  coeff.finish();
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
}


/** Element-local data transfer on a single leaf node of the basis tree
 *  Handles computations related to the finite element basis node
 */
template <class Node, class Container, class Basis>
class NodeDataTransfer
{
  using T = typename Container::value_type;
  using LocalView = typename Basis::LocalView;
  using Element = typename Node::Element;

  using LocalBasis = typename Node::FiniteElement::Traits::LocalBasisType;
  using LBRangeType = typename LocalBasis::Traits::RangeType;
  using LocalInterpolation = typename Node::FiniteElement::Traits::LocalBasisType;
  using LIDomainType = typename LocalInterpolation::Traits::DomainType;
  using LIRangeType = typename LocalInterpolation::Traits::RangeType;

public:
  using NodeElementData = std::vector<T>;

public:
  NodeDataTransfer() = default;

  /// To be called once before cacheLocal/restrictLocal are called within the preAdapt step
  void preAdaptInit(LocalView const& lv, Container const& coeff, Node const& node)
  {
    lv_ = &lv;
    node_ = &node;
    fatherNode_ = std::make_unique<Node>(node);
    constCoeff_ = &coeff;
253
254
  }

255
256
257
258
259
260
261
  /// \brief Cache data on the element bound to node_
  /**
   * This functions is used whenever the element does not vanish and thus the
   * data can trivially be transferred to the new element
   **/
  // [[expects: preAdaptInit to be called before]]
  void cacheLocal(NodeElementData& dofs) const
262
  {
263
    constCoeff_->gather(*lv_, *node_, dofs);
264
  }
265

266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
  /** \brief Evaluate data on the child element bound to node_ and interpolate onto
   *  father entity using the coordinate transformation trafo from father to child.
   *
   * Stores cached data in the NodeElementData argument. After grid adaption the
   * data is copied by \ref copyLocal or \ref prolongLocal to the target element
   * in the new grid.
   *
   * \param father      The father element to interpolate to
   * \param fatherDOFs  Container to store the interpolated DOFs
   * \param trafo       Coordinate transform from local coordinates in father to local
   *                    coordinates in child element
   * \param childDOFs   DOF values from the child element
   * \param init        The father element is visited for the first time
   **/
  // [[expects: preAdaptInit to be called before]]
  template <class Trafo>
  bool restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
                      NodeElementData const& childDOFs, bool init);


  /// To be called once before copyLocal/prolongLocal are called within the postAdapt step
  void postAdaptInit(LocalView const& lv, Container& coeff, Node const& node)
  {
    lv_ = &lv;
    node_ = &node;
    fatherNode_ = std::make_unique<Node>(node);
    coeff_ = &coeff;
  }
294

295
296
297
298
  /// \brief Copy already existing data to element bound to node_
  // [[expects: postAdaptInit to be called before]]
  void copyLocal(NodeElementData const& dofs) const
  {
299
    coeff_->scatter(*lv_, *node_, dofs, Assigner::assign{});
300
  }
301

302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
  /** \brief Interpolate data from father onto the child element bound to node_ using
   *  the transformation trafo from child to father
   *
   *  Stores the interpolated data from father to child in the container \ref coeff_.
   *
   * \param father      The father element
   * \param fatherDOFs  DOF values cached on the father element before adapt
   * \param trafo       Coordinate transform from local coordinates in child to local
   *                    coordinates in father element
   * \param init        Father element is visited for the first time
   **/
  // [[expects: postAdaptInit to be called before]]
  template <class Trafo>
  void prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
                    Trafo const& trafo, bool init);

private:
  LocalView const* lv_ = nullptr;
  Node const* node_ = nullptr;
  std::unique_ptr<Node> fatherNode_;
  Container const* constCoeff_ = nullptr;
  Container* coeff_ = nullptr;
  std::vector<bool> finishedDOFs_;
  NodeElementData fatherDOFsTemp_;
};


template <class N, class C, class B>
  template <class Trafo>
bool NodeDataTransfer<N,C,B>::
332
333
restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
              NodeElementData const& childDOFs, bool init)
334
335
336
337
338
339
340
341
342
343
344
{
  auto& fatherNode = *fatherNode_;
  std::size_t currentDOF = 0;
  if (init)
  {
    // TODO(FM): This is UB, replace with FE cache for father
    bindTree(fatherNode, father);
  }
  auto const& childNode = *node_;
  auto const& childFE = childNode.finiteElement();
  auto const& fatherFE = fatherNode.finiteElement();
345

346
347
348
349
  if (init) {
    finishedDOFs_.assign(fatherFE.size(), false);
    fatherDOFsTemp_.assign(fatherFE.size(), 0);
  }
350

351
352
  auto evalLeaf = [&](LIDomainType const& x) -> LIRangeType {
    if (!finishedDOFs_[currentDOF])
353
    {
354
355
356
357
358
359
360
      auto const& insideLocal = trafo(x);
      bool isInside = insideLocal.first;
      if (isInside)
      {
        auto const& local = insideLocal.second;
        thread_local std::vector<LBRangeType> shapeValues;
        childFE.localBasis().evaluateFunction(local, shapeValues);
361

362
        assert(childDOFs.size() == shapeValues.size());
363

364
365
366
        LIRangeType y(0);
        for (std::size_t i = 0; i < shapeValues.size(); ++i)
          y += shapeValues[i] * childDOFs[i];
367

368
369
370
371
        fatherDOFsTemp_[currentDOF] = y;
        finishedDOFs_[currentDOF++] = true;
        return y;
      }
372
    }
373
    return fatherDOFsTemp_[currentDOF++];
374
375
  };

376
377
  auto evalLeafFct = functionFromCallable<LIRangeType(LIDomainType)>(evalLeaf);
  fatherFE.localInterpolation().interpolate(evalLeafFct, fatherDOFs);
378

379
380
381
382
  // Return true if all father DOFs have been evaluated
  return std::accumulate(finishedDOFs_.begin(), finishedDOFs_.end(), true,
                          std::logical_and<bool>());
}
383
384


385
386
387
template <class N, class C, class B>
  template <class Trafo>
void NodeDataTransfer<N,C,B>::
388
389
prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
             Trafo const& trafo, bool init)
390
391
392
393
394
395
{
  auto& fatherNode = *fatherNode_;
  if (init)
  {
    // TODO(FM): This is UB, replace with FE cache for father
    bindTree(fatherNode, father);
396
  }
397
  auto const& childNode = *node_;
398

399
400
  // evaluate father in child coordinate x
  auto evalFather = [&](LIDomainType const& x) -> LIRangeType
401
  {
402
403
404
    thread_local std::vector<LBRangeType> shapeValues;
    fatherNode.finiteElement().localBasis().evaluateFunction(trafo(x), shapeValues);
    assert(shapeValues.size() == fatherDOFs.size());
405

406
407
408
    LIRangeType y(0);
    for (std::size_t i = 0; i < shapeValues.size(); ++i)
      y += shapeValues[i] * fatherDOFs[i];
409

410
411
    return y;
  };
412

413
414
  auto const& childFE = childNode.finiteElement();
  thread_local std::vector<T> childDOFs;
415
416
  auto evalFatherFct = functionFromCallable<LIRangeType(LIDomainType)>(evalFather);
  childFE.localInterpolation().interpolate(evalFatherFct, childDOFs);
417

418
  coeff_->scatter(*lv_, childNode, childDOFs, Assigner::assign{});
419
}
420
421

} // end namespace AMDiS