DOFMapping.hpp 7.67 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
#pragma once

#include <array>
#include <algorithm>
#include <iterator>
#include <numeric>
#include <vector>

#if HAVE_MPI
#include <dune/common/parallel/remoteindices.hh>
11
#include <amdis/common/parallel/Communicator.hpp>
12
13
14
15
16
17
18
#endif

#include <amdis/Environment.hpp>
#include <amdis/linearalgebra/AttributeSet.hpp>

namespace AMDiS
{
Praetorius, Simon's avatar
Praetorius, Simon committed
19
  /// \brief Fallback for \ref ParallelDofMapping in case there is only one mpi core.
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
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
86
87
88
89
90
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
147
148
149
150
151
152
153
154
155
  template <class IS, class GI = std::size_t>
  class SequentialDofMapping
  {
    using IndexSet = IS;

  public:
    using size_type = std::size_t;
    using DofIndex = size_type;
    using LocalIndex = size_type;
    using GlobalIndex = GI;

  public:
    SequentialDofMapping() = default;

    template <class Communication>
    SequentialDofMapping(Communication& c)
    {
      update(c);
    }

    /// How many DOFs are owned by my processor?
    size_type localSize() const
    {
      return localSize_;
    }

    /// Return the sequence of number of local indices for all processors
    std::array<size_type,1> localSizes() const
    {
      return {localSize_};
    }

    /// The total number of global DOFs.
    size_type globalSize() const
    {
      return globalSize_;
    }

    /// Return the sequence of starting points of the global indices for all processors
    std::array<GlobalIndex,1> globalStarts() const
    {
      return {0u};
    }

    /// Return the vector of global indices
    std::vector<GlobalIndex> const& globalIndices() const
    {
      return indices_;
    }

    /// Return number of ghost indices
    GlobalIndex ghostSize() const
    {
      return 0;
    }

    /// Return the vector of ghost indices
    std::array<GlobalIndex,0> ghostIndices() const
    {
      return {};
    }

    /// Map global index to local ghost index.
    LocalIndex globalToGhost(GlobalIndex const& n) const
    {
      assert(false && "There are no ghost indices in sequential dofmappings");
      return 0;
    }

    /// Map DOF index to local ghost index
    LocalIndex dofToGhost(DofIndex const& n) const
    {
      assert(false && "There are no ghost indices in sequential dofmappings");
      return 0;
    }

    ///	Global index of local index n.
    GlobalIndex global(LocalIndex const& n) const
    {
      return n;
    }

    /// Map global index to consecutive local owner index
    LocalIndex globalToLocal(GlobalIndex const& n) const
    {
      return n;
    }

    /// Map DOF index to consecutive local owner index
    LocalIndex dofToLocal(DofIndex const& n) const
    {
      return n;
    }


    /// DOF index n is owned by this processor
    bool owner(DofIndex const& n) const
    {
      assert(n < localSize());
      return true;
    }

    /// Global index n is owned by this processor
    bool globalOwner(GlobalIndex const& n) const
    {
      return globalOwner(0, n);
    }

    /// Global index n is owned by processor p
    bool globalOwner(int p, GlobalIndex const& n) const
    {
      assert(p == 0);
      assert(n < globalSize());
      return true;
    }

    /// Update the local to global mapping. Must be called before mapping local to global
    template <class Communication>
    void update(Communication& c)
    {
      localSize_ = c.indexSet().size();
      globalSize_ = c.indexSet().size();
      indices_.resize(globalSize_);
      std::iota(indices_.begin(), indices_.end(), size_type(0));
    }

    void debug() const {}

  private:
    size_type localSize_ = 0;
    size_type globalSize_ = 0;
    std::vector<GlobalIndex> indices_;
  };


#if HAVE_MPI
Praetorius, Simon's avatar
Praetorius, Simon committed
156
  /// \brief Mapping of local to global indices
157
158
159
160
161
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
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
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
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
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
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
  template <class PIS, class GI = std::size_t>
  class ParallelDofMapping
  {
    using ParallelIndexSet = PIS;
    using Attribute = typename AttributeSet<PIS>::type;
    using RemoteIndices = Dune::RemoteIndices<ParallelIndexSet>;

