MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RepartitionFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_REPARTITIONFACTORY_DEF_HPP
11#define MUELU_REPARTITIONFACTORY_DEF_HPP
12
13#include <algorithm>
14#include <iostream>
15#include <sstream>
16
17#include "MueLu_RepartitionFactory_decl.hpp" // TMP JG NOTE: before other includes, otherwise I cannot test the fwd declaration in _def
18
19#ifdef HAVE_MPI
20#include <Teuchos_DefaultMpiComm.hpp>
21#include <Teuchos_CommHelpers.hpp>
22#include <Teuchos_Details_MpiTypeTraits.hpp>
23
24#include <Xpetra_Map.hpp>
25#include <Xpetra_MapFactory.hpp>
26#include <Xpetra_MultiVectorFactory.hpp>
27#include <Xpetra_VectorFactory.hpp>
28#include <Xpetra_Import.hpp>
29#include <Xpetra_ImportFactory.hpp>
30#include <Xpetra_Export.hpp>
31#include <Xpetra_ExportFactory.hpp>
32#include <Xpetra_Matrix.hpp>
33#include <Xpetra_MatrixFactory.hpp>
34
35#include "MueLu_Utilities.hpp"
36
37#include "MueLu_CloneRepartitionInterface.hpp"
38
39#include "MueLu_Level.hpp"
40#include "MueLu_MasterList.hpp"
41#include "MueLu_Monitor.hpp"
42#include "MueLu_PerfUtils.hpp"
43
44#include "MueLu_RepartitionUtilities.hpp"
45
46namespace MueLu {
47
48template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
50
51template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
53
54template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56 RCP<ParameterList> validParamList = rcp(new ParameterList());
57
58#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
59 SET_VALID_ENTRY("repartition: print partition distribution");
60 SET_VALID_ENTRY("repartition: remap parts");
61 SET_VALID_ENTRY("repartition: remap num values");
62 SET_VALID_ENTRY("repartition: remap accept partition");
63 SET_VALID_ENTRY("repartition: node repartition level");
64 SET_VALID_ENTRY("repartition: save importer");
65#undef SET_VALID_ENTRY
66
67 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Factory of the matrix A");
68 validParamList->set<RCP<const FactoryBase> >("number of partitions", Teuchos::null, "Instance of RepartitionHeuristicFactory.");
69 validParamList->set<RCP<const FactoryBase> >("Partition", Teuchos::null, "Factory of the partition");
70
71 return validParamList;
72}
73
74template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
76 Input(currentLevel, "A");
77 Input(currentLevel, "number of partitions");
78 Input(currentLevel, "Partition");
79}
80
81template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 FactoryMonitor m(*this, "Build", currentLevel);
84
86
87 const Teuchos::ParameterList& pL = GetParameterList();
88 // Access parameters here to make sure that we set the parameter entry flag to "used" even in case of short-circuit evaluation.
89 // TODO (JG): I don't really know if we want to do this.
90 bool remapPartitions = pL.get<bool>("repartition: remap parts");
91
92 // TODO: We only need a CrsGraph. This class does not have to be templated on Scalar types.
93 RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel, "A");
94 if (A == Teuchos::null) {
95 Set<RCP<const Import> >(currentLevel, "Importer", Teuchos::null);
96 return;
97 }
98 RCP<const Map> rowMap = A->getRowMap();
99 GO indexBase = rowMap->getIndexBase();
100 Xpetra::UnderlyingLib lib = rowMap->lib();
101
102 RCP<const Teuchos::Comm<int> > origComm = rowMap->getComm();
103 RCP<const Teuchos::Comm<int> > comm = origComm;
104
105 int numProcs = comm->getSize();
106
107 RCP<const Teuchos::MpiComm<int> > tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
108 TEUCHOS_TEST_FOR_EXCEPTION(tmpic == Teuchos::null, Exceptions::RuntimeError, "Cannot cast base Teuchos::Comm to Teuchos::MpiComm object.");
109 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawMpiComm = tmpic->getRawMpiComm();
110
112 int numPartitions = Get<int>(currentLevel, "number of partitions");
113
114 // ======================================================================================================
115 // Construct decomposition vector
116 // ======================================================================================================
117 RCP<GOVector> decomposition = Get<RCP<GOVector> >(currentLevel, "Partition");
118
119 // check which factory provides "Partition"
120 if (remapPartitions == true && Teuchos::rcp_dynamic_cast<const CloneRepartitionInterface>(GetFactory("Partition")) != Teuchos::null) {
121 // if "Partition" is provided by a CloneRepartitionInterface class we have to switch of remapPartitions
122 // as we can assume the processor Ids in Partition to be the expected ones. If we would do remapping we
123 // would get different processors for the different blocks which screws up matrix-matrix multiplication.
124 remapPartitions = false;
125 }
126
127 // check special cases
128 if (numPartitions == 1) {
129 // Trivial case: decomposition is the trivial one, all zeros. We skip the call to Zoltan_Interface
130 // (this is mostly done to avoid extra output messages, as even if we didn't skip there is a shortcut
131 // in Zoltan[12]Interface).
132 // TODO: We can probably skip more work in this case (like building all extra data structures)
133 GetOStream(Runtime0) << "Only one partition: Skip call to the repartitioner." << std::endl;
134 } else if (numPartitions == -1) {
135 // No repartitioning necessary: decomposition should be Teuchos::null
136 GetOStream(Runtime0) << "No repartitioning necessary: partitions were left unchanged by the repartitioner" << std::endl;
137 Set<RCP<const Import> >(currentLevel, "Importer", Teuchos::null);
138 return;
139 }
140
141 // If we're doing node away, we need to be sure to get the mapping to the NodeComm's rank 0.
142 const int nodeRepartLevel = pL.get<int>("repartition: node repartition level");
143 if (currentLevel.GetLevelID() == nodeRepartLevel) {
144 // NodePartitionInterface returns the *ranks* of the guy who gets the info, not the *partition number*
145 // In a sense, we've already done remap here.
146
147 // FIXME: We need a low-comm import construction
148 remapPartitions = false;
149 }
150
151 // ======================================================================================================
152 // Remap if necessary
153 // ======================================================================================================
154 // From a user perspective, we want user to not care about remapping, thinking of it as only a performance feature.
155 // There are two problems, however.
156 // (1) Next level aggregation depends on the order of GIDs in the vector, if one uses "natural" or "random" orderings.
157 // This also means that remapping affects next level aggregation, despite the fact that the _set_ of GIDs for
158 // each partition is the same.
159 // (2) Even with the fixed order of GIDs, the remapping may influence the aggregation for the next-next level.
160 // Let us consider the following example. Lets assume that when we don't do remapping, processor 0 would have
161 // GIDs {0,1,2}, and processor 1 GIDs {3,4,5}, and if we do remapping processor 0 would contain {3,4,5} and
162 // processor 1 {0,1,2}. Now, when we run repartitioning algorithm on the next level (say Zoltan1 RCB), it may
163 // be dependent on whether whether it is [{0,1,2}, {3,4,5}] or [{3,4,5}, {0,1,2}]. Specifically, the tie-breaking
164 // algorithm can resolve these differently. For instance, running
165 // mpirun -np 5 ./MueLu_ScalingTestParamList.exe --xml=easy_sa.xml --nx=12 --ny=12 --nz=12
166 // with
167 // <ParameterList name="MueLu">
168 // <Parameter name="coarse: max size" type="int" value="1"/>
169 // <Parameter name="repartition: enable" type="bool" value="true"/>
170 // <Parameter name="repartition: min rows per proc" type="int" value="2"/>
171 // <ParameterList name="level 1">
172 // <Parameter name="repartition: remap parts" type="bool" value="false/true"/>
173 // </ParameterList>
174 // </ParameterList>
175 // produces different repartitioning for level 2.
176 // This different repartitioning may then escalate into different aggregation for the next level.
177 //
178 // We fix (1) by fixing the order of GIDs in a vector by sorting the resulting vector.
179 // Fixing (2) is more complicated.
180 // FIXME: Fixing (2) in Zoltan may not be enough, as we may use some arbitration in MueLu,
181 // for instance with CoupledAggregation. What we really need to do is to use the same order of processors containing
182 // the same order of GIDs. To achieve that, the newly created subcommunicator must be conforming with the order. For
183 // instance, if we have [{0,1,2}, {3,4,5}], we create a subcommunicator where processor 0 gets rank 0, and processor 1
184 // 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.
185 // This rank permutation requires help from Epetra/Tpetra, both of which have no such API in place.
186 // One should also be concerned that if we had such API in place, rank 0 in subcommunicator may no longer be rank 0 in
187 // MPI_COMM_WORLD, which may lead to issues for logging.
188 if (remapPartitions) {
189 SubFactoryMonitor m1(*this, "DeterminePartitionPlacement", currentLevel);
190
191 bool acceptPartition = pL.get<bool>("repartition: remap accept partition");
192 bool allSubdomainsAcceptPartitions;
193 int localNumAcceptPartition = acceptPartition;
194 int globalNumAcceptPartition;
195 MueLu_sumAll(comm, localNumAcceptPartition, globalNumAcceptPartition);
196 GetOStream(Statistics2) << "Number of ranks that accept partitions: " << globalNumAcceptPartition << std::endl;
197 if (globalNumAcceptPartition < numPartitions) {
198 GetOStream(Warnings0) << "Not enough ranks are willing to accept a partition, allowing partitions on all ranks." << std::endl;
199 acceptPartition = true;
200 allSubdomainsAcceptPartitions = true;
201 } else if (numPartitions > numProcs) {
202 // We are trying to repartition to a larger communicator.
203 allSubdomainsAcceptPartitions = true;
204 } else {
205 allSubdomainsAcceptPartitions = false;
206 }
207
208 DeterminePartitionPlacement(*A, *decomposition, numPartitions, acceptPartition, allSubdomainsAcceptPartitions);
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#ifdef HAVE_MUELU_DEBUG
220 int myRank = comm->getRank();
221 ArrayRCP<const GO> decompEntries;
222 if (decomposition->getLocalLength() > 0)
223 decompEntries = decomposition->getData(0);
224
225 // Test range of partition ids
226 int incorrectRank = -1;
227 for (int i = 0; i < decompEntries.size(); i++)
228 if (decompEntries[i] >= numProcs || decompEntries[i] < 0) {
229 incorrectRank = myRank;
230 break;
231 }
232
233 int incorrectGlobalRank = -1;
234 MueLu_maxAll(comm, incorrectRank, incorrectGlobalRank);
235 TEUCHOS_TEST_FOR_EXCEPTION(incorrectGlobalRank > -1, Exceptions::RuntimeError, "pid " + Teuchos::toString(incorrectGlobalRank) + " encountered a partition number is that out-of-range");
236#endif
237
238 // Step 1: Construct GIDs
239 Array<GO> myGIDs;
240 {
241 SubFactoryMonitor m1(*this, "Compute GIDs", currentLevel);
242 GO numLocalKept;
243 std::tie(myGIDs, numLocalKept) = Util::ConstructGIDs(decomposition);
244
245 if (IsPrint(Statistics2)) {
246 GO numGlobalKept, numGlobalRows = A->getGlobalNumRows();
247 MueLu_sumAll(comm, numLocalKept, numGlobalKept);
248 GetOStream(Statistics2) << "Unmoved rows: " << numGlobalKept << " / " << numGlobalRows << " (" << 100 * Teuchos::as<double>(numGlobalKept) / numGlobalRows << "%)" << std::endl;
249 }
250 }
251
252 // Step 2: Construct importer
253 RCP<Map> newRowMap;
254 {
255 SubFactoryMonitor m1(*this, "Map construction", currentLevel);
256 newRowMap = MapFactory ::Build(lib, rowMap->getGlobalNumElements(), myGIDs(), indexBase, origComm);
257 }
258
259 RCP<const Import> rowMapImporter;
260
261 RCP<const BlockedMap> blockedRowMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rowMap);
262
263 {
264 SubFactoryMonitor m1(*this, "Import construction", currentLevel);
265 // Generate a nonblocked rowmap if we need one
266 if (blockedRowMap.is_null())
267 rowMapImporter = ImportFactory::Build(rowMap, newRowMap);
268 else {
269 rowMapImporter = ImportFactory::Build(blockedRowMap->getMap(), newRowMap);
270 }
271 }
272
273 // If we're running BlockedCrs we should chop up the newRowMap into a newBlockedRowMap here (and do likewise for importers)
274 if (!blockedRowMap.is_null()) {
275 SubFactoryMonitor m1(*this, "Blocking newRowMap and Importer", currentLevel);
276 RCP<const BlockedMap> blockedTargetMap = MueLu::UtilitiesBase<Scalar, LocalOrdinal, GlobalOrdinal, Node>::GeneratedBlockedTargetMap(*blockedRowMap, *rowMapImporter);
277
278 // NOTE: This code qualifies as "correct but not particularly performant" If this needs to be sped up, we can probably read data from the existing importer to
279 // build sub-importers rather than generating new ones ex nihilo
280 size_t numBlocks = blockedRowMap->getNumMaps();
281 std::vector<RCP<const Import> > subImports(numBlocks);
282
283 for (size_t i = 0; i < numBlocks; i++) {
284 RCP<const Map> source = blockedRowMap->getMap(i);
285 RCP<const Map> target = blockedTargetMap->getMap(i);
286 subImports[i] = ImportFactory::Build(source, target);
287 }
288 Set(currentLevel, "SubImporters", subImports);
289 }
290
291 Set(currentLevel, "Importer", rowMapImporter);
292
293 // Importer saving
294 bool save_importer = pL.get<bool>("repartition: save importer");
295 if (save_importer) {
296 currentLevel.Set("Importer", rowMapImporter, NoFactory::get());
297 currentLevel.AddKeepFlag("Importer", NoFactory::get(), MueLu::Final);
298 currentLevel.RemoveKeepFlag("Importer", NoFactory::get(), MueLu::UserData); // FIXME: This is a hack
299 }
300 // ======================================================================================================
301 // Print some data
302 // ======================================================================================================
303 if (!rowMapImporter.is_null() && IsPrint(Statistics2)) {
304 // int oldRank = SetProcRankVerbose(rebalancedAc->getRowMap()->getComm()->getRank());
305 GetOStream(Statistics2) << PerfUtils::PrintImporterInfo(rowMapImporter, "Importer for rebalancing");
306 // SetProcRankVerbose(oldRank);
307 }
308 if (pL.get<bool>("repartition: print partition distribution") && IsPrint(Statistics2)) {
309 // Print the grid of processors
310 GetOStream(Statistics2) << "Partition distribution over cores (ownership is indicated by '+')" << std::endl;
311
312 char amActive = (myGIDs.size() ? 1 : 0);
313 std::vector<char> areActive(numProcs, 0);
314 MPI_Gather(&amActive, 1, MPI_CHAR, &areActive[0], 1, MPI_CHAR, 0, *rawMpiComm);
315
316 int rowWidth = std::min(Teuchos::as<int>(ceil(sqrt(numProcs))), 100);
317 for (int proc = 0; proc < numProcs; proc += rowWidth) {
318 for (int j = 0; j < rowWidth; j++)
319 if (proc + j < numProcs)
320 GetOStream(Statistics2) << (areActive[proc + j] ? "+" : ".");
321 else
322 GetOStream(Statistics2) << " ";
323
324 GetOStream(Statistics2) << " " << proc << ":" << std::min(proc + rowWidth, numProcs) - 1 << std::endl;
325 }
326 }
327
328} // Build
329
330//----------------------------------------------------------------------
331template <typename T, typename W>
332struct Triplet {
333 T i, j;
334 W v;
335};
336template <typename T, typename W>
337static bool compareTriplets(const Triplet<T, W>& a, const Triplet<T, W>& b) {
338 return (a.v > b.v); // descending order
339}
340
341template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
343 DeterminePartitionPlacement(const Matrix& A, GOVector& decomposition, GO numPartitions, bool willAcceptPartition, bool allSubdomainsAcceptPartitions) const {
344 RCP<const Map> rowMap = A.getRowMap();
345
346 RCP<const Teuchos::Comm<int> > comm = rowMap->getComm()->duplicate();
347 int numProcs = comm->getSize();
348
349 RCP<const Teuchos::MpiComm<int> > tmpic = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
350 TEUCHOS_TEST_FOR_EXCEPTION(tmpic == Teuchos::null, Exceptions::RuntimeError, "Cannot cast base Teuchos::Comm to Teuchos::MpiComm object.");
351 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > rawMpiComm = tmpic->getRawMpiComm();
352
353 const Teuchos::ParameterList& pL = GetParameterList();
354
355 // maxLocal is a constant which determins the number of largest edges which are being exchanged
356 // The idea is that we do not want to construct the full bipartite graph, but simply a subset of
357 // it, which requires less communication. By selecting largest local edges we hope to achieve
358 // similar results but at a lower cost.
359 const int maxLocal = pL.get<int>("repartition: remap num values");
360 const int dataSize = 2 * maxLocal;
361
362 ArrayRCP<GO> decompEntries;
363 if (decomposition.getLocalLength() > 0)
364 decompEntries = decomposition.getDataNonConst(0);
365
366 // Step 1: Sort local edges by weight
367 // Each edge of a bipartite graph corresponds to a triplet (i, j, v) where
368 // i: processor id that has some piece of part with part_id = j
369 // j: part id
370 // v: weight of the edge
371 // We set edge weights to be the total number of nonzeros in rows on this processor which
372 // correspond to this part_id. The idea is that when we redistribute matrix, this weight
373 // is a good approximation of the amount of data to move.
374 // We use two maps, original which maps a partition id of an edge to the corresponding weight,
375 // and a reverse one, which is necessary to sort by edges.
376 std::map<GO, GO> lEdges;
377 if (willAcceptPartition)
378 for (LO i = 0; i < decompEntries.size(); i++)
379 lEdges[decompEntries[i]] += A.getNumEntriesInLocalRow(i);
380
381 // Reverse map, so that edges are sorted by weight.
382 // This results in multimap, as we may have edges with the same weight
383 std::multimap<GO, GO> revlEdges;
384 for (typename std::map<GO, GO>::const_iterator it = lEdges.begin(); it != lEdges.end(); it++)
385 revlEdges.insert(std::make_pair(it->second, it->first));
386
387 // Both lData and gData are arrays of data which we communicate. The data is stored
388 // in pairs, so that data[2*i+0] is the part index, and data[2*i+1] is the corresponding edge weight.
389 // We do not store processor id in data, as we can compute that by looking on the offset in the gData.
390 Array<GO> lData(dataSize, -1), gData(numProcs * dataSize);
391 int numEdges = 0;
392 for (typename std::multimap<GO, GO>::reverse_iterator rit = revlEdges.rbegin(); rit != revlEdges.rend() && numEdges < maxLocal; rit++) {
393 lData[2 * numEdges + 0] = rit->second; // part id
394 lData[2 * numEdges + 1] = rit->first; // edge weight
395 numEdges++;
396 }
397
398 // Step 2: Gather most edges
399 // Each processors contributes maxLocal edges by providing maxLocal pairs <part id, weight>, which is of size dataSize
400 MPI_Datatype MpiType = Teuchos::Details::MpiTypeTraits<GO>::getType();
401 MPI_Allgather(static_cast<void*>(lData.getRawPtr()), dataSize, MpiType, static_cast<void*>(gData.getRawPtr()), dataSize, MpiType, *rawMpiComm);
402
403 // Step 3: Construct mapping
404
405 // Construct the set of triplets
406 Teuchos::Array<Triplet<int, int> > gEdges(numProcs * maxLocal);
407 Teuchos::Array<bool> procWillAcceptPartition(numProcs, allSubdomainsAcceptPartitions);
408 size_t k = 0;
409 for (LO i = 0; i < gData.size(); i += 2) {
410 int procNo = i / dataSize; // determine the processor by its offset (since every processor sends the same amount)
411 GO part = gData[i + 0];
412 GO weight = gData[i + 1];
413 if (part != -1) { // skip nonexistent edges
414 gEdges[k].i = procNo;
415 gEdges[k].j = part;
416 gEdges[k].v = weight;
417 procWillAcceptPartition[procNo] = true;
418 k++;
419 }
420 }
421 gEdges.resize(k);
422
423 // Sort edges by weight
424 // NOTE: compareTriplets is actually a reverse sort, so the edges weight is in decreasing order
425 std::sort(gEdges.begin(), gEdges.end(), compareTriplets<int, int>);
426
427 // Do matching
428 std::map<int, int> match;
429 Teuchos::Array<char> matchedRanks(numProcs, 0), matchedParts(numPartitions, 0);
430 int numMatched = 0;
431 for (typename Teuchos::Array<Triplet<int, int> >::const_iterator it = gEdges.begin(); it != gEdges.end(); it++) {
432 GO rank = it->i;
433 GO part = it->j;
434 if (matchedRanks[rank] == 0 && matchedParts[part] == 0) {
435 matchedRanks[rank] = 1;
436 matchedParts[part] = 1;
437 match[part] = rank;
438 numMatched++;
439 }
440 }
441 GetOStream(Statistics1) << "Number of unassigned partitions before cleanup stage: " << (numPartitions - numMatched) << " / " << numPartitions << std::endl;
442
443 // Step 4: Assign unassigned partitions if necessary.
444 // We do that through desperate matching for remaining partitions:
445 // We select the lowest rank that can still take a partition.
446 // The reason it is done this way is that we don't need any extra communication, as we don't
447 // need to know which parts are valid.
448 if (numPartitions - numMatched > 0) {
449 Teuchos::Array<char> partitionCounts(numPartitions, 0);
450 for (typename std::map<int, int>::const_iterator it = match.begin(); it != match.end(); it++)
451 partitionCounts[it->first] += 1;
452 for (int part = 0, matcher = 0; part < numPartitions; part++) {
453 if (partitionCounts[part] == 0) {
454 // Find first non-matched rank that accepts partitions
455 while (matchedRanks[matcher] || !procWillAcceptPartition[matcher])
456 matcher++;
457
458 match[part] = matcher++;
459 numMatched++;
460 }
461 }
462 }
463
464 TEUCHOS_TEST_FOR_EXCEPTION(numMatched != numPartitions, Exceptions::RuntimeError, "MueLu::RepartitionFactory::DeterminePartitionPlacement: Only " << numMatched << " partitions out of " << numPartitions << " got assigned to ranks.");
465
466 // Step 5: Permute entries in the decomposition vector
467 for (LO i = 0; i < decompEntries.size(); i++)
468 decompEntries[i] = match[decompEntries[i]];
469}
470
471} // namespace MueLu
472
473#endif // ifdef HAVE_MPI
474
475#endif // MUELU_REPARTITIONFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
#define MueLu_maxAll(rcpComm, in, out)
#define MueLu_sumAll(rcpComm, in, out)
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
void RemoveKeepFlag(const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All)
int GetLevelID() const
Return level number.
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
static const NoFactory * get()
static std::string PrintImporterInfo(RCP< const Import > importer, const std::string &msgTag)
virtual ~RepartitionFactory()
Destructor.
void Build(Level &currentLevel) const
Build an object with this factory.
void DeterminePartitionPlacement(const Matrix &A, GOVector &decomposition, GO numPartitions, bool willAcceptPartition=true, bool allSubdomainsAcceptPartitions=true) const
Determine which process should own each partition.
void DeclareInput(Level &currentLevel) const
Determines the data that RepartitionFactory needs, and the factories that generate that data.
RepartitionFactory()
Constructor.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static RCP< const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > > GeneratedBlockedTargetMap(const Xpetra::BlockedMap< LocalOrdinal, GlobalOrdinal, Node > &sourceBlockedMap, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &Importer)
Namespace for MueLu classes and methods.
@ Final
Keep data only for this run. Used to keep data useful for Hierarchy::Iterate(). Data will be deleted ...
@ UserData
User data are always kept. This flag is set automatically when Level::Set("data", data) is used....
static bool compareTriplets(const Triplet< T, W > &a, const Triplet< T, W > &b)
@ Warnings0
Important warning messages (one line)
@ Statistics2
Print even more statistics.
@ Statistics1
Print more statistics.
@ Runtime0
One-liner description of what is happening.