DataTransfer.inc.hpp 16.1 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
21
22
23
24
#include <dune/grid/common/mcmgmapper.hh>
#include <dune/grid/common/rangegenerators.hh>

#include <dune/functions/common/functionfromcallable.hh>
#include <dune/functions/functionspacebases/lagrangebasis.hh>

#include <amdis/Output.hpp>
25
#include <amdis/common/ConcurrentCache.hpp>
26
#include <amdis/typetree/Traversal.hpp>
27
#include <amdis/typetree/TreeContainer.hpp>
28
29
30
31

namespace AMDiS
{
  template <class Node, class Container, class Basis>
32
  class NodeDataTransfer;
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85


  /** Data Transfer implementation for a single grid
   *  Handles computations related to the geometric information of the grid and passes that to the
   *    underlying NodeDataTransfer classes
   */
  template <class Container, class Basis>
  class DataTransfer
      : public DataTransferInterface<Container>
  {
    using LocalView = typename Basis::LocalView;
    using Tree = typename LocalView::Tree;
    using GridView = typename Basis::GridView;
    using Grid = typename GridView::Grid;
    using Mapper = Dune::LeafMultipleCodimMultipleGeomTypeMapper<Grid>;

    using IdSet = typename Grid::LocalIdSet;
    using IdType = typename IdSet::IdType;

    using Element = typename GridView::template Codim<0>::Entity;
    using Geometry = typename Element::Geometry;
    using LocalCoordinate = typename Geometry::LocalCoordinate;
    using BoolCoordPair = std::pair<bool, LocalCoordinate>;

    // Hash function for cache container
    struct Hasher
    {
      std::size_t operator()(LocalCoordinate const& coord) const
      {
        std::size_t seed = 0;
        for (std::size_t i = 0; i < coord.size(); ++i)
          Dune::hash_combine(seed, coord[i]);
        return seed;
      }
    };

    using CacheImp = std::unordered_map<LocalCoordinate, BoolCoordPair, Hasher>;
    using ChildCache = ConcurrentCache<LocalCoordinate, BoolCoordPair, ConsecutivePolicy, CacheImp>;

    constexpr static int d = Geometry::coorddimension;
    using ctype = typename Geometry::ctype;

    // Returns the Node's NodeDataTransfer
    struct NDT
    {
      template<class Node>
      auto operator()(const Node& node)
      {
        using NDT = NodeDataTransfer<Node, Container, Basis>;
        return NDT{};
      }
    };

86
    using NodeDataTransferContainer = TYPEOF(makeTreeContainer<Tree, NDT>(std::declval<const Tree&>(), NDT()));
87
88
89
90
91
92
93
94
95
96
97
98
99

    // Returns the Node's NodeElementData
    struct NodeElementData
    {
      template<class Node>
      auto operator()(const Node& node)
      {
        using NED = typename NodeDataTransfer<Node, Container, Basis>
          ::NodeElementData;
        return NED{};
      }
    };

Praetorius, Simon's avatar
Praetorius, Simon committed
100
    using ElementData = TYPEOF(makeTreeContainer<Tree, NodeElementData>(std::declval<const Tree&>(), NodeElementData()));
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

  public:
    // Container with data that persists during grid adaptation
    using PersistentContainer = std::map<IdType, ElementData>;

  public:
    DataTransfer(Basis const& basis)
      : basis_(&basis)
      , mapper_(basis.gridView().grid(), Dune::mcmgElementLayout())
      , finished_(mapper_.size(), false)
      , nodeDataTransfer_(makeTreeContainer<Tree, NDT>(basis_->localView().tree(), NDT()))
    {}

    /** Saves data contained in coeff in the PersistentContainer
     *  To be called after grid.preAdapt() and before grid.adapt()
     */
    void preAdapt(Container const& coeff, bool mightCoarsen) override;

    /** Unpacks data from the PersistentContainer
     *  To be called after grid.adapt() and before grid.postAdapt()
     */
    void postAdapt(Container& coeff, bool refined) override;

  private:
    Basis const* basis_;
    PersistentContainer persistentContainer_;
    Mapper mapper_;
    std::vector<bool> finished_;
    NodeDataTransferContainer nodeDataTransfer_;
  };

132
133
  template <class C, class B>
  void DataTransfer<C,B>::preAdapt(C const& coeff, bool mightCoarsen)
134
135
136
137
138
  {
    auto gv = basis_->gridView();
    auto lv = basis_->localView();
    auto const& idSet = gv.grid().localIdSet();

139
    for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
140
141
142
143
144
145
146
147
148
149
150
151
      nodeDataTransfer_[tp].preAdaptInit(lv, coeff, node);
    });

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

      lv.bind(e);
      auto& treeContainer = it.first->second;
