10 #ifndef MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
11 #define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
14 #include <Teuchos_Tuple.hpp>
17 #include "Xpetra_MultiVectorFactory.hpp"
21 #include <Xpetra_MapFactory.hpp>
32 #include "MueLu_PerfUtils.hpp"
36 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
39 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
42 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
46 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
48 SET_VALID_ENTRY(
"repartition: explicit via new copy rebalance P and R");
52 #undef SET_VALID_ENTRY
56 RCP<validatorType> typeValidator =
rcp(
new validatorType(Teuchos::tuple<std::string>(
"Interpolation",
"Restriction")));
57 validParamList->
set(
"type",
"Interpolation",
"Type of the transfer operator that need to be rebalanced (Interpolation or Restriction)", typeValidator);
60 validParamList->
set<
RCP<const FactoryBase> >(
"P", null,
"Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
61 validParamList->
set<
RCP<const FactoryBase> >(
"R", null,
"Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
62 validParamList->
set<
RCP<const FactoryBase> >(
"Nullspace", null,
"Factory of the nullspace that need to be rebalanced (only used if type=Interpolation)");
63 validParamList->
set<
RCP<const FactoryBase> >(
"Coordinates", null,
"Factory of the coordinates that need to be rebalanced (only used if type=Interpolation)");
64 validParamList->
set<
RCP<const FactoryBase> >(
"Material", null,
"Factory of the material that need to be rebalanced (only used if type=Interpolation)");
65 validParamList->
set<
RCP<const FactoryBase> >(
"BlockNumber", null,
"Factory of the block ids that need to be rebalanced (only used if type=Interpolation)");
66 validParamList->
set<
RCP<const FactoryBase> >(
"Importer", null,
"Factory of the importer object used for the rebalancing");
67 validParamList->
set<
int>(
"write start", -1,
"First level at which coordinates should be written to file");
68 validParamList->
set<
int>(
"write end", -1,
"Last level at which coordinates should be written to file");
74 return validParamList;
77 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
81 if (pL.
get<std::string>(
"type") ==
"Interpolation") {
82 Input(coarseLevel,
"P");
83 if (pL.
get<
bool>(
"repartition: rebalance Nullspace"))
84 Input(coarseLevel,
"Nullspace");
86 Input(coarseLevel,
"Coordinates");
88 Input(coarseLevel,
"Material");
90 Input(coarseLevel,
"BlockNumber");
93 if (pL.
get<
bool>(
"transpose: use implicit") ==
false)
94 Input(coarseLevel,
"R");
97 Input(coarseLevel,
"Importer");
100 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
107 RCP<Matrix> originalP = Get<RCP<Matrix> >(coarseLevel,
"P");
110 if (originalP == Teuchos::null) {
111 Set(coarseLevel,
"P", originalP);
114 int implicit = !pL.
get<
bool>(
"repartition: rebalance P and R");
115 int reallyExplicit = pL.
get<
bool>(
"repartition: explicit via new copy rebalance P and R");
116 int writeStart = pL.
get<
int>(
"write start");
117 int writeEnd = pL.
get<
int>(
"write end");
119 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel,
"Coordinates")) {
120 std::string fileName =
"coordinates_level_0.m";
122 if (fineCoords != Teuchos::null)
126 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel,
"BlockNumber")) {
127 std::string fileName =
"BlockNumber_level_0.m";
129 if (fineBlockNumber != Teuchos::null)
141 params->
set(
"printLoadBalancingInfo",
true);
142 params->
set(
"printCommInfo",
true);
145 std::string transferType = pL.
get<std::string>(
"type");
146 if (transferType ==
"Interpolation") {
147 originalP = Get<RCP<Matrix> >(coarseLevel,
"P");
153 if (implicit || importer.
is_null()) {
154 GetOStream(
Runtime0) <<
"Using original prolongator" << std::endl;
155 Set(coarseLevel,
"P", originalP);
182 if (reallyExplicit) {
183 size_t totalMaxPerRow = 0;
184 ArrayRCP<size_t> nnzPerRow(originalP->getRowMap()->getLocalNumElements(), 0);
185 for (
size_t i = 0; i < originalP->getRowMap()->getLocalNumElements(); ++i) {
186 nnzPerRow[i] = originalP->getNumEntriesInLocalRow(i);
187 if (nnzPerRow[i] > totalMaxPerRow) totalMaxPerRow = nnzPerRow[i];
190 rebalancedP = MatrixFactory::Build(originalP->getRowMap(), totalMaxPerRow);
193 RCP<Import> trivialImporter = ImportFactory::Build(originalP->getRowMap(), originalP->getRowMap());
194 SubFactoryMonitor m2(*
this,
"Rebalancing prolongator -- import only", coarseLevel);
195 rebalancedP->doImport(*originalP, *trivialImporter,
Xpetra::INSERT);
197 rebalancedP->fillComplete(importer->getTargetMap(), originalP->getRangeMap());
200 rebalancedP = originalP;
205 TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error,
"Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
208 SubFactoryMonitor subM(*
this,
"Rebalancing prolongator -- fast map replacement", coarseLevel);
213 newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
215 rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
225 std::ostringstream oss;
226 oss <<
"P_" << coarseLevel.GetLevelID();
227 rebalancedP->setObjectLabel(oss.str());
229 Set(coarseLevel,
"P", rebalancedP);
237 if (IsAvailable(coarseLevel,
"Nullspace"))
238 Set(coarseLevel,
"Nullspace", Get<
RCP<MultiVector> >(coarseLevel,
"Nullspace"));
241 if (IsAvailable(coarseLevel,
"Coordinates"))
242 Set(coarseLevel,
"Coordinates", Get<
RCP<xdMV> >(coarseLevel,
"Coordinates"));
245 if (IsAvailable(coarseLevel,
"Material"))
246 Set(coarseLevel,
"Material", Get<
RCP<MultiVector> >(coarseLevel,
"Material"));
249 if (IsAvailable(coarseLevel,
"BlockNumber"))
257 IsAvailable(coarseLevel,
"Coordinates")) {
258 RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel,
"Coordinates");
263 LO nodeNumElts = coords->getMap()->getLocalNumElements();
266 LO myBlkSize = 0, blkSize = 0;
268 myBlkSize = importer->getSourceMap()->getLocalNumElements() / nodeNumElts;
269 MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
273 coordImporter = importer;
280 GO indexBase = origMap->getIndexBase();
283 LO numEntries = OEntries.
size() / blkSize;
285 for (LO i = 0; i < numEntries; i++)
286 Entries[i] = (OEntries[i * blkSize] - indexBase) / blkSize + indexBase;
288 RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
289 coordImporter = ImportFactory::Build(origMap, targetMap);
293 permutedCoords->doImport(*coords, *coordImporter,
Xpetra::INSERT);
295 if (pL.
isParameter(
"repartition: use subcommunicators") ==
true && pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
296 permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
298 if (permutedCoords->getMap() == Teuchos::null)
299 permutedCoords = Teuchos::null;
301 Set(coarseLevel,
"Coordinates", permutedCoords);
303 std::string fileName =
"rebalanced_coordinates_level_" +
toString(coarseLevel.GetLevelID()) +
".m";
304 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
308 if (IsAvailable(coarseLevel,
"Material")) {
309 RCP<MultiVector> material = Get<RCP<MultiVector> >(coarseLevel,
"Material");
314 RCP<MultiVector> permutedMaterial = MultiVectorFactory::Build(importer->getTargetMap(), material->getNumVectors());
317 if (pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
318 permutedMaterial->replaceMap(permutedMaterial->getMap()->removeEmptyProcesses());
320 if (permutedMaterial->getMap() == Teuchos::null)
321 permutedMaterial = Teuchos::null;
323 Set(coarseLevel,
"Material", permutedMaterial);
328 IsAvailable(coarseLevel,
"BlockNumber")) {
334 RCP<LocalOrdinalVector> permutedBlockNumber = LocalOrdinalVectorFactory::Build(importer->getTargetMap(),
false);
335 permutedBlockNumber->doImport(*BlockNumber, *importer,
Xpetra::INSERT);
337 if (pL.
isParameter(
"repartition: use subcommunicators") ==
true && pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
338 permutedBlockNumber->replaceMap(permutedBlockNumber->getMap()->removeEmptyProcesses());
340 if (permutedBlockNumber->getMap() == Teuchos::null)
341 permutedBlockNumber = Teuchos::null;
343 Set(coarseLevel,
"BlockNumber", permutedBlockNumber);
345 std::string fileName =
"rebalanced_BlockNumber_level_" +
toString(coarseLevel.GetLevelID()) +
".m";
346 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedBlockNumber->getMap() != Teuchos::null)
350 if (IsAvailable(coarseLevel,
"Nullspace")) {
351 RCP<MultiVector> nullspace = Get<RCP<MultiVector> >(coarseLevel,
"Nullspace");
356 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors());
357 permutedNullspace->doImport(*nullspace, *importer,
Xpetra::INSERT);
359 if (pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
360 permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
362 if (permutedNullspace->getMap() == Teuchos::null)
363 permutedNullspace = Teuchos::null;
365 Set(coarseLevel,
"Nullspace", permutedNullspace);
369 if (pL.
get<
bool>(
"transpose: use implicit") ==
false) {
370 RCP<Matrix> originalR = Get<RCP<Matrix> >(coarseLevel,
"R");
374 if (implicit || importer.
is_null()) {
375 GetOStream(
Runtime0) <<
"Using original restrictor" << std::endl;
376 Set(coarseLevel,
"R", originalR);
381 SubFactoryMonitor subM(*
this,
"Rebalancing restriction -- fusedImport", coarseLevel);
385 listLabel.
set(
"Timer Label",
"MueLu::RebalanceR-" +
Teuchos::toString(coarseLevel.GetLevelID()));
386 rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(),
Teuchos::rcp(&listLabel,
false));
389 std::ostringstream oss;
390 oss <<
"R_" << coarseLevel.GetLevelID();
391 rebalancedR->setObjectLabel(oss.str());
393 Set(coarseLevel,
"R", rebalancedR);
411 #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)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static const NoFactory * get()
virtual ~RebalanceTransferFactory()
Destructor.
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)
RebalanceTransferFactory()
Constructor.
int GetLevelID() const
Return level number.
std::string toString(const T &t)