MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RebalanceTransferFactory_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_REBALANCETRANSFERFACTORY_DEF_HPP
11#define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
12
13#include <sstream>
14#include <Teuchos_Tuple.hpp>
15
16#include "Xpetra_MultiVector.hpp"
17#include "Xpetra_MultiVectorFactory.hpp"
18#include "Xpetra_Vector.hpp"
19#include "Xpetra_VectorFactory.hpp"
20#include <Xpetra_Matrix.hpp>
21#include <Xpetra_MapFactory.hpp>
22#include <Xpetra_MatrixFactory.hpp>
23#include <Xpetra_Import.hpp>
24#include <Xpetra_ImportFactory.hpp>
25#include <Xpetra_IO.hpp>
26
28
29#include "MueLu_AmalgamationFactory.hpp"
30#include "MueLu_Level.hpp"
31#include "MueLu_MasterList.hpp"
32#include "MueLu_Monitor.hpp"
33#include "MueLu_PerfUtils.hpp"
34
35namespace MueLu {
36
37template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
39
40template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42
43template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45 RCP<ParameterList> validParamList = rcp(new ParameterList());
46
47#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
48 SET_VALID_ENTRY("repartition: rebalance P and R");
49 SET_VALID_ENTRY("repartition: explicit via new copy rebalance P and R");
50 SET_VALID_ENTRY("repartition: rebalance Nullspace");
51 SET_VALID_ENTRY("transpose: use implicit");
52 SET_VALID_ENTRY("repartition: use subcommunicators");
53 SET_VALID_ENTRY("repartition: send type");
54#undef SET_VALID_ENTRY
55
56 {
57 typedef Teuchos::StringValidator validatorType;
58 RCP<validatorType> typeValidator = rcp(new validatorType(Teuchos::tuple<std::string>("Interpolation", "Restriction")));
59 validParamList->set("type", "Interpolation", "Type of the transfer operator that need to be rebalanced (Interpolation or Restriction)", typeValidator);
60 }
61
62 validParamList->set<RCP<const FactoryBase> >("P", null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
63 validParamList->set<RCP<const FactoryBase> >("R", null, "Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
64 validParamList->set<RCP<const FactoryBase> >("Nullspace", null, "Factory of the nullspace that need to be rebalanced (only used if type=Interpolation)");
65 validParamList->set<RCP<const FactoryBase> >("Coordinates", null, "Factory of the coordinates that need to be rebalanced (only used if type=Interpolation)");
66 validParamList->set<RCP<const FactoryBase> >("Material", null, "Factory of the material that need to be rebalanced (only used if type=Interpolation)");
67 validParamList->set<RCP<const FactoryBase> >("BlockNumber", null, "Factory of the block ids that need to be rebalanced (only used if type=Interpolation)");
68 validParamList->set<RCP<const FactoryBase> >("Importer", null, "Factory of the importer object used for the rebalancing");
69 validParamList->set<int>("write start", -1, "First level at which coordinates should be written to file");
70 validParamList->set<int>("write end", -1, "Last level at which coordinates should be written to file");
71
72 // TODO validation: "P" parameter valid only for type="Interpolation" and "R" valid only for type="Restriction". Like so:
73 // if (paramList.isEntry("type") && paramList.get("type) == "Interpolation) {
74 // validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
75
76 return validParamList;
77}
78
79template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81 const ParameterList& pL = GetParameterList();
82
83 if (pL.get<std::string>("type") == "Interpolation") {
84 Input(coarseLevel, "P");
85 if (pL.get<bool>("repartition: rebalance Nullspace"))
86 Input(coarseLevel, "Nullspace");
87 if (pL.get<RCP<const FactoryBase> >("Coordinates") != Teuchos::null)
88 Input(coarseLevel, "Coordinates");
89 if (pL.get<RCP<const FactoryBase> >("Material") != Teuchos::null)
90 Input(coarseLevel, "Material");
91 if (pL.get<RCP<const FactoryBase> >("BlockNumber") != Teuchos::null)
92 Input(coarseLevel, "BlockNumber");
93
94 } else {
95 if (pL.get<bool>("transpose: use implicit") == false)
96 Input(coarseLevel, "R");
97 }
98
99 Input(coarseLevel, "Importer");
100}
101
102template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
104 FactoryMonitor m(*this, "Build", coarseLevel);
105 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO> xdMV;
106
107 const ParameterList& pL = GetParameterList();
108
109 RCP<Matrix> originalP = Get<RCP<Matrix> >(coarseLevel, "P");
110 // If we don't have a valid P (e.g., # global aggregates is 0), skip this rebalancing. This level will
111 // ultimately be removed in MueLu_Hierarchy_defs.h via a resize()
112 if (originalP == Teuchos::null) {
113 Set(coarseLevel, "P", originalP);
114 return;
115 }
116 int implicit = !pL.get<bool>("repartition: rebalance P and R");
117 int reallyExplicit = pL.get<bool>("repartition: explicit via new copy rebalance P and R");
118 int writeStart = pL.get<int>("write start");
119 int writeEnd = pL.get<int>("write end");
120
121 if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel, "Coordinates")) {
122 std::string fileName = "coordinates_level_0.m";
123 RCP<xdMV> fineCoords = fineLevel.Get<RCP<xdMV> >("Coordinates");
124 if (fineCoords != Teuchos::null)
125 Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *fineCoords);
126 }
127
128 if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel, "BlockNumber")) {
129 std::string fileName = "BlockNumber_level_0.m";
130 RCP<LocalOrdinalVector> fineBlockNumber = fineLevel.Get<RCP<LocalOrdinalVector> >("BlockNumber");
131 if (fineBlockNumber != Teuchos::null)
132 Xpetra::IO<SC, LO, GO, NO>::WriteLOMV(fileName, *fineBlockNumber);
133 }
134
135 RCP<const Import> importer = Get<RCP<const Import> >(coarseLevel, "Importer");
136 if (implicit) {
137 // Save the importer, we'll need it for solve
138 coarseLevel.Set("Importer", importer, NoFactory::get());
139 }
140
141 RCP<ParameterList> params = rcp(new ParameterList());
142 if (IsPrint(Statistics2)) {
143 params->set("printLoadBalancingInfo", true);
144 params->set("printCommInfo", true);
145 }
146
147 std::string transferType = pL.get<std::string>("type");
148 if (transferType == "Interpolation") {
149 originalP = Get<RCP<Matrix> >(coarseLevel, "P");
150
151 {
152 // This line must be after the Get call
153 if (implicit || importer.is_null()) {
154 GetOStream(Runtime0) << "Using original prolongator" << std::endl;
155 Set(coarseLevel, "P", originalP);
156
157 } else {
158 SubFactoryMonitor m1(*this, "Rebalancing prolongator", coarseLevel);
159 // There are two version of an explicit rebalanced P and R.
160 // The !reallyExplicit way, is sufficient for all MueLu purposes
161 // with the exception of the CombinePFactory that needs true domain
162 // and column maps.
163 // !reallyExplicit:
164 // Rather than calling fillComplete (which would entail creating a new
165 // column map), it's sufficient to replace the domain map and importer.
166 // Note that this potentially violates the assumption that in the
167 // column map, local IDs appear before any off-rank IDs.
168 //
169 // reallyExplicit:
170 // P transfers from coarse grid to the fine grid. Here, we change
171 // the domain map (coarse) of Paccording to the new partition. The
172 // range map (fine) is kept unchanged.
173 //
174 // The domain map of P must match the range map of R. To change the
175 // domain map of P, P needs to be fillCompleted again with the new
176 // domain map. To achieve this, P is copied into a new matrix that
177 // is not fill-completed. The doImport() operation is just used
178 // here to make a copy of P: the importer is trivial and there is
179 // no data movement involved. The reordering actually happens during
180 // fillComplete() with domainMap == importer->getTargetMap().
181
182 RCP<Matrix> rebalancedP;
183 if (reallyExplicit) {
184 size_t totalMaxPerRow = 0;
185 ArrayRCP<size_t> nnzPerRow(originalP->getRowMap()->getLocalNumElements(), 0);
186 for (size_t i = 0; i < originalP->getRowMap()->getLocalNumElements(); ++i) {
187 nnzPerRow[i] = originalP->getNumEntriesInLocalRow(i);
188 if (nnzPerRow[i] > totalMaxPerRow) totalMaxPerRow = nnzPerRow[i];
189 }
190
191 rebalancedP = MatrixFactory::Build(originalP->getRowMap(), totalMaxPerRow);
192
193 {
194 RCP<Import> trivialImporter = ImportFactory::Build(originalP->getRowMap(), originalP->getRowMap());
195
196 std::string monitorLabel{"Rebalancing prolongator -- import only"};
197 if (auto sendType = pL.get<std::string>("repartition: send type"); sendType != "") {
198 auto sendTypeList = Teuchos::rcp(new Teuchos::ParameterList());
199 sendTypeList->set("Send type", sendType);
200 trivialImporter->setDistributorParameters(sendTypeList);
201 monitorLabel += " [" + sendType + "]";
202 }
203
204 SubFactoryMonitor m2(*this, monitorLabel.c_str(), coarseLevel);
205 rebalancedP->doImport(*originalP, *trivialImporter, Xpetra::INSERT);
206 }
207 rebalancedP->fillComplete(importer->getTargetMap(), originalP->getRangeMap());
208
209 } else {
210 rebalancedP = originalP;
211 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(originalP);
212 TEUCHOS_TEST_FOR_EXCEPTION(crsOp == Teuchos::null, Exceptions::BadCast, "Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
213
214 RCP<CrsMatrix> rebalancedP2 = crsOp->getCrsMatrix();
215 TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error, "Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
216
217 {
218 SubFactoryMonitor subM(*this, "Rebalancing prolongator -- fast map replacement", coarseLevel);
219
220 RCP<const Import> newImporter;
221 {
222 SubFactoryMonitor subM2(*this, "Import construction", coarseLevel);
223 newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
224 }
225 rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
226 }
227 }
229 // TODO FIXME somehow we have to transfer the striding information of the permuted domain/range maps.
230 // That is probably something for an external permutation factory
231 // if (originalP->IsView("stridedMaps"))
232 // rebalancedP->CreateView("stridedMaps", originalP);
234 if (!rebalancedP.is_null()) {
235 std::ostringstream oss;
236 oss << "P_" << coarseLevel.GetLevelID();
237 rebalancedP->setObjectLabel(oss.str());
238 }
239 Set(coarseLevel, "P", rebalancedP);
240
241 if (IsPrint(Statistics2))
242 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedP, "P (rebalanced)", params);
243 }
244 }
245
246 if (importer.is_null()) {
247 if (IsAvailable(coarseLevel, "Nullspace"))
248 Set(coarseLevel, "Nullspace", Get<RCP<MultiVector> >(coarseLevel, "Nullspace"));
249
250 if (pL.isParameter("Coordinates") && pL.get<RCP<const FactoryBase> >("Coordinates") != Teuchos::null)
251 if (IsAvailable(coarseLevel, "Coordinates"))
252 Set(coarseLevel, "Coordinates", Get<RCP<xdMV> >(coarseLevel, "Coordinates"));
253
254 if (pL.isParameter("Material") && pL.get<RCP<const FactoryBase> >("Material") != Teuchos::null)
255 if (IsAvailable(coarseLevel, "Material"))
256 Set(coarseLevel, "Material", Get<RCP<MultiVector> >(coarseLevel, "Material"));
257
258 if (pL.isParameter("BlockNumber") && pL.get<RCP<const FactoryBase> >("BlockNumber") != Teuchos::null)
259 if (IsAvailable(coarseLevel, "BlockNumber"))
260 Set(coarseLevel, "BlockNumber", Get<RCP<LocalOrdinalVector> >(coarseLevel, "BlockNumber"));
261
262 return;
263 }
264
265 RCP<const Map> subCommMap;
266 bool subCommMapConstructed = false;
267
268 if (pL.isParameter("Coordinates") &&
269 pL.get<RCP<const FactoryBase> >("Coordinates") != Teuchos::null &&
270 IsAvailable(coarseLevel, "Coordinates")) {
271 RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel, "Coordinates");
272
273 // This line must be after the Get call
274 SubFactoryMonitor subM(*this, "Rebalancing coordinates", coarseLevel);
275
276 GO numElts = coords->getMap()->getGlobalNumElements();
277 LO blkSize = Teuchos::as<LO>(importer->getSourceMap()->getGlobalNumElements() / numElts);
278
279 RCP<const Import> coordImporter;
280
281 if (blkSize == 1) {
282 coordImporter = importer;
283 } else {
284 RCP<const Map> origMap = coords->getMap();
285 std::vector<size_t> stridingInfo{Teuchos::as<size_t>(blkSize)};
286 RCP<const Map> targetMap = StridedMapFactory::Build(origMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
287 importer->getTargetMap()->getLocalElementList(), origMap->getIndexBase(), stridingInfo, origMap->getComm());
288 RCP<const Map> targetVectorMap;
289
290 AmalgamationFactory<SC, LO, GO, NO>::AmalgamateMap(rcp_dynamic_cast<const StridedMap>(targetMap), targetVectorMap);
291 coordImporter = ImportFactory::Build(origMap, targetVectorMap);
292 }
293
294 RCP<xdMV> permutedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Build(coordImporter->getTargetMap(), coords->getNumVectors(), false);
295 permutedCoords->doImport(*coords, *coordImporter, Xpetra::INSERT);
296
297 if (pL.isParameter("repartition: use subcommunicators") && pL.get<bool>("repartition: use subcommunicators")) {
298 if (blkSize == 1) {
299 if (!subCommMapConstructed) {
300 subCommMap = permutedCoords->getMap()->removeEmptyProcesses();
301 subCommMapConstructed = true;
302 }
303 permutedCoords->replaceMap(subCommMap);
304 } else {
305 permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
306 }
307 }
308
309 if (permutedCoords->getMap() == Teuchos::null)
310 permutedCoords = Teuchos::null;
311
312 Set(coarseLevel, "Coordinates", permutedCoords);
313
314 std::string fileName = "rebalanced_coordinates_level_" + toString(coarseLevel.GetLevelID()) + ".m";
315 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
316 Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *permutedCoords);
317 }
318
319 if (IsAvailable(coarseLevel, "Material")) {
320 RCP<MultiVector> material = Get<RCP<MultiVector> >(coarseLevel, "Material");
321
322 // This line must be after the Get call
323 SubFactoryMonitor subM(*this, "Rebalancing material", coarseLevel);
324
325 GO numElts = material->getMap()->getGlobalNumElements();
326 LO blkSize = Teuchos::as<LO>(importer->getSourceMap()->getGlobalNumElements() / numElts);
327
328 RCP<const Import> materialImporter;
329
330 if (blkSize == 1) {
331 materialImporter = importer;
332 } else {
333 RCP<const Map> origMap = material->getMap();
334 std::vector<size_t> stridingInfo{Teuchos::as<size_t>(blkSize)};
335 RCP<const Map> targetMap = StridedMapFactory::Build(origMap->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
336 importer->getTargetMap()->getLocalElementList(), origMap->getIndexBase(), stridingInfo, origMap->getComm());
337 RCP<const Map> targetVectorMap;
338
339 AmalgamationFactory<SC, LO, GO, NO>::AmalgamateMap(rcp_dynamic_cast<const StridedMap>(targetMap), targetVectorMap);
340 materialImporter = ImportFactory::Build(origMap, targetVectorMap);
341 }
342
343 RCP<MultiVector> permutedMaterial = MultiVectorFactory::Build(materialImporter->getTargetMap(), material->getNumVectors(), false);
344 permutedMaterial->doImport(*material, *materialImporter, Xpetra::INSERT);
345
346 if (pL.get<bool>("repartition: use subcommunicators")) {
347 if (blkSize == 1) {
348 if (!subCommMapConstructed) {
349 subCommMap = permutedMaterial->getMap()->removeEmptyProcesses();
350 subCommMapConstructed = true;
351 }
352 permutedMaterial->replaceMap(subCommMap);
353 } else {
354 permutedMaterial->replaceMap(permutedMaterial->getMap()->removeEmptyProcesses());
355 }
356 }
357
358 if (permutedMaterial->getMap() == Teuchos::null)
359 permutedMaterial = Teuchos::null;
360
361 Set(coarseLevel, "Material", permutedMaterial);
362 }
363
364 if (pL.isParameter("BlockNumber") &&
365 pL.get<RCP<const FactoryBase> >("BlockNumber") != Teuchos::null &&
366 IsAvailable(coarseLevel, "BlockNumber")) {
367 RCP<LocalOrdinalVector> BlockNumber = Get<RCP<LocalOrdinalVector> >(coarseLevel, "BlockNumber");
368
369 // This line must be after the Get call
370 SubFactoryMonitor subM(*this, "Rebalancing BlockNumber", coarseLevel);
371
372 RCP<LocalOrdinalVector> permutedBlockNumber = LocalOrdinalVectorFactory::Build(importer->getTargetMap(), false);
373 permutedBlockNumber->doImport(*BlockNumber, *importer, Xpetra::INSERT);
374
375 if (pL.isParameter("repartition: use subcommunicators") && pL.get<bool>("repartition: use subcommunicators")) {
376 if (!subCommMapConstructed) {
377 subCommMap = permutedBlockNumber->getMap()->removeEmptyProcesses();
378 subCommMapConstructed = true;
379 }
380 permutedBlockNumber->replaceMap(subCommMap);
381 }
382
383 if (permutedBlockNumber->getMap() == Teuchos::null)
384 permutedBlockNumber = Teuchos::null;
385
386 Set(coarseLevel, "BlockNumber", permutedBlockNumber);
387
388 std::string fileName = "rebalanced_BlockNumber_level_" + toString(coarseLevel.GetLevelID()) + ".m";
389 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedBlockNumber->getMap() != Teuchos::null)
390 Xpetra::IO<SC, LO, GO, NO>::WriteLOMV(fileName, *permutedBlockNumber);
391 }
392
393 if (IsAvailable(coarseLevel, "Nullspace")) {
394 RCP<MultiVector> nullspace = Get<RCP<MultiVector> >(coarseLevel, "Nullspace");
395
396 // This line must be after the Get call
397 SubFactoryMonitor subM(*this, "Rebalancing nullspace", coarseLevel);
398
399 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors(), false);
400 permutedNullspace->doImport(*nullspace, *importer, Xpetra::INSERT);
401
402 if (pL.get<bool>("repartition: use subcommunicators")) {
403 if (!subCommMapConstructed) {
404 subCommMap = permutedNullspace->getMap()->removeEmptyProcesses();
405 subCommMapConstructed = true;
406 }
407 permutedNullspace->replaceMap(subCommMap);
408 }
409
410 if (permutedNullspace->getMap() == Teuchos::null)
411 permutedNullspace = Teuchos::null;
412
413 Set(coarseLevel, "Nullspace", permutedNullspace);
414 }
415
416 } else {
417 if (!pL.get<bool>("transpose: use implicit")) {
418 RCP<Matrix> originalR = Get<RCP<Matrix> >(coarseLevel, "R");
419
420 if (implicit || importer.is_null()) {
421 GetOStream(Runtime0) << "Using original restrictor" << std::endl;
422 Set(coarseLevel, "R", originalR);
423
424 } else {
425 SubFactoryMonitor m2(*this, "Rebalancing restrictor", coarseLevel);
426
427 RCP<Matrix> rebalancedR;
428 {
429 SubFactoryMonitor subM(*this, "Rebalancing restriction -- fusedImport", coarseLevel);
430
431 RCP<Map> dummy; // meaning: use originalR's domain map.
432 Teuchos::ParameterList listLabel;
433 listLabel.set("Timer Label", "MueLu::RebalanceR-" + Teuchos::toString(coarseLevel.GetLevelID()));
434 rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(), Teuchos::rcp(&listLabel, false));
435 }
436 if (!rebalancedR.is_null()) {
437 std::ostringstream oss;
438 oss << "R_" << coarseLevel.GetLevelID();
439 rebalancedR->setObjectLabel(oss.str());
440 }
441 Set(coarseLevel, "R", rebalancedR);
442
444 // TODO FIXME somehow we have to transfer the striding information of the permuted domain/range maps.
445 // That is probably something for an external permutation factory
446 // if (originalR->IsView("stridedMaps"))
447 // rebalancedR->CreateView("stridedMaps", originalR);
449
450 if (IsPrint(Statistics2))
451 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedR, "R (rebalanced)", params);
452 }
453 }
454 }
455}
456
457} // namespace MueLu
458
459#endif // MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
static void AmalgamateMap(const Map &sourceMap, const Matrix &A, RCP< const Map > &amalgamatedMap, Array< LO > &translation)
Method to create merged map for systems of PDEs.
Exception indicating invalid cast attempted.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
static const NoFactory * get()
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
virtual ~RebalanceTransferFactory()
Destructor.
RebalanceTransferFactory()
Constructor.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Runtime0
One-liner description of what is happening.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.