MueLu  Version of the Day
MueLu_RepartitionFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_REPARTITIONFACTORY_DEF_HPP
47 #define MUELU_REPARTITIONFACTORY_DEF_HPP
48 
49 #include <algorithm>
50 #include <iostream>
51 #include <sstream>
52 
53 #include "MueLu_RepartitionFactory_decl.hpp" // TMP JG NOTE: before other includes, otherwise I cannot test the fwd declaration in _def
54 
55 #ifdef HAVE_MPI
56 #include <Teuchos_DefaultMpiComm.hpp>
57 #include <Teuchos_CommHelpers.hpp>
58 
59 #include <Xpetra_Map.hpp>
60 #include <Xpetra_MapFactory.hpp>
61 #include <Xpetra_VectorFactory.hpp>
62 #include <Xpetra_Import.hpp>
63 #include <Xpetra_ImportFactory.hpp>
64 #include <Xpetra_Export.hpp>
65 #include <Xpetra_ExportFactory.hpp>
66 #include <Xpetra_Matrix.hpp>
67 #include <Xpetra_MatrixFactory.hpp>
68 
69 #include "MueLu_Utilities.hpp"
70 
71 #include "MueLu_CloneRepartitionInterface.hpp"
72 
73 #include "MueLu_Level.hpp"
74 #include "MueLu_MasterList.hpp"
75 #include "MueLu_Monitor.hpp"
76 
77 namespace MueLu {
78 
79  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81  RCP<ParameterList> validParamList = rcp(new ParameterList());
82 
83 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
84  SET_VALID_ENTRY("repartition: print partition distribution");
85  SET_VALID_ENTRY("repartition: remap parts");
86  SET_VALID_ENTRY("repartition: remap num values");
87 #undef SET_VALID_ENTRY
88 
89  validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Factory of the matrix A");
90  validParamList->set< RCP<const FactoryBase> >("number of partitions", Teuchos::null, "Instance of RepartitionHeuristicFactory.");
91  validParamList->set< RCP<const FactoryBase> >("Partition", Teuchos::null, "Factory of the partition");
92 
93  return validParamList;
94  }
95 
96  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98  Input(currentLevel, "A");
99  Input(currentLevel, "number of partitions");
100  Input(currentLevel, "Partition");
101  }
102 
103  template<class T> class MpiTypeTraits { public: static MPI_Datatype getType(); };
104  template<> class MpiTypeTraits<long> { public: static MPI_Datatype getType() { return MPI_LONG; } };
105  template<> class MpiTypeTraits<int> { public: static MPI_Datatype getType() { return MPI_INT; } };
106  template<> class MpiTypeTraits<short> { public: static MPI_Datatype getType() { return MPI_SHORT; } };
107  template<> class MpiTypeTraits<unsigned> { public: static MPI_Datatype getType() { return MPI_UNSIGNED; } };
108  template<> class MpiTypeTraits<long long> { public: static MPI_Datatype getType() { return MPI_LONG_LONG; } };
109 
110  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
112  FactoryMonitor m(*this, "Build", currentLevel);
113 
114  const Teuchos::ParameterList & pL = GetParameterList();
115  // Access parameters here to make sure that we set the parameter entry flag to "used" even in case of short-circuit evaluation.
116  // TODO (JG): I don't really know if we want to do this.
117  bool remapPartitions = pL.get<bool> ("repartition: remap parts");
118 
119  // TODO: We only need a CrsGraph. This class does not have to be templated on Scalar types.
120  RCP<Matrix> A = Get< RCP<Matrix> >(currentLevel, "A");
121  RCP<const Map> rowMap = A->getRowMap();
122  GO indexBase = rowMap->getIndexBase();
123  Xpetra::UnderlyingLib lib = rowMap->lib();
124 
125  // NOTE: Teuchos::MPIComm::duplicate() calls MPI_Bcast inside, so this is
126  // a synchronization point. However, as we do MueLu_sumAll afterwards anyway, it
127  // does not matter.
128  RCP<const Teuchos::Comm<int> > origComm = rowMap->getComm();
129  RCP<const Teuchos::Comm<int> > comm = origComm->duplicate();
130 
131  int myRank = comm->getRank();
132  int numProcs = comm->getSize();
133 
134  RCP<const Teuchos::MpiComm<int> > tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
135  TEUCHOS_TEST_FOR_EXCEPTION(tmpic == Teuchos::null, Exceptions::RuntimeError, "Cannot cast base Teuchos::Comm to Teuchos::MpiComm object.");
136  RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawMpiComm = tmpic->getRawMpiComm();
137 
139  int numPartitions = Get<int>(currentLevel, "number of partitions");
140 
141  // ======================================================================================================
142  // Construct decomposition vector
143  // ======================================================================================================
144  RCP<GOVector> decomposition = Get<RCP<GOVector> >(currentLevel, "Partition");
145 
146  // check which factory provides "Partition"
147  if(remapPartitions == true && Teuchos::rcp_dynamic_cast<const CloneRepartitionInterface>(GetFactory("Partition")) != Teuchos::null) {
148  // if "Partition" is provided by a CloneRepartitionInterface class we have to switch of remapPartitions
149  // as we can assume the processor Ids in Partition to be the expected ones. If we would do remapping we
150  // would get different processors for the different blocks which screws up matrix-matrix multiplication.
151  remapPartitions = false;
152  }
153 
154  // check special cases
155  if (numPartitions == 1) {
156  // Trivial case: decomposition is the trivial one, all zeros. We skip the call to Zoltan_Interface
157  // (this is mostly done to avoid extra output messages, as even if we didn't skip there is a shortcut
158  // in Zoltan[12]Interface).
159  // TODO: We can probably skip more work in this case (like building all extra data structures)
160  GetOStream(Warnings0) << "Only one partition: Skip call to the repartitioner." << std::endl;
161  } else if (numPartitions == -1) {
162  // No repartitioning necessary: decomposition should be Teuchos::null
163  GetOStream(Warnings0) << "No repartitioning necessary: partitions were left unchanged by the repartitioner" << std::endl;
164  Set<RCP<const Import> >(currentLevel, "Importer", Teuchos::null);
165  return;
166  }
167 
168  // ======================================================================================================
169  // Remap if necessary
170  // ======================================================================================================
171  // From a user perspective, we want user to not care about remapping, thinking of it as only a performance feature.
172  // There are two problems, however.
173  // (1) Next level aggregation depends on the order of GIDs in the vector, if one uses "natural" or "random" orderings.
174  // This also means that remapping affects next level aggregation, despite the fact that the _set_ of GIDs for
175  // each partition is the same.
176  // (2) Even with the fixed order of GIDs, the remapping may influence the aggregation for the next-next level.
177  // Let us consider the following example. Lets assume that when we don't do remapping, processor 0 would have
178  // GIDs {0,1,2}, and processor 1 GIDs {3,4,5}, and if we do remapping processor 0 would contain {3,4,5} and
179  // processor 1 {0,1,2}. Now, when we run repartitioning algorithm on the next level (say Zoltan1 RCB), it may
180  // be dependent on whether whether it is [{0,1,2}, {3,4,5}] or [{3,4,5}, {0,1,2}]. Specifically, the tie-breaking
181  // algorithm can resolve these differently. For instance, running
182  // mpirun -np 5 ./MueLu_ScalingTestParamList.exe --xml=easy_sa.xml --nx=12 --ny=12 --nz=12
183  // with
184  // <ParameterList name="MueLu">
185  // <Parameter name="coarse: max size" type="int" value="1"/>
186  // <Parameter name="repartition: enable" type="bool" value="true"/>
187  // <Parameter name="repartition: min rows per proc" type="int" value="2"/>
188  // <ParameterList name="level 1">
189  // <Parameter name="repartition: remap parts" type="bool" value="false/true"/>
190  // </ParameterList>
191  // </ParameterList>
192  // produces different repartitioning for level 2.
193  // This different repartitioning may then escalate into different aggregation for the next level.
194  //
195  // We fix (1) by fixing the order of GIDs in a vector by sorting the resulting vector.
196  // Fixing (2) is more complicated.
197  // FIXME: Fixing (2) in Zoltan may not be enough, as we may use some arbitration in MueLu,
198  // for instance with CoupledAggregation. What we really need to do is to use the same order of processors containing
199  // the same order of GIDs. To achieve that, the newly created subcommunicator must be conforming with the order. For
200  // instance, if we have [{0,1,2}, {3,4,5}], we create a subcommunicator where processor 0 gets rank 0, and processor 1
201  // gets rank 1. If, on the other hand, we have [{3,4,5}, {0,1,2}], we assign rank 1 to processor 0, and rank 0 to processor 1.
202  // This rank permutation requires help from Epetra/Tpetra, both of which have no such API in place.
203  // One should also be concerned that if we had such API in place, rank 0 in subcommunicator may no longer be rank 0 in
204  // MPI_COMM_WORLD, which may lead to issues for logging.
205  if (remapPartitions) {
206  SubFactoryMonitor m1(*this, "DeterminePartitionPlacement", currentLevel);
207 
208  DeterminePartitionPlacement(*A, *decomposition, numPartitions);
209  }
210 
211  // ======================================================================================================
212  // Construct importer
213  // ======================================================================================================
214  // At this point, the following is true:
215  // * Each processors owns 0 or 1 partitions
216  // * If a processor owns a partition, that partition number is equal to the processor rank
217  // * The decomposition vector contains the partitions ids that the corresponding GID belongs to
218 
219  ArrayRCP<const GO> decompEntries;
220  if (decomposition->getLocalLength() > 0)
221  decompEntries = decomposition->getData(0);
222 
223 #ifdef HAVE_MUELU_DEBUG
224  // Test range of partition ids
225  int incorrectRank = -1;
226  for (int i = 0; i < decompEntries.size(); i++)
227  if (decompEntries[i] >= numProcs || decompEntries[i] < 0) {
228  incorrectRank = myRank;
229  break;
230  }
231 
232  int incorrectGlobalRank = -1;
233  MueLu_maxAll(comm, incorrectRank, incorrectGlobalRank);
234  TEUCHOS_TEST_FOR_EXCEPTION(incorrectGlobalRank >- 1, Exceptions::RuntimeError, "pid " + Teuchos::toString(incorrectGlobalRank) + " encountered a partition number is that out-of-range");
235 #endif
236 
237  Array<GO> myGIDs;
238  myGIDs.reserve(decomposition->getLocalLength());
239 
240  // Step 0: Construct mapping
241  // part number -> GIDs I own which belong to this part
242  // NOTE: my own part GIDs are not part of the map
243  typedef std::map<GO, Array<GO> > map_type;
244  map_type sendMap;
245  for (LO i = 0; i < decompEntries.size(); i++) {
246  GO id = decompEntries[i];
247  GO GID = rowMap->getGlobalElement(i);
248 
249  if (id == myRank)
250  myGIDs .push_back(GID);
251  else
252  sendMap[id].push_back(GID);
253  }
254  decompEntries = Teuchos::null;
255 
256  if (IsPrint(Statistics2)) {
257  GO numLocalKept = myGIDs.size(), numGlobalKept, numGlobalRows = A->getGlobalNumRows();
258  MueLu_sumAll(comm,numLocalKept, numGlobalKept);
259  GetOStream(Statistics2) << "Unmoved rows: " << numGlobalKept << " / " << numGlobalRows << " (" << 100*Teuchos::as<double>(numGlobalKept)/numGlobalRows << "%)" << std::endl;
260  }
261 
262  int numSend = sendMap.size(), numRecv;
263 
264  // Arrayify map keys
265  Array<GO> myParts(numSend), myPart(1);
266  int cnt = 0;
267  myPart[0] = myRank;
268  for (typename map_type::const_iterator it = sendMap.begin(); it != sendMap.end(); it++)
269  myParts[cnt++] = it->first;
270 
271  // Step 1: Find out how many processors send me data
272  // partsIndexBase starts from zero, as the processors ids start from zero
273  GO partsIndexBase = 0;
274  RCP<Map> partsIHave = MapFactory ::Build(lib, Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), myParts(), partsIndexBase, comm);
275  RCP<Map> partsIOwn = MapFactory ::Build(lib, numProcs, myPart(), partsIndexBase, comm);
276  RCP<Export> partsExport = ExportFactory::Build(partsIHave, partsIOwn);
277 
278  RCP<GOVector> partsISend = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(partsIHave);
279  RCP<GOVector> numPartsIRecv = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(partsIOwn);
280  if (numSend) {
281  ArrayRCP<GO> partsISendData = partsISend->getDataNonConst(0);
282  for (int i = 0; i < numSend; i++)
283  partsISendData[i] = 1;
284  }
285  (numPartsIRecv->getDataNonConst(0))[0] = 0;
286 
287  numPartsIRecv->doExport(*partsISend, *partsExport, Xpetra::ADD);
288  numRecv = (numPartsIRecv->getData(0))[0];
289 
290  // Step 2: Get my GIDs from everybody else
291  MPI_Datatype MpiType = MpiTypeTraits<GO>::getType();
292  int msgTag = 12345; // TODO: use Comm::dup for all internal messaging
293 
294  // Post sends
295  Array<MPI_Request> sendReqs(numSend);
296  cnt = 0;
297  for (typename map_type::iterator it = sendMap.begin(); it != sendMap.end(); it++)
298  MPI_Isend(static_cast<void*>(it->second.getRawPtr()), it->second.size(), MpiType, Teuchos::as<GO>(it->first), msgTag, *rawMpiComm, &sendReqs[cnt++]);
299 
300  map_type recvMap;
301  size_t totalGIDs = myGIDs.size();
302  for (int i = 0; i < numRecv; i++) {
303  MPI_Status status;
304  MPI_Probe(MPI_ANY_SOURCE, msgTag, *rawMpiComm, &status);
305 
306  // Get rank and number of elements from status
307  int fromRank = status.MPI_SOURCE, count;
308  MPI_Get_count(&status, MpiType, &count);
309 
310  recvMap[fromRank].resize(count);
311  MPI_Recv(static_cast<void*>(recvMap[fromRank].getRawPtr()), count, MpiType, fromRank, msgTag, *rawMpiComm, &status);
312 
313  totalGIDs += count;
314  }
315 
316  // Do waits on send requests
317  if (numSend) {
318  Array<MPI_Status> sendStatuses(numSend);
319  MPI_Waitall(numSend, sendReqs.getRawPtr(), sendStatuses.getRawPtr());
320  }
321 
322  // Merge GIDs
323  myGIDs.reserve(totalGIDs);
324  for (typename map_type::const_iterator it = recvMap.begin(); it != recvMap.end(); it++) {
325  int offset = myGIDs.size(), len = it->second.size();
326  if (len) {
327  myGIDs.resize(offset + len);
328  memcpy(myGIDs.getRawPtr() + offset, it->second.getRawPtr(), len*sizeof(GO));
329  }
330  }
331  // NOTE 2: The general sorting algorithm could be sped up by using the knowledge that original myGIDs and all received chunks
332  // (i.e. it->second) are sorted. Therefore, a merge sort would work well in this situation.
333  std::sort(myGIDs.begin(), myGIDs.end());
334 
335  // Step 3: Construct importer
336  RCP<Map> newRowMap = MapFactory ::Build(lib, rowMap->getGlobalNumElements(), myGIDs(), indexBase, origComm);
337  RCP<const Import> rowMapImporter;
338  {
339  SubFactoryMonitor m1(*this, "Import construction", currentLevel);
340  rowMapImporter = ImportFactory::Build(rowMap, newRowMap);
341  }
342 
343  Set(currentLevel, "Importer", rowMapImporter);
344 
345  // ======================================================================================================
346  // Print some data
347  // ======================================================================================================
348  if (pL.get<bool>("repartition: print partition distribution") && IsPrint(Statistics2)) {
349  // Print the grid of processors
350  GetOStream(Statistics2) << "Partition distribution over cores (ownership is indicated by '+')" << std::endl;
351 
352  char amActive = (myGIDs.size() ? 1 : 0);
353  std::vector<char> areActive(numProcs, 0);
354  MPI_Gather(&amActive, 1, MPI_CHAR, &areActive[0], 1, MPI_CHAR, 0, *rawMpiComm);
355 
356  int rowWidth = std::min(Teuchos::as<int>(ceil(sqrt(numProcs))), 100);
357  for (int proc = 0; proc < numProcs; proc += rowWidth) {
358  for (int j = 0; j < rowWidth; j++)
359  if (proc + j < numProcs)
360  GetOStream(Statistics2) << (areActive[proc + j] ? "+" : ".");
361  else
362  GetOStream(Statistics2) << " ";
363 
364  GetOStream(Statistics2) << " " << proc << ":" << std::min(proc + rowWidth, numProcs) - 1 << std::endl;
365  }
366  }
367 
368  } // Build
369 
370  //----------------------------------------------------------------------
371  template<typename T, typename W>
372  struct Triplet {
373  T i, j;
374  W v;
375  };
376  template<typename T, typename W>
377  static bool compareTriplets(const Triplet<T,W>& a, const Triplet<T,W>& b) {
378  return (a.v > b.v); // descending order
379  }
380 
381  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
383  DeterminePartitionPlacement(const Matrix& A, GOVector& decomposition, GO numPartitions) const {
384  RCP<const Map> rowMap = A.getRowMap();
385 
386  RCP<const Teuchos::Comm<int> > comm = rowMap->getComm()->duplicate();
387  int numProcs = comm->getSize();
388 
389  RCP<const Teuchos::MpiComm<int> > tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
390  TEUCHOS_TEST_FOR_EXCEPTION(tmpic == Teuchos::null, Exceptions::RuntimeError, "Cannot cast base Teuchos::Comm to Teuchos::MpiComm object.");
391  RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawMpiComm = tmpic->getRawMpiComm();
392 
393  const Teuchos::ParameterList& pL = GetParameterList();
394 
395  // maxLocal is a constant which determins the number of largest edges which are being exchanged
396  // The idea is that we do not want to construct the full bipartite graph, but simply a subset of
397  // it, which requires less communication. By selecting largest local edges we hope to achieve
398  // similar results but at a lower cost.
399  const int maxLocal = pL.get<int>("repartition: remap num values");
400  const int dataSize = 2*maxLocal;
401 
402  ArrayRCP<GO> decompEntries;
403  if (decomposition.getLocalLength() > 0)
404  decompEntries = decomposition.getDataNonConst(0);
405 
406  // Step 1: Sort local edges by weight
407  // Each edge of a bipartite graph corresponds to a triplet (i, j, v) where
408  // i: processor id that has some piece of part with part_id = j
409  // j: part id
410  // v: weight of the edge
411  // We set edge weights to be the total number of nonzeros in rows on this processor which
412  // correspond to this part_id. The idea is that when we redistribute matrix, this weight
413  // is a good approximation of the amount of data to move.
414  // We use two maps, original which maps a partition id of an edge to the corresponding weight,
415  // and a reverse one, which is necessary to sort by edges.
416  std::map<GO,GO> lEdges;
417  for (LO i = 0; i < decompEntries.size(); i++)
418  lEdges[decompEntries[i]] += A.getNumEntriesInLocalRow(i);
419 
420  // Reverse map, so that edges are sorted by weight.
421  // This results in multimap, as we may have edges with the same weight
422  std::multimap<GO,GO> revlEdges;
423  for (typename std::map<GO,GO>::const_iterator it = lEdges.begin(); it != lEdges.end(); it++)
424  revlEdges.insert(std::make_pair(it->second, it->first));
425 
426  // Both lData and gData are arrays of data which we communicate. The data is stored
427  // in pairs, so that data[2*i+0] is the part index, and data[2*i+1] is the corresponding edge weight.
428  // We do not store processor id in data, as we can compute that by looking on the offset in the gData.
429  Array<GO> lData(dataSize, -1), gData(numProcs * dataSize);
430  int numEdges = 0;
431  for (typename std::multimap<GO,GO>::reverse_iterator rit = revlEdges.rbegin(); rit != revlEdges.rend() && numEdges < maxLocal; rit++) {
432  lData[2*numEdges+0] = rit->second; // part id
433  lData[2*numEdges+1] = rit->first; // edge weight
434  numEdges++;
435  }
436 
437  // Step 2: Gather most edges
438  // Each processors contributes maxLocal edges by providing maxLocal pairs <part id, weight>, which is of size dataSize
439  MPI_Datatype MpiType = MpiTypeTraits<GO>::getType();
440  MPI_Allgather(static_cast<void*>(lData.getRawPtr()), dataSize, MpiType, static_cast<void*>(gData.getRawPtr()), dataSize, MpiType, *rawMpiComm);
441 
442  // Step 3: Construct mapping
443 
444  // Construct the set of triplets
445  std::vector<Triplet<int,int> > gEdges(numProcs * maxLocal);
446  size_t k = 0;
447  for (LO i = 0; i < gData.size(); i += 2) {
448  GO part = gData[i+0];
449  GO weight = gData[i+1];
450  if (part != -1) { // skip nonexistent edges
451  gEdges[k].i = i/dataSize; // determine the processor by its offset (since every processor sends the same amount)
452  gEdges[k].j = part;
453  gEdges[k].v = weight;
454  k++;
455  }
456  }
457  gEdges.resize(k);
458 
459  // Sort edges by weight
460  // NOTE: compareTriplets is actually a reverse sort, so the edges weight is in decreasing order
461  std::sort(gEdges.begin(), gEdges.end(), compareTriplets<int,int>);
462 
463  // Do matching
464  std::map<int,int> match;
465  std::vector<char> matchedRanks(numProcs, 0), matchedParts(numProcs, 0);
466  int numMatched = 0;
467  for (typename std::vector<Triplet<int,int> >::const_iterator it = gEdges.begin(); it != gEdges.end(); it++) {
468  GO rank = it->i;
469  GO part = it->j;
470  if (matchedRanks[rank] == 0 && matchedParts[part] == 0) {
471  matchedRanks[rank] = 1;
472  matchedParts[part] = 1;
473  match[part] = rank;
474  numMatched++;
475  }
476  }
477  GetOStream(Statistics1) << "Number of unassigned paritions before cleanup stage: " << (numPartitions - numMatched) << " / " << numPartitions << std::endl;
478 
479  // Step 4: Assign unassigned partitions
480  // We do that through random matching for remaining partitions. Not all part numbers are valid, but valid parts are a subset of [0, numProcs).
481  // The reason it is done this way is that we don't need any extra communication, as we don't need to know which parts are valid.
482  for (int part = 0, matcher = 0; part < numProcs; part++)
483  if (match.count(part) == 0) {
484  // Find first non-matched rank
485  while (matchedRanks[matcher])
486  matcher++;
487 
488  match[part] = matcher++;
489  }
490 
491  // Step 5: Permute entries in the decomposition vector
492  for (LO i = 0; i < decompEntries.size(); i++)
493  decompEntries[i] = match[decompEntries[i]];
494  }
495 
496 } // namespace MueLu
497 
498 #endif //ifdef HAVE_MPI
499 
500 #endif // MUELU_REPARTITIONFACTORY_DEF_HPP
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Important warning messages (one line)
#define MueLu_sumAll(rcpComm, in, out)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
#define MueLu_maxAll(rcpComm, in, out)
GlobalOrdinal GO
Timer to be used in factories. Similar to Monitor but with additional timers.
static MPI_Datatype getType()
Print more statistics.
LocalOrdinal LO
Namespace for MueLu classes and methods.
void DeterminePartitionPlacement(const Matrix &A, GOVector &decomposition, GO numPartitions) const
Determine which process should own each partition.
#define SET_VALID_ENTRY(name)
Print even more statistics.
void Build(Level &currentLevel) const
Build an object with this factory.
static bool compareTriplets(const Triplet< T, W > &a, const Triplet< T, W > &b)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
void DeclareInput(Level &currentLevel) const
Determines the data that RepartitionFactory needs, and the factories that generate that data...
static RCP< Vector > Build(const Teuchos::RCP< const Map > &map, bool zeroOut=true)
Exception throws to report errors in the internal logical of the program.