continuousdatacollector.hh 4.24 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
#pragma once

3
#include <numeric>
4
#include "unstructureddatacollector.hh"
Praetorius, Simon's avatar
Praetorius, Simon committed
5

6
7
#include <dune/grid/utility/globalindexset.hh>

8
namespace Dune
Praetorius, Simon's avatar
Praetorius, Simon committed
9
10
11
{

/// Implementation of \ref DataCollector for linear cells, with continuous data.
12
template <class GridView, class Partition>
13
class ContinuousDataCollector
14
    : public UnstructuredDataCollectorInterface<GridView, ContinuousDataCollector<GridView,Partition>, Partition>
Praetorius, Simon's avatar
Praetorius, Simon committed
15
{
16
  using Self = ContinuousDataCollector;
17
18
19
20
21
  using Super = UnstructuredDataCollectorInterface<GridView, Self, Partition>;

public:
  using Super::dim;
  using Super::partition;
Praetorius, Simon's avatar
Praetorius, Simon committed
22
23

public:
24
  ContinuousDataCollector (GridView const& gridView)
Praetorius, Simon's avatar
Praetorius, Simon committed
25
26
27
    : Super(gridView)
  {}

28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
  /// Collect the vertex indices
  void updateImpl ()
  {
    numPoints_ = 0;
    indexMap_.resize(gridView_.size(dim));
    auto const& indexSet = gridView_.indexSet();
    for (auto const& vertex : vertices(gridView_, partition))
      indexMap_[indexSet.index(vertex)] = std::int64_t(numPoints_++);

    if (gridView_.comm().size() > 1) {
      auto&& e = elements(gridView_, partition);
      numCells_ = std::distance(std::begin(e), std::end(e));
    } else {
      numCells_ = gridView_.size(0);
    }
  }

Praetorius, Simon's avatar
Praetorius, Simon committed
45
46
47
  /// Return number of grid vertices
  std::uint64_t numPointsImpl () const
  {
48
    return numPoints_;
Praetorius, Simon's avatar
Praetorius, Simon committed
49
50
51
52
53
54
  }

  /// Return the coordinates of all grid vertices in the order given by the indexSet
  template <class T>
  std::vector<T> pointsImpl () const
  {
55
56
57
    std::vector<T> data;
    data.reserve(numPoints_ * 3);
    for (auto const& vertex : vertices(gridView_, partition)) {
Praetorius, Simon's avatar
Praetorius, Simon committed
58
59
      auto v = vertex.geometry().center();
      for (std::size_t j = 0; j < v.size(); ++j)
60
        data.emplace_back(v[j]);
Praetorius, Simon's avatar
Praetorius, Simon committed
61
      for (std::size_t j = v.size(); j < 3u; ++j)
62
        data.emplace_back(0);
Praetorius, Simon's avatar
Praetorius, Simon committed
63
64
65
66
    }
    return data;
  }

67
68
69
  /// Return a vector of global unique ids of the points
  std::vector<std::uint64_t> pointIdsImpl () const
  {
70
71
    std::vector<std::uint64_t> data;
    data.reserve(numPoints_);
72
    GlobalIndexSet<GridView> globalIndexSet(gridView_, dim);
73
74
    for (auto const& vertex : vertices(gridView_, partition)) {
      data.emplace_back(globalIndexSet.index(vertex));
75
76
77
78
    }
    return data;
  }

79
80
81
  /// Return number of grid cells
  std::uint64_t numCellsImpl () const
  {
82
    return numCells_;
83
84
  }

Praetorius, Simon's avatar
Praetorius, Simon committed
85
86
87
88
89
90
91
92
93
94
95
96
  /// Return the types, offsets and connectivity of the cells, using the same connectivity as
  /// given by the grid.
  Cells cellsImpl () const
  {
    auto const& indexSet = gridView_.indexSet();
    auto types = indexSet.types(0);
    int maxVertices = std::accumulate(types.begin(), types.end(), 1, [](int m, GeometryType t) {
      auto refElem = referenceElement<double,dim>(t);
      return std::max(m, refElem.size(dim));
    });

    Cells cells;
97
98
99
    cells.connectivity.reserve(numCells_ * maxVertices);
    cells.offsets.reserve(numCells_);
    cells.types.reserve(numCells_);
Praetorius, Simon's avatar
Praetorius, Simon committed
100
101

    std::int64_t old_o = 0;
102
    for (auto const& c : elements(gridView_, partition)) {
Praetorius, Simon's avatar
Praetorius, Simon committed
103
      Vtk::CellType cellType(c.type());
104
      for (unsigned int j = 0; j < c.subEntities(dim); ++j)
105
        cells.connectivity.emplace_back(indexMap_[indexSet.subIndex(c,cellType.permutation(j),dim)]);
Praetorius, Simon's avatar
Praetorius, Simon committed
106
107
108
109
110
111
112
113
114
115
      cells.offsets.push_back(old_o += c.subEntities(dim));
      cells.types.push_back(cellType.type());
    }
    return cells;
  }

  /// Evaluate the `fct` at the corners of the elements
  template <class T, class GlobalFunction>
  std::vector<T> pointDataImpl (GlobalFunction const& fct) const
  {
116
    std::vector<T> data(numPoints_ * fct.ncomps());
Praetorius, Simon's avatar
Praetorius, Simon committed
117
118
    auto const& indexSet = gridView_.indexSet();
    auto localFct = localFunction(fct);
119
    for (auto const& e : elements(gridView_, partition)) {
Praetorius, Simon's avatar
Praetorius, Simon committed
120
121
122
      localFct.bind(e);
      Vtk::CellType cellType{e.type()};
      auto refElem = referenceElement(e.geometry());
123
      for (unsigned int j = 0; j < e.subEntities(dim); ++j) {
124
        std::size_t idx = fct.ncomps() * indexMap_[indexSet.subIndex(e,cellType.permutation(j),dim)];
Praetorius, Simon's avatar
Praetorius, Simon committed
125
126
127
128
129
130
131
        for (int comp = 0; comp < fct.ncomps(); ++comp)
          data[idx + comp] = T(localFct.evaluate(comp, refElem.position(cellType.permutation(j),dim)));
      }
      localFct.unbind();
    }
    return data;
  }
132
133
134

protected:
  using Super::gridView_;
135
136
137
  std::uint64_t numPoints_ = 0;
  std::uint64_t numCells_ = 0;
  std::vector<std::int64_t> indexMap_;
Praetorius, Simon's avatar
Praetorius, Simon committed
138
139
};

140
} // end namespace Dune