  public:
    using size_type = std::size_t;
    using DofIndex = size_type;
    using LocalIndex = size_type;
    using GlobalIndex = GI;

  public:
    ParallelDofMapping() = default;

    template <class Communication>
    ParallelDofMapping(Communication& c)
    {
      update(c);
    }

    /// How many DOFs are owned by my processor?
    size_type localSize() const
    {
      return localSize_;
    }

    /// Return the sequence of number of local indices for all processors
    std::vector<size_type> const& localSizes() const
    {
      return sizes_;
    }

    /// The total number of global DOFs.
    size_type globalSize() const
    {
      return globalSize_;
    }

    /// Return the sequence of starting points of the global indices for all processors
    std::vector<GlobalIndex> const& globalStarts() const
    {
      return starts_;
    }

    /// Return vector of global indices
    std::vector<GlobalIndex> const& globalIndices() const
    {
      return globalIndices_;
    }

    /// Return number of ghost indices
    size_type ghostSize() const
    {
      return ghostSize_;
    }

    /// Return vector of global ghost indices
    std::vector<GlobalIndex> const& ghostIndices() const
    {
      return ghostIndices_;
    }

    /// Map global index to local ghost index. NOTE: expensive
    LocalIndex globalToGhost(GlobalIndex const& n) const
    {
      auto it = std::find(ghostIndices_.begin(), ghostIndices_.end(), n);
      assert(it != ghostIndices_.end());
      return std::distance(ghostIndices_.begin(), it);
    }

    /// Map DOF index to local ghost index
    LocalIndex dofToGhost(DofIndex const& n) const
    {
      assert(!owner(n));
      assert(n < ghostLocalIndices_.size());
      assert(ghostLocalIndices_[n] < ghostSize_);

      return ghostLocalIndices_[n];
    }

    /// Map DOF index to global index
    GlobalIndex global(DofIndex const& n) const
    {
      assert(n < globalIndices_.size());
      return globalIndices_[n];
    }


    /// Map global index to consecutive local owner index
    LocalIndex globalToLocal(GlobalIndex const& n) const
    {
      return n - starts_[Environment::mpiRank()];
    }

    /// Map DOF index to consecutive local owner index
    LocalIndex dofToLocal(DofIndex const& n) const
    {
      assert(n < globalIndices_.size());
      return globalToLocal(globalIndices_[n]);
    }

    /// DOF index n is owned by this processor
    bool owner(DofIndex const& n) const
    {
      assert(n < owner_.size());
      return owner_[n];
    }

    /// Global index n is owned by this processor
    bool globalOwner(GlobalIndex const& n) const
    {
      return globalOwner(Environment::mpiRank(), n);
    }

    /// Global index n is owned by processor p
    bool globalOwner(int p, GlobalIndex const& n) const
    {
      assert(p < Environment::mpiSize());
      return n >= starts_[p] && n < starts_[p+1];
    }

    /// Update the local to global mapping. Must be called before mapping local to global
    template <class Communication>
    void update(Communication& c);

    void debug() const;

  private:
    void reset()
    {
      sizes_.clear();
      starts_.clear();
      localSize_ = 0;
      globalSize_ = 0;
      ghostSize_ = 0;

      globalIndices_.clear();
      ghostIndices_.clear();
      ghostLocalIndices_.clear();
      owner_.clear();
    }

  private:
    std::vector<size_type> sizes_;
    std::vector<GlobalIndex> starts_;
    size_type localSize_;
    size_type globalSize_;
    size_type ghostSize_;

    std::vector<GlobalIndex> globalIndices_; // indexed by LocalIndex
    std::vector<GlobalIndex> ghostIndices_;
    std::vector<LocalIndex> ghostLocalIndices_;
    std::vector<bool> owner_;

    const Mpi::Tag tag_{7513};
  };

  template <class PIS,class GI>
  using DofMapping = ParallelDofMapping<PIS,GI>;
#else
  template <class IS, class GI>
  using DofMapping = SequentialDofMapping<IS,GI>;
#endif

} // end namespace AMDiS

#include <amdis/linearalgebra/DOFMapping.inc.hpp>