46 #ifndef MUELU_REBALANCETRANSFERFACTORY_DEF_HPP 47 #define MUELU_REBALANCETRANSFERFACTORY_DEF_HPP 49 #include <Teuchos_Tuple.hpp> 67 #include "MueLu_PerfUtils.hpp" 71 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 RCP<ParameterList> validParamList = rcp(
new ParameterList());
75 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name)) 79 #undef SET_VALID_ENTRY 82 typedef Teuchos::StringToIntegralParameterEntryValidator<int> validatorType;
83 RCP<validatorType> typeValidator = rcp (
new validatorType(Teuchos::tuple<std::string>(
"Interpolation",
"Restriction"),
"type"));
84 validParamList->set(
"type",
"Interpolation",
"Type of the transfer operator that need to be rebalanced (Interpolation or Restriction)", typeValidator);
87 validParamList->set< RCP<const FactoryBase> >(
"P", null,
"Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
88 validParamList->set< RCP<const FactoryBase> >(
"R", null,
"Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
89 validParamList->set< RCP<const FactoryBase> >(
"Nullspace", null,
"Factory of the nullspace that need to be rebalanced (only used if type=Restriction)");
90 validParamList->set< RCP<const FactoryBase> >(
"Coordinates", null,
"Factory of the coordinates that need to be rebalanced (only used if type=Restriction)");
91 validParamList->set< RCP<const FactoryBase> >(
"Importer", null,
"Factory of the importer object used for the rebalancing");
92 validParamList->set<
int > (
"write start", -1,
"First level at which coordinates should be written to file");
93 validParamList->set<
int > (
"write end", -1,
"Last level at which coordinates should be written to file");
99 return validParamList;
102 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
104 const ParameterList& pL = GetParameterList();
106 if (pL.get<std::string>(
"type") ==
"Interpolation") {
107 Input(coarseLevel,
"P");
108 Input(coarseLevel,
"Nullspace");
109 if (pL.get< RCP<const FactoryBase> >(
"Coordinates") != Teuchos::null)
110 Input(coarseLevel,
"Coordinates");
113 if (pL.get<
bool>(
"transpose: use implicit") ==
false)
114 Input(coarseLevel,
"R");
117 Input(coarseLevel,
"Importer");
120 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 const ParameterList& pL = GetParameterList();
127 int implicit = !pL.get<
bool>(
"repartition: rebalance P and R");
128 int writeStart = pL.get<
int> (
"write start");
129 int writeEnd = pL.get<
int> (
"write end");
131 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel,
"Coordinates")) {
132 std::string fileName =
"coordinates_level_0.m";
133 RCP<xdMV> fineCoords = fineLevel.
Get< RCP<xdMV> >(
"Coordinates");
134 if (fineCoords != Teuchos::null)
138 RCP<const Import> importer = Get<RCP<const Import> >(coarseLevel,
"Importer");
144 RCP<ParameterList> params = rcp(
new ParameterList());;
145 params->set(
"printLoadBalancingInfo",
true);
146 params->set(
"printCommInfo",
true);
148 std::string transferType = pL.get<std::string>(
"type");
149 if (transferType ==
"Interpolation") {
150 RCP<Matrix> originalP = Get< RCP<Matrix> >(coarseLevel,
"P");
156 if (implicit || importer.is_null()) {
157 GetOStream(
Runtime0) <<
"Using original prolongator" << std::endl;
158 Set(coarseLevel,
"P", originalP);
176 RCP<Matrix> rebalancedP = originalP;
177 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<
const CrsMatrixWrap>(originalP);
178 TEUCHOS_TEST_FOR_EXCEPTION(crsOp == Teuchos::null,
Exceptions::BadCast,
"Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
180 RCP<CrsMatrix> rebalancedP2 = crsOp->getCrsMatrix();
181 TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error,
"Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
184 SubFactoryMonitor subM(*
this,
"Rebalancing prolongator -- fast map replacement", coarseLevel);
186 RCP<const Import> newImporter;
189 newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
191 rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
201 Set(coarseLevel,
"P", rebalancedP);
208 if (importer.is_null()) {
209 if (IsAvailable(coarseLevel,
"Nullspace"))
210 Set(coarseLevel,
"Nullspace", Get<RCP<MultiVector> >(coarseLevel,
"Nullspace"));
212 if (pL.isParameter(
"Coordinates") && pL.get< RCP<const FactoryBase> >(
"Coordinates") != Teuchos::null)
213 if (IsAvailable(coarseLevel,
"Coordinates"))
214 Set(coarseLevel,
"Coordinates", Get< RCP<xdMV> >(coarseLevel,
"Coordinates"));
219 if (pL.isParameter(
"Coordinates") &&
220 pL.get< RCP<const FactoryBase> >(
"Coordinates") != Teuchos::null &&
221 IsAvailable(coarseLevel,
"Coordinates")) {
222 RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel,
"Coordinates");
227 LO nodeNumElts = coords->getMap()->getNodeNumElements();
230 LO myBlkSize = 0, blkSize = 0;
232 myBlkSize = importer->getSourceMap()->getNodeNumElements() / nodeNumElts;
233 MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
235 RCP<const Import> coordImporter;
237 coordImporter = importer;
243 RCP<const Map> origMap = coords->getMap();
244 GO indexBase = origMap->getIndexBase();
246 ArrayView<const GO> OEntries = importer->getTargetMap()->getNodeElementList();
247 LO numEntries = OEntries.size()/blkSize;
248 ArrayRCP<GO> Entries(numEntries);
249 for (
LO i = 0; i < numEntries; i++)
250 Entries[i] = (OEntries[i*blkSize]-indexBase)/blkSize + indexBase;
252 RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
253 coordImporter = ImportFactory::Build(origMap, targetMap);
257 permutedCoords->doImport(*coords, *coordImporter,
Xpetra::INSERT);
259 if (pL.get<
bool>(
"repartition: use subcommunicators") ==
true)
260 permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
262 Set(coarseLevel,
"Coordinates", permutedCoords);
264 std::string fileName =
"rebalanced_coordinates_level_" +
toString(coarseLevel.GetLevelID()) +
".m";
265 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
269 if (IsAvailable(coarseLevel,
"Nullspace")) {
270 RCP<MultiVector> nullspace = Get< RCP<MultiVector> >(coarseLevel,
"Nullspace");
275 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors());
276 permutedNullspace->doImport(*nullspace, *importer,
Xpetra::INSERT);
278 if (pL.get<
bool>(
"repartition: use subcommunicators") ==
true)
279 permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
281 Set(coarseLevel,
"Nullspace", permutedNullspace);
285 if (pL.get<
bool>(
"transpose: use implicit") ==
false) {
286 RCP<Matrix> originalR = Get< RCP<Matrix> >(coarseLevel,
"R");
290 if (implicit || importer.is_null()) {
291 GetOStream(
Runtime0) <<
"Using original restrictor" << std::endl;
292 Set(coarseLevel,
"R", originalR);
295 RCP<Matrix> rebalancedR;
297 SubFactoryMonitor subM(*
this,
"Rebalancing restriction -- fusedImport", coarseLevel);
300 Teuchos::ParameterList listLabel;
301 listLabel.set(
"Timer Label",
"MueLu::RebalanceR-" +
Teuchos::toString(coarseLevel.GetLevelID()));
302 rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(),Teuchos::rcp(&listLabel,
false));
304 Set(coarseLevel,
"R", rebalancedR);
322 #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.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
#define MueLu_maxAll(rcpComm, in, out)
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 Monitor but with additional timers.
static void Write(const std::string &fileName, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &M)
One-liner description of what is happening.
Namespace for MueLu classes and methods.
static const NoFactory * get()
int GetLevelID() const
Return level number.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Class that holds all level-specific information.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
#define SET_VALID_ENTRY(name)
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)