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