globalindices.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_GLOBAL_INDICES_HH
28#define EWOMS_GLOBAL_INDICES_HH
29
30#include <dune/grid/common/datahandleif.hh>
31#include <dune/istl/bcrsmatrix.hh>
32#include <dune/istl/scalarproducts.hh>
33#include <dune/istl/operators.hh>
34
35#include <algorithm>
36#include <set>
37#include <map>
38#include <iostream>
39#include <tuple>
40
41#if HAVE_MPI
42#include <mpi.h>
43#endif
44
45#include "overlaptypes.hh"
46
47namespace Opm {
48namespace Linear {
54template <class ForeignOverlap>
56{
57 GlobalIndices(const GlobalIndices& ) = delete;
58
59 using GlobalToDomesticMap = std::map<Index, Index>;
60 using DomesticToGlobalMap = std::map<Index, Index>;
61
62public:
63 GlobalIndices(const ForeignOverlap& foreignOverlap)
64 : foreignOverlap_(foreignOverlap)
65 {
66 myRank_ = 0;
67 mpiSize_ = 1;
68
69#if HAVE_MPI
70 {
71 int tmp;
72 MPI_Comm_rank(MPI_COMM_WORLD, &tmp);
73 myRank_ = static_cast<ProcessRank>(tmp);
74 MPI_Comm_size(MPI_COMM_WORLD, &tmp);
75 mpiSize_ = static_cast<size_t>(tmp);
76 }
77#endif
78
79 // calculate the domestic overlap (i.e. all overlap indices in
80 // foreign processes which the current process overlaps.)
81 // This requires communication via MPI.
83 }
84
88 Index domesticToGlobal(Index domesticIdx) const
89 {
90 assert(domesticToGlobal_.find(domesticIdx) != domesticToGlobal_.end());
91
92 return domesticToGlobal_.find(domesticIdx)->second;
93 }
94
98 Index globalToDomestic(Index globalIdx) const
99 {
100 const auto& tmp = globalToDomestic_.find(globalIdx);
101
102 if (tmp == globalToDomestic_.end())
103 return -1;
104
105 return tmp->second;
106 }
107
112 size_t numLocal() const
113 { return foreignOverlap_.numLocal(); }
114
121 size_t numDomestic() const
122 { return numDomestic_; }
123
127 void addIndex(Index domesticIdx, Index globalIdx)
128 {
129 domesticToGlobal_[domesticIdx] = globalIdx;
130 globalToDomestic_[globalIdx] = domesticIdx;
132
133 assert(domesticToGlobal_.size() == globalToDomestic_.size());
134 }
135
139 void sendBorderIndex([[maybe_unused]] ProcessRank peerRank,
140 [[maybe_unused]] Index domesticIdx,
141 [[maybe_unused]] Index peerLocalIdx)
142 {
143#if HAVE_MPI
144 PeerIndexGlobalIndex sendBuf;
145 sendBuf.peerIdx = peerLocalIdx;
146 sendBuf.globalIdx = domesticToGlobal(domesticIdx);
147 MPI_Send(&sendBuf, // buff
148 sizeof(PeerIndexGlobalIndex), // count
149 MPI_BYTE, // data type
150 static_cast<int>(peerRank), // peer process
151 0, // tag
152 MPI_COMM_WORLD); // communicator
153#endif
154 }
155
160 void receiveBorderIndex([[maybe_unused]] ProcessRank peerRank)
161 {
162#if HAVE_MPI
163 PeerIndexGlobalIndex recvBuf;
164 MPI_Recv(&recvBuf, // buff
165 sizeof(PeerIndexGlobalIndex), // count
166 MPI_BYTE, // data type
167 static_cast<int>(peerRank), // peer process
168 0, // tag
169 MPI_COMM_WORLD, // communicator
170 MPI_STATUS_IGNORE); // status
171
172 Index domesticIdx = foreignOverlap_.nativeToLocal(recvBuf.peerIdx);
173 if (domesticIdx >= 0) {
174 Index globalIdx = recvBuf.globalIdx;
175 addIndex(domesticIdx, globalIdx);
176 }
177#endif // HAVE_MPI
178 }
179
183 bool hasGlobalIndex(Index globalIdx) const
184 { return globalToDomestic_.find(globalIdx) != globalToDomestic_.end(); }
185
190 void print() const
191 {
192 std::cout << "(domestic index, global index, domestic->global->domestic)"
193 << " list for rank " << myRank_ << "\n";
194
195 for (size_t domIdx = 0; domIdx < domesticToGlobal_.size(); ++domIdx)
196 std::cout << "(" << domIdx << ", " << domesticToGlobal(domIdx)
197 << ", " << globalToDomestic(domesticToGlobal(domIdx)) << ") ";
198 std::cout << "\n" << std::flush;
199 }
200
201protected:
202 // retrieve the offset for the indices where we are master in the
203 // global index list
205 {
206#if HAVE_MPI
207 numDomestic_ = 0;
208#else
209 numDomestic_ = foreignOverlap_.numLocal();
210#endif
211
212#if HAVE_MPI
213 if (myRank_ == 0) {
214 // the first rank starts at index zero
215 domesticOffset_ = 0;
216 }
217 else {
218 // all other ranks retrieve their offset from the next
219 // lower rank
220 MPI_Recv(&domesticOffset_, // buffer
221 1, // count
222 MPI_INT, // data type
223 static_cast<int>(myRank_ - 1), // peer rank
224 0, // tag
225 MPI_COMM_WORLD, // communicator
226 MPI_STATUS_IGNORE);
227 }
228
229 // create maps for all indices for which the current process
230 // is the master
231 int numMaster = 0;
232 for (unsigned i = 0; i < foreignOverlap_.numLocal(); ++i) {
233 if (!foreignOverlap_.iAmMasterOf(static_cast<Index>(i)))
234 continue;
235
236 addIndex(static_cast<Index>(i),
237 static_cast<Index>(domesticOffset_ + numMaster));
238 ++numMaster;
239 }
240
241 if (myRank_ < mpiSize_ - 1) {
242 // send the domestic offset plus the number of master
243 // indices to the process which is one rank higher
244 int tmp = domesticOffset_ + numMaster;
245 MPI_Send(&tmp, // buff
246 1, // count
247 MPI_INT, // data type
248 static_cast<int>(myRank_ + 1), // peer rank
249 0, // tag
250 MPI_COMM_WORLD); // communicator
251 }
252
253 typename PeerSet::const_iterator peerIt;
254 typename PeerSet::const_iterator peerEndIt = peerSet_().end();
255 // receive the border indices from the lower ranks
256 peerIt = peerSet_().begin();
257 for (; peerIt != peerEndIt; ++peerIt) {
258 if (*peerIt < myRank_)
259 receiveBorderFrom_(*peerIt);
260 }
261
262 // send the border indices to the higher ranks
263 peerIt = peerSet_().begin();
264 for (; peerIt != peerEndIt; ++peerIt) {
265 if (*peerIt > myRank_)
266 sendBorderTo_(*peerIt);
267 }
268
269 // receive the border indices from the higher ranks
270 peerIt = peerSet_().begin();
271 for (; peerIt != peerEndIt; ++peerIt) {
272 if (*peerIt > myRank_)
273 receiveBorderFrom_(*peerIt);
274 }
275
276 // send the border indices to the lower ranks
277 peerIt = peerSet_().begin();
278 for (; peerIt != peerEndIt; ++peerIt) {
279 if (*peerIt < myRank_)
280 sendBorderTo_(*peerIt);
281 }
282#endif // HAVE_MPI
283 }
284
285 void sendBorderTo_([[maybe_unused]] ProcessRank peerRank)
286 {
287#if HAVE_MPI
288 // send (local index on myRank, global index) pairs to the
289 // peers
290 BorderList::const_iterator borderIt = borderList_().begin();
291 BorderList::const_iterator borderEndIt = borderList_().end();
292 for (; borderIt != borderEndIt; ++borderIt) {
293 ProcessRank borderPeer = borderIt->peerRank;
294 BorderDistance borderDistance = borderIt->borderDistance;
295 if (borderPeer != peerRank || borderDistance != 0)
296 continue;
297
298 Index localIdx = foreignOverlap_.nativeToLocal(borderIt->localIdx);
299 Index peerIdx = borderIt->peerIdx;
300 assert(localIdx >= 0);
301 if (foreignOverlap_.iAmMasterOf(localIdx)) {
302 sendBorderIndex(borderPeer, localIdx, peerIdx);
303 }
304 }
305#endif // HAVE_MPI
306 }
307
308 void receiveBorderFrom_([[maybe_unused]] ProcessRank peerRank)
309 {
310#if HAVE_MPI
311 // retrieve the global indices for which we are not master
312 // from the processes with lower rank
313 BorderList::const_iterator borderIt = borderList_().begin();
314 BorderList::const_iterator borderEndIt = borderList_().end();
315 for (; borderIt != borderEndIt; ++borderIt) {
316 ProcessRank borderPeer = borderIt->peerRank;
317 BorderDistance borderDistance = borderIt->borderDistance;
318 if (borderPeer != peerRank || borderDistance != 0)
319 continue;
320
321 Index nativeIdx = borderIt->localIdx;
322 Index localIdx = foreignOverlap_.nativeToLocal(nativeIdx);
323 if (localIdx >= 0 && foreignOverlap_.masterRank(localIdx) == borderPeer)
324 receiveBorderIndex(borderPeer);
325 }
326#endif // HAVE_MPI
327 }
328
329 const PeerSet& peerSet_() const
330 { return foreignOverlap_.peerSet(); }
331
332 const BorderList& borderList_() const
333 { return foreignOverlap_.borderList(); }
334
336 size_t mpiSize_;
337
340 const ForeignOverlap& foreignOverlap_;
341
342 GlobalToDomesticMap globalToDomestic_;
343 DomesticToGlobalMap domesticToGlobal_;
344};
345
346} // namespace Linear
347} // namespace Opm
348
349#endif
This class maps domestic row indices to and from "global" indices which is used to construct an algeb...
Definition: globalindices.hh:56
void buildGlobalIndices_()
Definition: globalindices.hh:204
GlobalToDomesticMap globalToDomestic_
Definition: globalindices.hh:342
DomesticToGlobalMap domesticToGlobal_
Definition: globalindices.hh:343
int domesticOffset_
Definition: globalindices.hh:338
ProcessRank myRank_
Definition: globalindices.hh:335
void print() const
Prints the global indices of all domestic indices for debugging purposes.
Definition: globalindices.hh:190
const ForeignOverlap & foreignOverlap_
Definition: globalindices.hh:340
void addIndex(Index domesticIdx, Index globalIdx)
Add an index to the domestic<->global mapping.
Definition: globalindices.hh:127
size_t numDomestic() const
Returns the number domestic indices.
Definition: globalindices.hh:121
Index domesticToGlobal(Index domesticIdx) const
Converts a domestic index to a global one.
Definition: globalindices.hh:88
GlobalIndices(const ForeignOverlap &foreignOverlap)
Definition: globalindices.hh:63
size_t mpiSize_
Definition: globalindices.hh:336
void sendBorderIndex(ProcessRank peerRank, Index domesticIdx, Index peerLocalIdx)
Send a border index to a remote process.
Definition: globalindices.hh:139
void receiveBorderFrom_(ProcessRank peerRank)
Definition: globalindices.hh:308
void receiveBorderIndex(ProcessRank peerRank)
Receive an index on the border from a remote process and add it the translation maps.
Definition: globalindices.hh:160
Index globalToDomestic(Index globalIdx) const
Converts a global index to a domestic one.
Definition: globalindices.hh:98
size_t numDomestic_
Definition: globalindices.hh:339
void sendBorderTo_(ProcessRank peerRank)
Definition: globalindices.hh:285
const BorderList & borderList_() const
Definition: globalindices.hh:332
const PeerSet & peerSet_() const
Definition: globalindices.hh:329
bool hasGlobalIndex(Index globalIdx) const
Return true iff a given global index already exists.
Definition: globalindices.hh:183
size_t numLocal() const
Returns the number of indices which are in the interior or on the border of the current rank.
Definition: globalindices.hh:112
A set of process ranks.
Definition: overlaptypes.hh:149
unsigned BorderDistance
The type representing the distance of an index to the border.
Definition: overlaptypes.hh:54
unsigned ProcessRank
The type of the rank of a process.
Definition: overlaptypes.hh:49
int Index
The type of an index of a degree of freedom.
Definition: overlaptypes.hh:44
std::list< BorderIndex > BorderList
This class managages a list of indices which are on the border of a process' partition of the grid.
Definition: overlaptypes.hh:120
Definition: blackoilboundaryratevector.hh:37
This files provides several data structures for storing tuples of indices of remote and/or local proc...
This structure stores a local index on a peer process and a global index.
Definition: overlaptypes.hh:70
Index globalIdx
Definition: overlaptypes.hh:72
Index peerIdx
Definition: overlaptypes.hh:71