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

3
#include <numeric>
4
#include <vector>
Praetorius, Simon's avatar
Praetorius, Simon committed
5

6
#include <dune/geometry/referenceelements.hh>
7
#include <dune/grid/common/partitionset.hh>
8
#include <dune/grid/utility/globalindexset.hh>
Stenger, Florian's avatar
Stenger, Florian committed
9
#include <dune/vtk/types.hh>
10
11
12

#include "unstructureddatacollector.hh"

13
namespace Dune
Praetorius, Simon's avatar
Praetorius, Simon committed
14
{
Stenger, Florian's avatar
Stenger, Florian committed
15
16
17
  namespace Vtk
  {
    /// Implementation of \ref DataCollector for linear cells, with continuous data.
18
    template <class GridView, class Partition = Partitions::InteriorBorder>
Stenger, Florian's avatar
Stenger, Florian committed
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
    class ContinuousDataCollector
        : public UnstructuredDataCollectorInterface<GridView, ContinuousDataCollector<GridView,Partition>, Partition>
    {
      using Self = ContinuousDataCollector;
      using Super = UnstructuredDataCollectorInterface<GridView, Self, Partition>;

    public:
      using Super::dim;
      using Super::partition;

    public:
      ContinuousDataCollector (GridView const& gridView)
        : Super(gridView)
      {}

      /// 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
50

Stenger, Florian's avatar
Stenger, Florian committed
51
52
53
54
55
      /// Return number of grid vertices
      std::uint64_t numPointsImpl () const
      {
        return numPoints_;
      }
Praetorius, Simon's avatar
Praetorius, Simon committed
56

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

Stenger, Florian's avatar
Stenger, Florian committed
73
74
75
76
77
78
79
80
81
82
83
      /// Return a vector of global unique ids of the points
      std::vector<std::uint64_t> pointIdsImpl () const
      {
        std::vector<std::uint64_t> data;
        data.reserve(numPoints_);
        GlobalIndexSet<GridView> globalIndexSet(gridView_, dim);
        for (auto const& vertex : vertices(gridView_, partition)) {
          data.emplace_back(globalIndexSet.index(vertex));
        }
        return data;
      }
84

Stenger, Florian's avatar
Stenger, Florian committed
85
86
87
88
      /// Return number of grid cells
      std::uint64_t numCellsImpl () const
      {
        return numCells_;
Praetorius, Simon's avatar
Praetorius, Simon committed
89
90
      }

Stenger, Florian's avatar
Stenger, Florian committed
91
92
93
94
95
96
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
      /// 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;
        cells.connectivity.reserve(numCells_ * maxVertices);
        cells.offsets.reserve(numCells_);
        cells.types.reserve(numCells_);

        std::int64_t old_o = 0;
        for (auto const& c : elements(gridView_, partition)) {
          Vtk::CellType cellType(c.type());
          for (unsigned int j = 0; j < c.subEntities(dim); ++j)
            cells.connectivity.emplace_back(indexMap_[indexSet.subIndex(c,cellType.permutation(j),dim)]);
          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
      {
        std::vector<T> data(numPoints_ * fct.ncomps());
        auto const& indexSet = gridView_.indexSet();
        auto localFct = localFunction(fct);
        for (auto const& e : elements(gridView_, partition)) {
          localFct.bind(e);
          Vtk::CellType cellType{e.type()};
          auto refElem = referenceElement(e.geometry());
          for (unsigned int j = 0; j < e.subEntities(dim); ++j) {
            std::size_t idx = fct.ncomps() * indexMap_[indexSet.subIndex(e,cellType.permutation(j),dim)];
            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;
      }

    protected:
      using Super::gridView_;
      std::uint64_t numPoints_ = 0;
      std::uint64_t numCells_ = 0;
      std::vector<std::int64_t> indexMap_;
    };

  } // end namespace Vtk
147
} // end namespace Dune