46 #ifndef MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
47 #define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
50 #include <Teuchos_Tuple.hpp>
53 #include "Xpetra_MultiVectorFactory.hpp"
57 #include <Xpetra_MapFactory.hpp>
68 #include "MueLu_PerfUtils.hpp"
72 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
78 SET_VALID_ENTRY(
"repartition: explicit via new copy rebalance P and R");
82 #undef SET_VALID_ENTRY
86 RCP<validatorType> typeValidator =
rcp(
new validatorType(Teuchos::tuple<std::string>(
"Interpolation",
"Restriction"),
"type"));
87 validParamList->
set(
"type",
"Interpolation",
"Type of the transfer operator that need to be rebalanced (Interpolation or Restriction)", typeValidator);
90 validParamList->
set<
RCP<const FactoryBase> >(
"P", null,
"Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
91 validParamList->
set<
RCP<const FactoryBase> >(
"R", null,
"Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
92 validParamList->
set<
RCP<const FactoryBase> >(
"Nullspace", null,
"Factory of the nullspace that need to be rebalanced (only used if type=Interpolation)");
93 validParamList->
set<
RCP<const FactoryBase> >(
"Coordinates", null,
"Factory of the coordinates that need to be rebalanced (only used if type=Interpolation)");
94 validParamList->
set<
RCP<const FactoryBase> >(
"BlockNumber", null,
"Factory of the block ids that need to be rebalanced (only used if type=Interpolation)");
95 validParamList->
set<
RCP<const FactoryBase> >(
"Importer", null,
"Factory of the importer object used for the rebalancing");
96 validParamList->
set<
int>(
"write start", -1,
"First level at which coordinates should be written to file");
97 validParamList->
set<
int>(
"write end", -1,
"Last level at which coordinates should be written to file");
103 return validParamList;
106 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 if (pL.
get<std::string>(
"type") ==
"Interpolation") {
111 Input(coarseLevel,
"P");
112 if (pL.
get<
bool>(
"repartition: rebalance Nullspace"))
113 Input(coarseLevel,
"Nullspace");
115 Input(coarseLevel,
"Coordinates");
117 Input(coarseLevel,
"BlockNumber");
120 if (pL.
get<
bool>(
"transpose: use implicit") ==
false)
121 Input(coarseLevel,
"R");
124 Input(coarseLevel,
"Importer");
127 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 RCP<Matrix> originalP = Get<RCP<Matrix> >(coarseLevel,
"P");
137 if (originalP == Teuchos::null) {
138 Set(coarseLevel,
"P", originalP);
141 int implicit = !pL.
get<
bool>(
"repartition: rebalance P and R");
142 int reallyExplicit = pL.
get<
bool>(
"repartition: explicit via new copy rebalance P and R");
143 int writeStart = pL.
get<
int>(
"write start");
144 int writeEnd = pL.
get<
int>(
"write end");
146 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel,
"Coordinates")) {
147 std::string fileName =
"coordinates_level_0.m";
149 if (fineCoords != Teuchos::null)
153 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel,
"BlockNumber")) {
154 std::string fileName =
"BlockNumber_level_0.m";
156 if (fineBlockNumber != Teuchos::null)
168 params->
set(
"printLoadBalancingInfo",
true);
169 params->
set(
"printCommInfo",
true);
172 std::string transferType = pL.
get<std::string>(
"type");
173 if (transferType ==
"Interpolation") {
174 originalP = Get<RCP<Matrix> >(coarseLevel,
"P");
180 if (implicit || importer.
is_null()) {
181 GetOStream(
Runtime0) <<
"Using original prolongator" << std::endl;
182 Set(coarseLevel,
"P", originalP);
209 if (reallyExplicit) {
210 size_t totalMaxPerRow = 0;
211 ArrayRCP<size_t> nnzPerRow(originalP->getRowMap()->getLocalNumElements(), 0);
212 for (
size_t i = 0; i < originalP->getRowMap()->getLocalNumElements(); ++i) {
213 nnzPerRow[i] = originalP->getNumEntriesInLocalRow(i);
214 if (nnzPerRow[i] > totalMaxPerRow) totalMaxPerRow = nnzPerRow[i];
217 rebalancedP = MatrixFactory::Build(originalP->getRowMap(), totalMaxPerRow);
220 RCP<Import> trivialImporter = ImportFactory::Build(originalP->getRowMap(), originalP->getRowMap());
221 SubFactoryMonitor m2(*
this,
"Rebalancing prolongator -- import only", coarseLevel);
222 rebalancedP->doImport(*originalP, *trivialImporter,
Xpetra::INSERT);
224 rebalancedP->fillComplete(importer->getTargetMap(), originalP->getRangeMap());
227 rebalancedP = originalP;
232 TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error,
"Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
235 SubFactoryMonitor subM(*
this,
"Rebalancing prolongator -- fast map replacement", coarseLevel);
240 newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
242 rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
252 std::ostringstream oss;
253 oss <<
"P_" << coarseLevel.GetLevelID();
254 rebalancedP->setObjectLabel(oss.str());
256 Set(coarseLevel,
"P", rebalancedP);
264 if (IsAvailable(coarseLevel,
"Nullspace"))
265 Set(coarseLevel,
"Nullspace", Get<
RCP<MultiVector> >(coarseLevel,
"Nullspace"));
268 if (IsAvailable(coarseLevel,
"Coordinates"))
269 Set(coarseLevel,
"Coordinates", Get<
RCP<xdMV> >(coarseLevel,
"Coordinates"));
272 if (IsAvailable(coarseLevel,
"BlockNumber"))
280 IsAvailable(coarseLevel,
"Coordinates")) {
281 RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel,
"Coordinates");
286 LO nodeNumElts = coords->getMap()->getLocalNumElements();
289 LO myBlkSize = 0, blkSize = 0;
291 myBlkSize = importer->getSourceMap()->getLocalNumElements() / nodeNumElts;
292 MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
296 coordImporter = importer;
303 GO indexBase = origMap->getIndexBase();
306 LO numEntries = OEntries.
size() / blkSize;
308 for (LO i = 0; i < numEntries; i++)
309 Entries[i] = (OEntries[i * blkSize] - indexBase) / blkSize + indexBase;
311 RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
312 coordImporter = ImportFactory::Build(origMap, targetMap);
316 permutedCoords->doImport(*coords, *coordImporter,
Xpetra::INSERT);
318 if (pL.
isParameter(
"repartition: use subcommunicators") ==
true && pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
319 permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
321 if (permutedCoords->getMap() == Teuchos::null)
322 permutedCoords = Teuchos::null;
324 Set(coarseLevel,
"Coordinates", permutedCoords);
326 std::string fileName =
"rebalanced_coordinates_level_" +
toString(coarseLevel.GetLevelID()) +
".m";
327 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
333 IsAvailable(coarseLevel,
"BlockNumber")) {
339 RCP<LocalOrdinalVector> permutedBlockNumber = LocalOrdinalVectorFactory::Build(importer->getTargetMap(),
false);
340 permutedBlockNumber->doImport(*BlockNumber, *importer,
Xpetra::INSERT);
342 if (pL.
isParameter(
"repartition: use subcommunicators") ==
true && pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
343 permutedBlockNumber->replaceMap(permutedBlockNumber->getMap()->removeEmptyProcesses());
345 if (permutedBlockNumber->getMap() == Teuchos::null)
346 permutedBlockNumber = Teuchos::null;
348 Set(coarseLevel,
"BlockNumber", permutedBlockNumber);
350 std::string fileName =
"rebalanced_BlockNumber_level_" +
toString(coarseLevel.GetLevelID()) +
".m";
351 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedBlockNumber->getMap() != Teuchos::null)
355 if (IsAvailable(coarseLevel,
"Nullspace")) {
356 RCP<MultiVector> nullspace = Get<RCP<MultiVector> >(coarseLevel,
"Nullspace");
361 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors());
362 permutedNullspace->doImport(*nullspace, *importer,
Xpetra::INSERT);
364 if (pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
365 permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
367 if (permutedNullspace->getMap() == Teuchos::null)
368 permutedNullspace = Teuchos::null;
370 Set(coarseLevel,
"Nullspace", permutedNullspace);
374 if (pL.
get<
bool>(
"transpose: use implicit") ==
false) {
375 RCP<Matrix> originalR = Get<RCP<Matrix> >(coarseLevel,
"R");
379 if (implicit || importer.
is_null()) {
380 GetOStream(
Runtime0) <<
"Using original restrictor" << std::endl;
381 Set(coarseLevel,
"R", originalR);
386 SubFactoryMonitor subM(*
this,
"Rebalancing restriction -- fusedImport", coarseLevel);
390 listLabel.
set(
"Timer Label",
"MueLu::RebalanceR-" +
Teuchos::toString(coarseLevel.GetLevelID()));
391 rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(),
Teuchos::rcp(&listLabel,
false));
394 std::ostringstream oss;
395 oss <<
"R_" << coarseLevel.GetLevelID();
396 rebalancedR->setObjectLabel(oss.str());
398 Set(coarseLevel,
"R", rebalancedR);
416 #endif // MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
Exception indicating invalid cast attempted.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level->Get< RCP<Matrix> >("A", factory) if factory == NULL => use default factory.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
#define MueLu_maxAll(rcpComm, in, out)
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static void Write(const std::string &fileName, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &M)
One-liner description of what is happening.
static const NoFactory * get()
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Print even more statistics.
bool isParameter(const std::string &name) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
#define SET_VALID_ENTRY(name)
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
int GetLevelID() const
Return level number.
std::string toString(const T &t)