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