152
      for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
153
154
155
156
157
158
159
        nodeDataTransfer_[tp].cacheLocal(treeContainer[tp]);
      });
    }

    // Interpolate from possibly vanishing elements
    if (mightCoarsen) {
      auto maxLevel = gv.grid().maxLevel();
160
161
      using std::sqrt;
      ctype const checkInsideTolerance = sqrt(std::numeric_limits<ctype>::epsilon());
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
      for (auto const& e : elements(gv))
      {
        auto father = e;
        while (father.mightVanish() && father.hasFather())
        {
          father = father.father();
          auto it = persistentContainer_.emplace(idSet.id(father),
            makeTreeContainer<Tree, NodeElementData>(lv.tree(), NodeElementData()));
          if (it.second) {
            auto& treeContainer = it.first->second;
            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()) {
                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 const& ref = Dune::referenceElement(childGeo);
                auto const& refTypeId = ref.type().id();

                // 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
194
                  bool isInside = Dune::Geo::Impl::checkInside(refTypeId, d, local, checkInsideTolerance);
195
196
197
198
199
200
201
202
203
204
                  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;
205
                for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
                  restrictLocalCompleted &=
                    nodeDataTransfer_[tp].restrictLocal(father, treeContainer[tp], xInChildCached,
                                                        childContainer[tp], init);
                });
                init = false;
              }
            }
            // test if restrictLocal was completed on all nodes
            assert(restrictLocalCompleted);
          }
        }
      }
    }
  }

221
222
  template <class C, class B>
  void DataTransfer<C,B>::postAdapt(C& coeff, bool refined)
223
224
225
226
227
  {
    coeff.resize(*basis_);
    auto gv = basis_->gridView();
    auto lv = basis_->localView();
    auto const& idSet = gv.grid().localIdSet();
228
    for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
      nodeDataTransfer_[tp].postAdaptInit(lv, coeff, node);
    });

    mapper_.update();
    finished_.clear();
    finished_.resize(mapper_.size(), false);
    for (const auto& e : elements(gv))
    {
      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;
247
        for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
          nodeDataTransfer_[tp].copyLocal(treeContainer[tp]);
        });
        finished_[index] = true;
      }
      // Data needs to be interpolated
      else {
        auto father = e;
        while (father.hasFather() && father.isNew())
          father = father.father();

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

        auto it = persistentContainer_.find(idSet.id(father));
        assert(it != persistentContainer_.end());
        auto const& treeContainer = it->second;

        auto hItEnd = father.hend(maxLevel);
        for (auto hIt = father.hbegin(maxLevel); hIt != hItEnd; ++hIt) {
          if (hIt->isLeaf()) {
            auto const& child = *hIt;
            lv.bind(child);

            // coordinate transform from child to father element
            auto xInFather = [&fatherGeo, childGeo = child.geometry()]
              (LocalCoordinate const& x) -> LocalCoordinate
            {
              return fatherGeo.local(childGeo.global(x));
            };

279
            for_each_leaf_node(lv.tree(), [&](auto const& node, auto const& tp) {
280
281
282
283
284
285
286
287
288
289
290
291
292
293
              nodeDataTransfer_[tp].prolongLocal(father, treeContainer[tp], xInFather, init);
            });

            finished_[mapper_.index(child)] = true;
            init = false;
          }
        }
      }
    }
  }

  /** Element-local data transfer on a single leaf node of the basis tree
   *  Handles computations related to the finite element basis node
   */
