46 #ifndef MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
47 #define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
50 #include <Teuchos_Tuple.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))
81 #undef SET_VALID_ENTRY
85 RCP<validatorType> typeValidator =
rcp (
new validatorType(Teuchos::tuple<std::string>(
"Interpolation",
"Restriction"),
"type"));
86 validParamList->
set(
"type",
"Interpolation",
"Type of the transfer operator that need to be rebalanced (Interpolation or Restriction)", typeValidator);
89 validParamList->
set<
RCP<const FactoryBase> >(
"P", null,
"Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
90 validParamList->
set<
RCP<const FactoryBase> >(
"R", null,
"Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
91 validParamList->
set<
RCP<const FactoryBase> >(
"Nullspace", null,
"Factory of the nullspace that need to be rebalanced (only used if type=Interpolation)");
92 validParamList->
set<
RCP<const FactoryBase> >(
"Coordinates", null,
"Factory of the coordinates that need to be rebalanced (only used if type=Interpolation)");
93 validParamList->
set<
RCP<const FactoryBase> >(
"Importer", null,
"Factory of the importer object used for the rebalancing");
94 validParamList->
set<
int > (
"write start", -1,
"First level at which coordinates should be written to file");
95 validParamList->
set<
int > (
"write end", -1,
"Last level at which coordinates should be written to file");
101 return validParamList;
104 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
108 if (pL.
get<std::string>(
"type") ==
"Interpolation") {
109 Input(coarseLevel,
"P");
110 if (pL.
get<
bool>(
"repartition: rebalance Nullspace"))
111 Input(coarseLevel,
"Nullspace");
113 Input(coarseLevel,
"Coordinates");
116 if (pL.
get<
bool>(
"transpose: use implicit") ==
false)
117 Input(coarseLevel,
"R");
120 Input(coarseLevel,
"Importer");
123 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
130 int implicit = !pL.
get<
bool>(
"repartition: rebalance P and R");
131 int writeStart = pL.
get<
int> (
"write start");
132 int writeEnd = pL.
get<
int> (
"write end");
134 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel,
"Coordinates")) {
135 std::string fileName =
"coordinates_level_0.m";
137 if (fineCoords != Teuchos::null)
149 params->
set(
"printLoadBalancingInfo",
true);
150 params->
set(
"printCommInfo",
true);
153 std::string transferType = pL.
get<std::string>(
"type");
154 if (transferType ==
"Interpolation") {
155 RCP<Matrix> originalP = Get< RCP<Matrix> >(coarseLevel,
"P");
161 if (implicit || importer.
is_null()) {
162 GetOStream(
Runtime0) <<
"Using original prolongator" << std::endl;
163 Set(coarseLevel,
"P", originalP);
186 TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error,
"Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
189 SubFactoryMonitor subM(*
this,
"Rebalancing prolongator -- fast map replacement", coarseLevel);
194 newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
196 rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
205 if(!rebalancedP.
is_null()) {std::ostringstream oss; oss <<
"P_" << coarseLevel.GetLevelID(); rebalancedP->setObjectLabel(oss.str());}
206 Set(coarseLevel,
"P", rebalancedP);
214 if (IsAvailable(coarseLevel,
"Nullspace"))
215 Set(coarseLevel,
"Nullspace", Get<
RCP<MultiVector> >(coarseLevel,
"Nullspace"));
218 if (IsAvailable(coarseLevel,
"Coordinates"))
219 Set(coarseLevel,
"Coordinates", Get<
RCP<xdMV> >(coarseLevel,
"Coordinates"));
226 IsAvailable(coarseLevel,
"Coordinates")) {
227 RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel,
"Coordinates");
232 LO nodeNumElts = coords->getMap()->getNodeNumElements();
235 LO myBlkSize = 0, blkSize = 0;
237 myBlkSize = importer->getSourceMap()->getNodeNumElements() / nodeNumElts;
238 MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
242 coordImporter = importer;
249 GO indexBase = origMap->getIndexBase();
252 LO numEntries = OEntries.
size()/blkSize;
254 for (LO i = 0; i < numEntries; i++)
255 Entries[i] = (OEntries[i*blkSize]-indexBase)/blkSize + indexBase;
257 RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
258 coordImporter = ImportFactory::Build(origMap, targetMap);
262 permutedCoords->doImport(*coords, *coordImporter,
Xpetra::INSERT);
264 if (pL.
isParameter(
"repartition: use subcommunicators") ==
true && pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
265 permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
267 if (permutedCoords->getMap() == Teuchos::null)
268 permutedCoords = Teuchos::null;
270 Set(coarseLevel,
"Coordinates", permutedCoords);
272 std::string fileName =
"rebalanced_coordinates_level_" +
toString(coarseLevel.GetLevelID()) +
".m";
273 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
277 if (IsAvailable(coarseLevel,
"Nullspace")) {
278 RCP<MultiVector> nullspace = Get< RCP<MultiVector> >(coarseLevel,
"Nullspace");
283 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors());
284 permutedNullspace->doImport(*nullspace, *importer,
Xpetra::INSERT);
286 if (pL.
get<
bool>(
"repartition: use subcommunicators") ==
true)
287 permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
289 if (permutedNullspace->getMap() == Teuchos::null)
290 permutedNullspace = Teuchos::null;
292 Set(coarseLevel,
"Nullspace", permutedNullspace);
296 if (pL.
get<
bool>(
"transpose: use implicit") ==
false) {
297 RCP<Matrix> originalR = Get< RCP<Matrix> >(coarseLevel,
"R");
301 if (implicit || importer.
is_null()) {
302 GetOStream(
Runtime0) <<
"Using original restrictor" << std::endl;
303 Set(coarseLevel,
"R", originalR);
308 SubFactoryMonitor subM(*
this,
"Rebalancing restriction -- fusedImport", coarseLevel);
312 listLabel.
set(
"Timer Label",
"MueLu::RebalanceR-" +
Teuchos::toString(coarseLevel.GetLevelID()));
313 rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(),
Teuchos::rcp(&listLabel,
false));
315 if(!rebalancedR.
is_null()) {std::ostringstream oss; oss <<
"R_" << coarseLevel.GetLevelID(); rebalancedR->setObjectLabel(oss.str());}
316 Set(coarseLevel,
"R", rebalancedR);
334 #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)
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)