foreignoverlapfrombcrsmatrix.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_FOREIGN_OVERLAP_FROM_BCRS_MATRIX_HH
28#define EWOMS_FOREIGN_OVERLAP_FROM_BCRS_MATRIX_HH
29
30#include "overlaptypes.hh"
31#include "blacklist.hh"
32
34
35#include <dune/grid/common/datahandleif.hh>
36#include <dune/istl/bcrsmatrix.hh>
37#include <dune/istl/scalarproducts.hh>
38#include <dune/istl/operators.hh>
39
40#include <algorithm>
41#include <iostream>
42#include <map>
43#include <vector>
44
45#if HAVE_MPI
46#include <mpi.h>
47#endif // HAVE_MPI
48
49namespace Opm {
50namespace Linear {
51
60{
61public:
62 // overlaps should never be copied!
64
69 template <class BCRSMatrix>
70 ForeignOverlapFromBCRSMatrix(const BCRSMatrix& A,
72 const BlackList& blackList,
73 unsigned overlapSize)
75 {
77
78 myRank_ = 0;
79#if HAVE_MPI
80 {
81 int tmp;
82 MPI_Comm_rank(MPI_COMM_WORLD, &tmp);
83 myRank_ = static_cast<ProcessRank>(tmp);
84 }
85#endif
86 numNative_ = A.N();
87
88 // Computes the local <-> native index maps
90
91 // calculate the set of local indices on the border (beware:
92 // _not_ the native ones)
93 auto it = borderList.begin();
94 const auto& endIt = borderList.end();
95 for (; it != endIt; ++it) {
96 Index localIdx = nativeToLocal(it->localIdx);
97 if (localIdx < 0)
98 continue;
99
100 localBorderIndices_.insert(localIdx);
101 }
102
103 // compute the set of processes which are neighbors of the
104 // local process ...
106 // ... and the initial set of processes which we will have to
107 // communicate with. We must always communicate with our
108 // neighbors, but depending on the size of the overlap region,
109 // we might have to communicate with additional processes as
110 // well (these will be added later).
112
113 // Create an initial seed list of indices which are in the
114 // overlap.
115 SeedList initialSeedList;
116 initialSeedList.update(borderList);
117
118 // calculate the minimum distance from the border of the
119 // initial seed list
120 unsigned minBorderDist = overlapSize;
121 auto borderIt = borderList.begin();
122 const auto& borderEndIt = borderList.end();
123 for (; borderIt != borderEndIt; ++borderIt) {
124 minBorderDist = std::min(minBorderDist, borderIt->borderDistance);
125 }
126
127 // calculate the foreign overlap for the local partition,
128 // i.e. find the distance of each row from the seed set.
130 extendForeignOverlap_(A, initialSeedList, minBorderDist, overlapSize);
131
132 // computes the process with the lowest rank for all local
133 // indices.
135
136 // group foreign overlap by peer process rank
138 }
139
143 unsigned overlapSize() const
144 { return overlapSize_; }
145
149 bool isBorder(Index localIdx) const
150 { return localBorderIndices_.count(localIdx) > 0; }
151
156 bool isBorderWith(Index localIdx, ProcessRank peerRank) const
157 {
158 const auto& indexOverlap = foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)];
159 const auto& borderDistIt = indexOverlap.find(peerRank);
160 if (borderDistIt == indexOverlap.end())
161 return false;
162
163 // border distance of the index needs to be 0
164 return borderDistIt->second == 0;
165 }
166
172 { return masterRank_[static_cast<unsigned>(localIdx)]; }
173
182 bool iAmMasterOf(Index localIdx) const
183 { return masterRank_[static_cast<unsigned>(localIdx)] == myRank_; }
184
189 const BorderList& borderList() const
190 { return borderList_; }
191
198 {
199 assert(foreignOverlapByRank_.find(peerRank) != foreignOverlapByRank_.end());
200 return foreignOverlapByRank_.find(peerRank)->second;
201 }
202
207 const std::map<ProcessRank, BorderDistance> &
209 {
210 assert(isLocal(localIdx));
211 return foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)];
212 }
213
217 bool peerHasIndex(ProcessRank peerRank, Index localIdx) const
218 {
219 const auto& idxOverlap = foreignOverlapByLocalIndex_[localIdx];
220 return idxOverlap.find(peerRank) != idxOverlap.end();
221 }
222
227 size_t numFront(ProcessRank peerRank) const
228 {
229 const auto& peerOverlap = foreignOverlapByRank_.find(peerRank)->second;
230
231 size_t n = 0;
232 auto it = peerOverlap.begin();
233 const auto& endIt = peerOverlap.end();
234 for (; it != endIt; ++it) {
235 if (it->borderDistance == overlapSize_)
236 ++n;
237 }
238 return n;
239 }
240
245 bool isFrontFor(ProcessRank peerRank, Index localIdx) const
246 {
247 const auto& idxOverlap = foreignOverlapByLocalIndex_[localIdx];
248
249 auto it = idxOverlap.find(peerRank);
250 if (it == idxOverlap.end())
251 return false; // index is not in overlap
252
253 return it->second == overlapSize_;
254 }
255
260 const PeerSet& peerSet() const
261 { return peerSet_; }
262
268 { return neighborPeerSet_; }
269
273 size_t numNative() const
274 { return numNative_; }
275
279 size_t numLocal() const
280 { return numLocal_; }
281
285 bool isLocal(Index domesticIdx) const
286 { return static_cast<unsigned>(domesticIdx) < numLocal(); }
287
294 Index nativeToLocal(Index nativeIdx) const
295 { return nativeToLocalIndices_[static_cast<unsigned>(nativeIdx)]; }
296
300 Index localToNative(Index localIdx) const
301 {
302 assert(localIdx < static_cast<Index>(localToNativeIndices_.size()));
303 return localToNativeIndices_[static_cast<unsigned>(localIdx)];
304 }
305
309 const BlackList& blackList() const
310 { return blackList_; }
311
316 size_t numPeers(Index localIdx) const
317 { return foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].size(); }
318
323 bool isInOverlap(Index localIdx) const
324 { return foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].size() > 0; }
325
329 void print() const
330 {
331 auto it = foreignOverlapByRank_.begin();
332 const auto& endIt = foreignOverlapByRank_.end();
333 for (; it != endIt; ++it) {
334 std::cout << "Overlap rows(distance) for rank " << it->first << ": ";
335
336 auto rowIt = it->second.begin();
337 const auto& rowEndIt = it->second.end();
338 for (; rowIt != rowEndIt; ++rowIt)
339 std::cout << rowIt->index << "(" << rowIt->borderDistance << ") ";
340 std::cout << "\n" << std::flush;
341 }
342 }
343
344protected:
345 // extend the foreign overlaps by 'overlapSize' levels. this uses
346 // a greedy algorithm which extends the region by one level and
347 // then calls itself recursively...
348 template <class BCRSMatrix>
349 void extendForeignOverlap_(const BCRSMatrix& A,
350 SeedList& seedList,
351 BorderDistance borderDistance,
353 {
354 // communicate the non-neigbor overlap indices
355 addNonNeighborOverlapIndices_(A, seedList, borderDistance);
356
357 // add all processes in the seed rows of the current overlap level
358 auto seedIt = seedList.begin();
359 const auto& seedEndIt = seedList.end();
360 for (; seedIt != seedEndIt; ++seedIt) {
361 Index localIdx = nativeToLocal(seedIt->index);
362 ProcessRank peerRank = seedIt->peerRank;
363 unsigned distance = borderDistance;
364 if (localIdx < 0)
365 continue;
366 if (foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].count(peerRank) == 0)
367 foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)][peerRank] = distance;
368 }
369
370 // if we have reached the maximum overlap distance, i.e. we're
371 // finished and break the recursion
372 if (borderDistance >= overlapSize)
373 return;
374
375 // find the seed list for the next overlap level using the
376 // seed set for the current level
377 SeedList nextSeedList;
378 seedIt = seedList.begin();
379 for (; seedIt != seedEndIt; ++seedIt) {
380 Index nativeRowIdx = seedIt->index;
381 if (nativeToLocal(nativeRowIdx) < 0)
382 continue; // ignore blacklisted indices
383 ProcessRank peerRank = seedIt->peerRank;
384
385 // find all column indices in the row. The indices of the
386 // columns are the additional indices of the overlap which
387 // we would like to add
388 using ColIterator = typename BCRSMatrix::ConstColIterator;
389 ColIterator colIt = A[static_cast<unsigned>(nativeRowIdx)].begin();
390 ColIterator colEndIt = A[static_cast<unsigned>(nativeRowIdx)].end();
391 for (; colIt != colEndIt; ++colIt) {
392 Index nativeColIdx = static_cast<Index>(colIt.index());
393 Index localColIdx = nativeToLocal(nativeColIdx);
394
395 // ignore if the native index is not a local one
396 if (localColIdx < 0)
397 continue;
398 // if the process is already is in the overlap of the
399 // column index, ignore this column index!
400 else if (foreignOverlapByLocalIndex_[static_cast<unsigned>(localColIdx)].count(peerRank) > 0)
401 continue;
402
403 // check whether the new index is already in the overlap
404 bool hasIndex = false;
405 typename SeedList::iterator sIt = nextSeedList.begin();
406 typename SeedList::iterator sEndIt = nextSeedList.end();
407 for (; sIt != sEndIt; ++sIt) {
408 if (sIt->index == nativeColIdx && sIt->peerRank == peerRank) {
409 hasIndex = true;
410 break;
411 }
412 }
413 if (hasIndex)
414 continue; // we already have this index
415
416 // add the current processes to the seed list for the
417 // next overlap level
418 IndexRankDist newTuple;
419 newTuple.index = nativeColIdx;
420 newTuple.peerRank = peerRank;
421 newTuple.borderDistance = seedIt->borderDistance + 1;
422 nextSeedList.push_back(newTuple);
423 }
424 }
425
426 // clear the old seed list to save some memory
427 seedList.clear();
428
429 // Perform the same excercise for the next overlap distance
430 extendForeignOverlap_(A, nextSeedList, borderDistance + 1, overlapSize);
431 }
432
433 // Computes the local <-> native index maps
435 {
436 // create the native <-> local maps
437 Index localIdx = 0;
438 for (unsigned nativeIdx = 0; nativeIdx < numNative_;) {
439 if (!blackList_.hasIndex(static_cast<Index>(nativeIdx))) {
440 localToNativeIndices_.push_back(static_cast<Index>(nativeIdx));
441 nativeToLocalIndices_.push_back(static_cast<Index>(localIdx));
442 ++nativeIdx;
443 ++localIdx;
444 }
445 else {
446 nativeToLocalIndices_.push_back(-1);
447 ++nativeIdx;
448 }
449 }
450
452 }
453
454 Index localToPeerIdx_(Index localIdx, ProcessRank peerRank) const
455 {
456 auto it = borderList_.begin();
457 const auto& endIt = borderList_.end();
458 for (; it != endIt; ++it) {
459 if (it->localIdx == localIdx && it->peerRank == peerRank)
460 return it->peerIdx;
461 }
462
463 return -1;
464 }
465
466 template <class BCRSMatrix>
467 void addNonNeighborOverlapIndices_(const BCRSMatrix&,
468 [[maybe_unused]] SeedList& seedList,
469 [[maybe_unused]] BorderDistance borderDist)
470 {
471 // TODO: this probably does not work! (the matrix A is unused, but it is needed
472 // from a logical POV.)
473#if HAVE_MPI
474 // first, create the buffers which will contain the number of
475 // border indices relevant for a neighbor peer
476 std::map<ProcessRank, std::vector<BorderIndex> > borderIndices;
477
478 // get all indices in the border which have borderDist as
479 // their distance to the closest border of their local process
480 auto it = seedList.begin();
481 const auto& endIt = seedList.end();
482 for (; it != endIt; ++it) {
483 Index localIdx = nativeToLocal(it->index);
484 if (!isBorder(localIdx))
485 continue;
486 BorderIndex borderHandle;
487 borderHandle.localIdx = localIdx;
488 borderHandle.peerRank = it->peerRank;
489 borderHandle.borderDistance = it->borderDistance;
490
491 // add the border index to all the neighboring peers
492 auto neighborIt = foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].begin();
493 const auto& neighborEndIt = foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].end();
494 for (; neighborIt != neighborEndIt; ++neighborIt) {
495 if (neighborIt->second != 0)
496 // not a border index for the neighbor
497 continue;
498 else if (neighborIt->first == borderHandle.peerRank)
499 // don't communicate the indices which are owned
500 // by the peer to itself
501 continue;
502
503 Index peerIdx = localToPeerIdx_(localIdx, neighborIt->first);
504 if (peerIdx < 0)
505 // the index is on the border, but is not on the border
506 // with the considered neighboring process. Ignore it!
507 continue;
508 borderHandle.peerIdx = peerIdx;
509 borderIndices[neighborIt->first].push_back(borderHandle);
510 }
511 }
512
513 // now borderIndices contains the lists of indices which we
514 // would like to send to each neighbor. Let's create the MPI
515 // buffers.
516 std::map<ProcessRank, Opm::MpiBuffer<unsigned> > numIndicesSendBufs;
517 std::map<ProcessRank, Opm::MpiBuffer<BorderIndex> > indicesSendBufs;
518 auto peerIt = neighborPeerSet().begin();
519 const auto& peerEndIt = neighborPeerSet().end();
520 for (; peerIt != peerEndIt; ++peerIt) {
521 ProcessRank peerRank = *peerIt;
522 size_t numIndices = borderIndices[peerRank].size();
523 numIndicesSendBufs[peerRank].resize(1);
524 numIndicesSendBufs[peerRank][0] = static_cast<unsigned>(numIndices);
525
526 const auto& peerBorderIndices = borderIndices[peerRank];
527 indicesSendBufs[peerRank].resize(numIndices);
528
529 auto tmpIt = peerBorderIndices.begin();
530 const auto& tmpEndIt = peerBorderIndices.end();
531 size_t i = 0;
532 for (; tmpIt != tmpEndIt; ++tmpIt, ++i) {
533 indicesSendBufs[peerRank][i] = *tmpIt;
534 }
535 }
536
537 // now, send all these nice buffers to our neighbors
538 peerIt = neighborPeerSet().begin();
539 for (; peerIt != peerEndIt; ++peerIt) {
540 ProcessRank neighborPeer = *peerIt;
541 numIndicesSendBufs[neighborPeer].send(neighborPeer);
542 indicesSendBufs[neighborPeer].send(neighborPeer);
543 }
544
545 // receive all data from the neighbors
546 std::map<ProcessRank, MpiBuffer<unsigned> > numIndicesRcvBufs;
547 std::map<ProcessRank, MpiBuffer<BorderIndex> > indicesRcvBufs;
548 peerIt = neighborPeerSet().begin();
549 for (; peerIt != peerEndIt; ++peerIt) {
550 ProcessRank neighborPeer = *peerIt;
551 auto& numIndicesRcvBuf = numIndicesRcvBufs[neighborPeer];
552 auto& indicesRcvBuf = indicesRcvBufs[neighborPeer];
553
554 numIndicesRcvBuf.resize(1);
555 numIndicesRcvBuf.receive(neighborPeer);
556 unsigned numIndices = numIndicesRcvBufs[neighborPeer][0];
557 indicesRcvBuf.resize(numIndices);
558 indicesRcvBuf.receive(neighborPeer);
559
560 // filter out all indices which are already in the peer
561 // processes' overlap and add them to the seed list. also
562 // extend the set of peer processes.
563 for (unsigned i = 0; i < numIndices; ++i) {
564 // swap the local and the peer indices, because they were
565 // created with the point view of the sender
566 std::swap(indicesRcvBuf[i].localIdx, indicesRcvBuf[i].peerIdx);
567
568 ProcessRank peerRank = indicesRcvBuf[i].peerRank;
569 // Index peerIdx = indicesRcvBuf[i].peerIdx;
570 Index localIdx = indicesRcvBuf[i].localIdx;
571
572 // check if the index is already in the overlap for
573 // the peer
574 const auto& distIt = foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].find(peerRank);
575 if (distIt != foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].end())
576 continue;
577
578 // make sure the index is not already in the seed list
579 bool inSeedList = false;
580 auto seedIt = seedList.begin();
581 const auto& seedEndIt = seedList.end();
582 for (; seedIt != seedEndIt; ++seedIt) {
583 if (seedIt->index == localIdx && seedIt->peerRank == peerRank) {
584 inSeedList = true;
585 break;
586 }
587 }
588 if (inSeedList)
589 continue;
590
591 IndexRankDist seedEntry;
592 seedEntry.index = localIdx;
593 seedEntry.peerRank = peerRank;
594 seedEntry.borderDistance = borderDist;
595 seedList.push_back(seedEntry);
596
597 // update the peer set
598 peerSet_.insert(peerRank);
599 }
600 }
601
602 // make sure all data was send
603 peerIt = neighborPeerSet().begin();
604 for (; peerIt != peerEndIt; ++peerIt) {
605 ProcessRank neighborPeer = *peerIt;
606 numIndicesSendBufs[neighborPeer].wait();
607 indicesSendBufs[neighborPeer].wait();
608 }
609#endif // HAVE_MPI
610 }
611
612 // given a list of border indices and provided that
613 // borderListToSeedList_() was already called, calculate the
614 // master process of each local index.
616 {
617 // determine the minimum rank for all indices
618 masterRank_.resize(numLocal_);
619 for (unsigned localIdx = 0; localIdx < numLocal_; ++localIdx) {
620 unsigned masterRank = myRank_;
621 if (isBorder(static_cast<Index>(localIdx))) {
622 // if the local index is a border index, loop over all ranks
623 // for which this index is also a border index. the lowest
624 // rank wins!
625 auto it = foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].begin();
626 const auto& endIt = foreignOverlapByLocalIndex_[static_cast<unsigned>(localIdx)].end();
627 for (; it != endIt; ++it) {
628 if (it->second == 0) {
629 // if the border distance is zero, the rank with the
630 // minimum
631 masterRank = std::min<ProcessRank>(masterRank, it->first);
632 }
633 }
634 }
635 masterRank_[static_cast<unsigned>(localIdx)] = masterRank;
636 }
637 }
638
639 // assuming that the foreign overlap has been created for each
640 // local index, this method groups the foreign overlap by peer
641 // process rank
643 {
644 // loop over all indices which are in the overlap of some
645 // process
647 for (unsigned localIdx = 0; localIdx < numLocal; ++localIdx) {
648 // loop over the list of processes for the current index
649 auto it = foreignOverlapByLocalIndex_[localIdx].begin();
650 const auto& endIt = foreignOverlapByLocalIndex_[localIdx].end();
651 size_t nRanks = foreignOverlapByLocalIndex_[localIdx].size();
652 for (; it != endIt; ++it) {
654 tmp.index = static_cast<Index>(localIdx);
655 tmp.borderDistance = it->second;
656 tmp.numPeers = static_cast<unsigned>(nRanks);
657 foreignOverlapByRank_[it->first].push_back(tmp);
658 }
659 }
660 }
661
662 // set of processes with which we have to communicate
664
665 // set of processes which are direct neighbors of us
667
668 // the list of indices on the border
670
671 // the set of indices which should not be considered
673
674 // local indices are the native indices sans the black listed ones
675 std::vector<Index> nativeToLocalIndices_;
676 std::vector<Index> localToNativeIndices_;
677
678 // an array which contains the rank of the master process for each
679 // index
680 std::vector<ProcessRank> masterRank_;
681
682 // set of all local indices which are on the border of some remote
683 // process
684 std::set<Index> localBorderIndices_;
685
686 // stores the set of process ranks which are in the overlap for a
687 // given row index "owned" by the current rank. The second value
688 // store the distance from the nearest process border.
690
691 // stores a list of foreign overlap indices for each rank
693
694 // size of the overlap region
695 unsigned overlapSize_;
696
697 // number of local indices
698 size_t numLocal_;
699
700 // number of native indices
702
703 // the MPI rank of the local process
705};
706
707} // namespace Linear
708} // namespace Opm
709
710#endif
Expresses which degrees of freedom are blacklisted for the parallel linear solvers and which domestic...
Definition: blacklist.hh:49
bool hasIndex(Index nativeIdx) const
Definition: blacklist.hh:63
This class creates and manages the foreign overlap given an initial list of border indices and a BCRS...
Definition: foreignoverlapfrombcrsmatrix.hh:60
std::vector< Index > localToNativeIndices_
Definition: foreignoverlapfrombcrsmatrix.hh:676
Index localToPeerIdx_(Index localIdx, ProcessRank peerRank) const
Definition: foreignoverlapfrombcrsmatrix.hh:454
const PeerSet & peerSet() const
Return the set of process ranks which share an overlap with the current process.
Definition: foreignoverlapfrombcrsmatrix.hh:260
size_t numLocal() const
Returns the number of local indices.
Definition: foreignoverlapfrombcrsmatrix.hh:279
size_t numNative_
Definition: foreignoverlapfrombcrsmatrix.hh:701
bool isInOverlap(Index localIdx) const
Returns true if a given local index is in the foreign overlap of any rank.
Definition: foreignoverlapfrombcrsmatrix.hh:323
void computeMasterRanks_()
Definition: foreignoverlapfrombcrsmatrix.hh:615
const OverlapWithPeer & foreignOverlapWithPeer(ProcessRank peerRank) const
Return the list of (local indices, border distance, number of processes) triples which are in the ove...
Definition: foreignoverlapfrombcrsmatrix.hh:197
void extendForeignOverlap_(const BCRSMatrix &A, SeedList &seedList, BorderDistance borderDistance, BorderDistance overlapSize)
Definition: foreignoverlapfrombcrsmatrix.hh:349
std::vector< Index > nativeToLocalIndices_
Definition: foreignoverlapfrombcrsmatrix.hh:675
void print() const
Print the foreign overlap for debugging purposes.
Definition: foreignoverlapfrombcrsmatrix.hh:329
unsigned overlapSize() const
Returns the size of the overlap region.
Definition: foreignoverlapfrombcrsmatrix.hh:143
bool isBorder(Index localIdx) const
Returns true iff a local index is a border index.
Definition: foreignoverlapfrombcrsmatrix.hh:149
std::set< Index > localBorderIndices_
Definition: foreignoverlapfrombcrsmatrix.hh:684
Index nativeToLocal(Index nativeIdx) const
Convert a native index to a local one.
Definition: foreignoverlapfrombcrsmatrix.hh:294
const PeerSet & neighborPeerSet() const
Return the set of process ranks which share a border index with the current process.
Definition: foreignoverlapfrombcrsmatrix.hh:267
ForeignOverlapFromBCRSMatrix(const BCRSMatrix &A, const BorderList &borderList, const BlackList &blackList, unsigned overlapSize)
Constructs the foreign overlap given a BCRS matrix and an initial list of border indices.
Definition: foreignoverlapfrombcrsmatrix.hh:70
PeerSet neighborPeerSet_
Definition: foreignoverlapfrombcrsmatrix.hh:666
bool peerHasIndex(ProcessRank peerRank, Index localIdx) const
Returns true iff a local index is seen by a peer rank.
Definition: foreignoverlapfrombcrsmatrix.hh:217
PeerSet peerSet_
Definition: foreignoverlapfrombcrsmatrix.hh:663
size_t numNative() const
Returns the number of native indices.
Definition: foreignoverlapfrombcrsmatrix.hh:273
void createLocalIndices_()
Definition: foreignoverlapfrombcrsmatrix.hh:434
OverlapByIndex foreignOverlapByLocalIndex_
Definition: foreignoverlapfrombcrsmatrix.hh:689
ProcessRank myRank_
Definition: foreignoverlapfrombcrsmatrix.hh:704
const std::map< ProcessRank, BorderDistance > & foreignOverlapByLocalIndex(Index localIdx) const
Return the map of (peer rank, border distance) for a given local index.
Definition: foreignoverlapfrombcrsmatrix.hh:208
bool isFrontFor(ProcessRank peerRank, Index localIdx) const
Returns whether a given local index is on the front of a given peer rank.
Definition: foreignoverlapfrombcrsmatrix.hh:245
const BlackList & blackList_
Definition: foreignoverlapfrombcrsmatrix.hh:672
unsigned overlapSize_
Definition: foreignoverlapfrombcrsmatrix.hh:695
ForeignOverlapFromBCRSMatrix(const ForeignOverlapFromBCRSMatrix &)=delete
const BorderList & borderList() const
Returns the list of indices which intersect the process border.
Definition: foreignoverlapfrombcrsmatrix.hh:189
size_t numFront(ProcessRank peerRank) const
Returns the number of front indices of a peer process in the local partition.
Definition: foreignoverlapfrombcrsmatrix.hh:227
const BorderList & borderList_
Definition: foreignoverlapfrombcrsmatrix.hh:669
Index localToNative(Index localIdx) const
Convert a local index to a native one.
Definition: foreignoverlapfrombcrsmatrix.hh:300
void groupForeignOverlapByRank_()
Definition: foreignoverlapfrombcrsmatrix.hh:642
size_t numLocal_
Definition: foreignoverlapfrombcrsmatrix.hh:698
size_t numPeers(Index localIdx) const
Return the number of peer ranks for which a given local index is visible.
Definition: foreignoverlapfrombcrsmatrix.hh:316
ProcessRank masterRank(Index localIdx) const
Return the rank of the master process of an index.
Definition: foreignoverlapfrombcrsmatrix.hh:171
const BlackList & blackList() const
Returns the object which represents the black-listed native indices.
Definition: foreignoverlapfrombcrsmatrix.hh:309
bool isBorderWith(Index localIdx, ProcessRank peerRank) const
Returns true iff a local index is a border index shared with a given peer process.
Definition: foreignoverlapfrombcrsmatrix.hh:156
bool iAmMasterOf(Index localIdx) const
Return true if the current rank is the "master" of an index.
Definition: foreignoverlapfrombcrsmatrix.hh:182
std::vector< ProcessRank > masterRank_
Definition: foreignoverlapfrombcrsmatrix.hh:680
void addNonNeighborOverlapIndices_(const BCRSMatrix &, SeedList &seedList, BorderDistance borderDist)
Definition: foreignoverlapfrombcrsmatrix.hh:467
bool isLocal(Index domesticIdx) const
Returns true iff a domestic index is local.
Definition: foreignoverlapfrombcrsmatrix.hh:285
OverlapByRank foreignOverlapByRank_
Definition: foreignoverlapfrombcrsmatrix.hh:692
A set of process ranks.
Definition: overlaptypes.hh:149
void update(const BorderList &borderList)
Definition: overlaptypes.hh:151
The list of indices which are on the process boundary.
Definition: overlaptypes.hh:126
void update(const BorderList &borderList)
Definition: overlaptypes.hh:128
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
std::vector< std::map< ProcessRank, BorderDistance > > OverlapByIndex
Maps each index to a list of processes .
Definition: overlaptypes.hh:176
std::vector< IndexDistanceNpeers > OverlapWithPeer
The list of indices which overlap with a peer rank.
Definition: overlaptypes.hh:165
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
std::map< ProcessRank, OverlapWithPeer > OverlapByRank
A type mapping the process rank to the list of indices shared with this peer.
Definition: overlaptypes.hh:171
Definition: blackoilboundaryratevector.hh:37
This files provides several data structures for storing tuples of indices of remote and/or local proc...
A single index intersecting with the process boundary.
Definition: overlaptypes.hh:102
Index localIdx
Index of the entity for the local process.
Definition: overlaptypes.hh:104
BorderDistance borderDistance
Distance to the process border for the peer (in hops)
Definition: overlaptypes.hh:113
ProcessRank peerRank
Rank of the peer process.
Definition: overlaptypes.hh:110
Index peerIdx
Index of the entity for the peer process.
Definition: overlaptypes.hh:107
This structure stores an index, a process rank, and the number of processes which "see" the degree of...
Definition: overlaptypes.hh:92
BorderDistance borderDistance
Definition: overlaptypes.hh:94
unsigned numPeers
Definition: overlaptypes.hh:95
Index index
Definition: overlaptypes.hh:93
This structure stores an index, a process rank, and the distance of the degree of freedom to the proc...
Definition: overlaptypes.hh:80
BorderDistance borderDistance
Definition: overlaptypes.hh:83
ProcessRank peerRank
Definition: overlaptypes.hh:82
Index index
Definition: overlaptypes.hh:81