294
295
  template <class Node, class Container, class Basis>
  class NodeDataTransfer
296
297
298
  {
    using T = typename Container::value_type;
    using LocalView = typename Basis::LocalView;
299
    using Element = typename Node::Element;
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
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373

    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;
      auto nodeCopy = node;
      fatherNode_ = std::make_unique<Node>(std::move(nodeCopy));
      constCoeff_ = &coeff;
    }

    /// 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;
      auto nodeCopy = node;
      fatherNode_ = std::make_unique<Node>(std::move(nodeCopy));
      coeff_ = &coeff;
    }

    /// Cache data on the element bound to node_
    void cacheLocal(NodeElementData& dofs) const
    {
      auto const& fe = node_->finiteElement();
      auto feSize = fe.size();
      dofs.resize(feSize);
      for (std::size_t i = 0; i < feSize; ++i)
        dofs[i] = (*constCoeff_)[lv_->index(node_->localIndex(i))];
    }

    /** Evaluate data on the child element bound to node_ and interpolate onto father entity
     *    using the coordinate transformation trafo from father to child
     */
    template <class Trafo>
    bool restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
                       NodeElementData const& childDOFs, bool init);

    /// Copy already existing data to element bound to node_
    void copyLocal(NodeElementData const& dofs) const
    {
      for (std::size_t i = 0; i < dofs.size(); ++i)
        (*coeff_)[lv_->index(node_->localIndex(i))] = dofs[i];
    }

    /** Interpolate data from father onto the child element bound to node_ using the transformation
     *    trafo from child to father
     */
    template <class Trafo>
    void prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
                      Trafo const& trafo, bool init);

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

374
  template <class Node, class C, class B>
375
    template <class Trafo>
376
  bool NodeDataTransfer<Node,C,B>::
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
    restrictLocal(Element const& father, NodeElementData& fatherDOFs, Trafo const& trafo,
                  NodeElementData const& childDOFs, bool init)
  {
    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();
    if (init)
    {
      finishedDOFs_.assign(fatherFE.size(), false);
      fatherDOFsTemp_.assign(fatherFE.size(), 0);
    }

    auto evalLeaf = [&](LIDomainType const& x) -> LIRangeType {
      if (!finishedDOFs_[currentDOF])
      {
        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);

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

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

          fatherDOFsTemp_[currentDOF] = y;
          finishedDOFs_[currentDOF++] = true;
          return y;
        }
      }
      return fatherDOFsTemp_[currentDOF++];
    };

    using FFC
      = Dune::Functions::FunctionFromCallable<LIRangeType(LIDomainType), decltype(evalLeaf)>;
    fatherFE.localInterpolation().interpolate(FFC(evalLeaf), fatherDOFs);

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

430
  template <class Node, class C, class B>
431
    template <class Trafo>
432
  void NodeDataTransfer<Node,C,B>::
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
    prolongLocal(Element const& father, NodeElementData const& fatherDOFs,
                 Trafo const& trafo, bool init)
  {
    auto& fatherNode = *fatherNode_;
    if (init)
    {
      // TODO(FM): This is UB, replace with FE cache for father
      bindTree(fatherNode, father);
    }
    auto const& childNode = *node_;

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

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

      return y;
    };

    auto const& childFE = childNode.finiteElement();
    thread_local std::vector<T> childDOFs;
    using FFC = Dune::Functions
      ::FunctionFromCallable<LIRangeType(LIDomainType), decltype(evalFather)>;
    childFE.localInterpolation().interpolate(FFC(evalFather), childDOFs);

    for (std::size_t i = 0; i < childDOFs.size(); ++i)
      (*coeff_)[lv_->index(childNode.localIndex(i))] = childDOFs[i];
  }

} // end namespace AMDiS