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_Level.hpp"
30#include "MueLu_MasterList.hpp"
31#include "MueLu_Monitor.hpp"
32#include "MueLu_PerfUtils.hpp"
33
34namespace MueLu {
35
36template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38
39template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
41
42template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
44 RCP<ParameterList> validParamList = rcp(new ParameterList());
45
46#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
47 SET_VALID_ENTRY("repartition: rebalance P and R");
48 SET_VALID_ENTRY("repartition: explicit via new copy rebalance P and R");
49 SET_VALID_ENTRY("repartition: rebalance Nullspace");
50 SET_VALID_ENTRY("transpose: use implicit");
51 SET_VALID_ENTRY("repartition: use subcommunicators");
52#undef SET_VALID_ENTRY
53
54 {
55 typedef Teuchos::StringValidator validatorType;
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);
58 }
59
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");
69
70 // TODO validation: "P" parameter valid only for type="Interpolation" and "R" valid only for type="Restriction". Like so:
71 // if (paramList.isEntry("type") && paramList.get("type) == "Interpolation) {
72 // validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
73
74 return validParamList;
75}
76
77template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79 const ParameterList& pL = GetParameterList();
80
81 if (pL.get<std::string>("type") == "Interpolation") {
82 Input(coarseLevel, "P");
83 if (pL.get<bool>("repartition: rebalance Nullspace"))
84 Input(coarseLevel, "Nullspace");
85 if (pL.get<RCP<const FactoryBase> >("Coordinates") != Teuchos::null)
86 Input(coarseLevel, "Coordinates");
87 if (pL.get<RCP<const FactoryBase> >("Material") != Teuchos::null)
88 Input(coarseLevel, "Material");
89 if (pL.get<RCP<const FactoryBase> >("BlockNumber") != Teuchos::null)
90 Input(coarseLevel, "BlockNumber");
91
92 } else {
93 if (pL.get<bool>("transpose: use implicit") == false)
94 Input(coarseLevel, "R");
95 }
96
97 Input(coarseLevel, "Importer");
98}
99
100template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
102 FactoryMonitor m(*this, "Build", coarseLevel);
103 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO> xdMV;
104
105 const ParameterList& pL = GetParameterList();
106
107 RCP<Matrix> originalP = Get<RCP<Matrix> >(coarseLevel, "P");
108 // If we don't have a valid P (e.g., # global aggregates is 0), skip this rebalancing. This level will
109 // ultimately be removed in MueLu_Hierarchy_defs.h via a resize()
110 if (originalP == Teuchos::null) {
111 Set(coarseLevel, "P", originalP);
112 return;
113 }
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");
118
119 if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel, "Coordinates")) {
120 std::string fileName = "coordinates_level_0.m";
121 RCP<xdMV> fineCoords = fineLevel.Get<RCP<xdMV> >("Coordinates");
122 if (fineCoords != Teuchos::null)
123 Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *fineCoords);
124 }
125
126 if (writeStart == 0 && fineLevel.GetLevelID() == 0 && writeStart <= writeEnd && IsAvailable(fineLevel, "BlockNumber")) {
127 std::string fileName = "BlockNumber_level_0.m";
128 RCP<LocalOrdinalVector> fineBlockNumber = fineLevel.Get<RCP<LocalOrdinalVector> >("BlockNumber");
129 if (fineBlockNumber != Teuchos::null)
130 Xpetra::IO<SC, LO, GO, NO>::WriteLOMV(fileName, *fineBlockNumber);
131 }
132
133 RCP<const Import> importer = Get<RCP<const Import> >(coarseLevel, "Importer");
134 if (implicit) {
135 // Save the importer, we'll need it for solve
136 coarseLevel.Set("Importer", importer, NoFactory::get());
137 }
138
139 RCP<ParameterList> params = rcp(new ParameterList());
140 if (IsPrint(Statistics2)) {
141 params->set("printLoadBalancingInfo", true);
142 params->set("printCommInfo", true);
143 }
144
145 std::string transferType = pL.get<std::string>("type");
146 if (transferType == "Interpolation") {
147 originalP = Get<RCP<Matrix> >(coarseLevel, "P");
148
149 {
150 // This line must be after the Get call
151 SubFactoryMonitor m1(*this, "Rebalancing prolongator", coarseLevel);
152
153 if (implicit || importer.is_null()) {
154 GetOStream(Runtime0) << "Using original prolongator" << std::endl;
155 Set(coarseLevel, "P", originalP);
156
157 } else {
158 // There are two version of an explicit rebalanced P and R.
159 // The !reallyExplicit way, is sufficient for all MueLu purposes
160 // with the exception of the CombinePFactory that needs true domain
161 // and column maps.
162 // !reallyExplicit:
163 // Rather than calling fillComplete (which would entail creating a new
164 // column map), it's sufficient to replace the domain map and importer.
165 // Note that this potentially violates the assumption that in the
166 // column map, local IDs appear before any off-rank IDs.
167 //
168 // reallyExplicit:
169 // P transfers from coarse grid to the fine grid. Here, we change
170 // the domain map (coarse) of Paccording to the new partition. The
171 // range map (fine) is kept unchanged.
172 //
173 // The domain map of P must match the range map of R. To change the
174 // domain map of P, P needs to be fillCompleted again with the new
175 // domain map. To achieve this, P is copied into a new matrix that
176 // is not fill-completed. The doImport() operation is just used
177 // here to make a copy of P: the importer is trivial and there is
178 // no data movement involved. The reordering actually happens during
179 // fillComplete() with domainMap == importer->getTargetMap().
180
181 RCP<Matrix> rebalancedP;
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];
188 }
189
190 rebalancedP = MatrixFactory::Build(originalP->getRowMap(), totalMaxPerRow);
191
192 {
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);
196 }
197 rebalancedP->fillComplete(importer->getTargetMap(), originalP->getRangeMap());
198
199 } else {
200 rebalancedP = originalP;
201 RCP<const CrsMatrixWrap> crsOp = rcp_dynamic_cast<const CrsMatrixWrap>(originalP);
202 TEUCHOS_TEST_FOR_EXCEPTION(crsOp == Teuchos::null, Exceptions::BadCast, "Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
203
204 RCP<CrsMatrix> rebalancedP2 = crsOp->getCrsMatrix();
205 TEUCHOS_TEST_FOR_EXCEPTION(rebalancedP2 == Teuchos::null, std::runtime_error, "Xpetra::CrsMatrixWrap doesn't have a CrsMatrix");
206
207 {
208 SubFactoryMonitor subM(*this, "Rebalancing prolongator -- fast map replacement", coarseLevel);
209
210 RCP<const Import> newImporter;
211 {
212 SubFactoryMonitor subM2(*this, "Import construction", coarseLevel);
213 newImporter = ImportFactory::Build(importer->getTargetMap(), rebalancedP->getColMap());
214 }
215 rebalancedP2->replaceDomainMapAndImporter(importer->getTargetMap(), newImporter);
216 }
217 }
219 // TODO FIXME somehow we have to transfer the striding information of the permuted domain/range maps.
220 // That is probably something for an external permutation factory
221 // if (originalP->IsView("stridedMaps"))
222 // rebalancedP->CreateView("stridedMaps", originalP);
224 if (!rebalancedP.is_null()) {
225 std::ostringstream oss;
226 oss << "P_" << coarseLevel.GetLevelID();
227 rebalancedP->setObjectLabel(oss.str());
228 }
229 Set(coarseLevel, "P", rebalancedP);
230
231 if (IsPrint(Statistics2))
232 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedP, "P (rebalanced)", params);
233 }
234 }
235
236 if (importer.is_null()) {
237 if (IsAvailable(coarseLevel, "Nullspace"))
238 Set(coarseLevel, "Nullspace", Get<RCP<MultiVector> >(coarseLevel, "Nullspace"));
239
240 if (pL.isParameter("Coordinates") && pL.get<RCP<const FactoryBase> >("Coordinates") != Teuchos::null)
241 if (IsAvailable(coarseLevel, "Coordinates"))
242 Set(coarseLevel, "Coordinates", Get<RCP<xdMV> >(coarseLevel, "Coordinates"));
243
244 if (pL.isParameter("Material") && pL.get<RCP<const FactoryBase> >("Material") != Teuchos::null)
245 if (IsAvailable(coarseLevel, "Material"))
246 Set(coarseLevel, "Material", Get<RCP<MultiVector> >(coarseLevel, "Material"));
247
248 if (pL.isParameter("BlockNumber") && pL.get<RCP<const FactoryBase> >("BlockNumber") != Teuchos::null)
249 if (IsAvailable(coarseLevel, "BlockNumber"))
250 Set(coarseLevel, "BlockNumber", Get<RCP<LocalOrdinalVector> >(coarseLevel, "BlockNumber"));
251
252 return;
253 }
254
255 if (pL.isParameter("Coordinates") &&
256 pL.get<RCP<const FactoryBase> >("Coordinates") != Teuchos::null &&
257 IsAvailable(coarseLevel, "Coordinates")) {
258 RCP<xdMV> coords = Get<RCP<xdMV> >(coarseLevel, "Coordinates");
259
260 // This line must be after the Get call
261 SubFactoryMonitor subM(*this, "Rebalancing coordinates", coarseLevel);
262
263 LO nodeNumElts = coords->getMap()->getLocalNumElements();
264
265 // If a process has no matrix rows, then we can't calculate blocksize using the formula below.
266 LO myBlkSize = 0, blkSize = 0;
267 if (nodeNumElts > 0)
268 myBlkSize = importer->getSourceMap()->getLocalNumElements() / nodeNumElts;
269 MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
270
271 RCP<const Import> coordImporter;
272 if (blkSize == 1) {
273 coordImporter = importer;
274
275 } else {
276 // NOTE: there is an implicit assumption here: we assume that dof any node are enumerated consequently
277 // Proper fix would require using decomposition similar to how we construct importer in the
278 // RepartitionFactory
279 RCP<const Map> origMap = coords->getMap();
280 GO indexBase = origMap->getIndexBase();
281
282 ArrayView<const GO> OEntries = importer->getTargetMap()->getLocalElementList();
283 LO numEntries = OEntries.size() / blkSize;
284 ArrayRCP<GO> Entries(numEntries);
285 for (LO i = 0; i < numEntries; i++)
286 Entries[i] = (OEntries[i * blkSize] - indexBase) / blkSize + indexBase;
287
288 RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
289 coordImporter = ImportFactory::Build(origMap, targetMap);
290 }
291
292 RCP<xdMV> permutedCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Build(coordImporter->getTargetMap(), coords->getNumVectors());
293 permutedCoords->doImport(*coords, *coordImporter, Xpetra::INSERT);
294
295 if (pL.isParameter("repartition: use subcommunicators") == true && pL.get<bool>("repartition: use subcommunicators") == true)
296 permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
297
298 if (permutedCoords->getMap() == Teuchos::null)
299 permutedCoords = Teuchos::null;
300
301 Set(coarseLevel, "Coordinates", permutedCoords);
302
303 std::string fileName = "rebalanced_coordinates_level_" + toString(coarseLevel.GetLevelID()) + ".m";
304 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedCoords->getMap() != Teuchos::null)
305 Xpetra::IO<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>::Write(fileName, *permutedCoords);
306 }
307
308 if (IsAvailable(coarseLevel, "Material")) {
309 RCP<MultiVector> material = Get<RCP<MultiVector> >(coarseLevel, "Material");
310
311 // This line must be after the Get call
312 SubFactoryMonitor subM(*this, "Rebalancing material", coarseLevel);
313
314 RCP<MultiVector> permutedMaterial = MultiVectorFactory::Build(importer->getTargetMap(), material->getNumVectors());
315 permutedMaterial->doImport(*material, *importer, Xpetra::INSERT);
316
317 if (pL.get<bool>("repartition: use subcommunicators") == true)
318 permutedMaterial->replaceMap(permutedMaterial->getMap()->removeEmptyProcesses());
319
320 if (permutedMaterial->getMap() == Teuchos::null)
321 permutedMaterial = Teuchos::null;
322
323 Set(coarseLevel, "Material", permutedMaterial);
324 }
325
326 if (pL.isParameter("BlockNumber") &&
327 pL.get<RCP<const FactoryBase> >("BlockNumber") != Teuchos::null &&
328 IsAvailable(coarseLevel, "BlockNumber")) {
329 RCP<LocalOrdinalVector> BlockNumber = Get<RCP<LocalOrdinalVector> >(coarseLevel, "BlockNumber");
330
331 // This line must be after the Get call
332 SubFactoryMonitor subM(*this, "Rebalancing BlockNumber", coarseLevel);
333
334 RCP<LocalOrdinalVector> permutedBlockNumber = LocalOrdinalVectorFactory::Build(importer->getTargetMap(), false);
335 permutedBlockNumber->doImport(*BlockNumber, *importer, Xpetra::INSERT);
336
337 if (pL.isParameter("repartition: use subcommunicators") == true && pL.get<bool>("repartition: use subcommunicators") == true)
338 permutedBlockNumber->replaceMap(permutedBlockNumber->getMap()->removeEmptyProcesses());
339
340 if (permutedBlockNumber->getMap() == Teuchos::null)
341 permutedBlockNumber = Teuchos::null;
342
343 Set(coarseLevel, "BlockNumber", permutedBlockNumber);
344
345 std::string fileName = "rebalanced_BlockNumber_level_" + toString(coarseLevel.GetLevelID()) + ".m";
346 if (writeStart <= coarseLevel.GetLevelID() && coarseLevel.GetLevelID() <= writeEnd && permutedBlockNumber->getMap() != Teuchos::null)
347 Xpetra::IO<SC, LO, GO, NO>::WriteLOMV(fileName, *permutedBlockNumber);
348 }
349
350 if (IsAvailable(coarseLevel, "Nullspace")) {
351 RCP<MultiVector> nullspace = Get<RCP<MultiVector> >(coarseLevel, "Nullspace");
352
353 // This line must be after the Get call
354 SubFactoryMonitor subM(*this, "Rebalancing nullspace", coarseLevel);
355
356 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(importer->getTargetMap(), nullspace->getNumVectors());
357 permutedNullspace->doImport(*nullspace, *importer, Xpetra::INSERT);
358
359 if (pL.get<bool>("repartition: use subcommunicators") == true)
360 permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
361
362 if (permutedNullspace->getMap() == Teuchos::null)
363 permutedNullspace = Teuchos::null;
364
365 Set(coarseLevel, "Nullspace", permutedNullspace);
366 }
367
368 } else {
369 if (pL.get<bool>("transpose: use implicit") == false) {
370 RCP<Matrix> originalR = Get<RCP<Matrix> >(coarseLevel, "R");
371
372 SubFactoryMonitor m2(*this, "Rebalancing restrictor", coarseLevel);
373
374 if (implicit || importer.is_null()) {
375 GetOStream(Runtime0) << "Using original restrictor" << std::endl;
376 Set(coarseLevel, "R", originalR);
377
378 } else {
379 RCP<Matrix> rebalancedR;
380 {
381 SubFactoryMonitor subM(*this, "Rebalancing restriction -- fusedImport", coarseLevel);
382
383 RCP<Map> dummy; // meaning: use originalR's domain map.
384 Teuchos::ParameterList listLabel;
385 listLabel.set("Timer Label", "MueLu::RebalanceR-" + Teuchos::toString(coarseLevel.GetLevelID()));
386 rebalancedR = MatrixFactory::Build(originalR, *importer, dummy, importer->getTargetMap(), Teuchos::rcp(&listLabel, false));
387 }
388 if (!rebalancedR.is_null()) {
389 std::ostringstream oss;
390 oss << "R_" << coarseLevel.GetLevelID();
391 rebalancedR->setObjectLabel(oss.str());
392 }
393 Set(coarseLevel, "R", rebalancedR);
394
396 // TODO FIXME somehow we have to transfer the striding information of the permuted domain/range maps.
397 // That is probably something for an external permutation factory
398 // if (originalR->IsView("stridedMaps"))
399 // rebalancedR->CreateView("stridedMaps", originalR);
401
402 if (IsPrint(Statistics2))
403 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedR, "R (rebalanced)", params);
404 }
405 }
406 }
407}
408
409} // namespace MueLu
410
411#endif // MUELU_REBALANCETRANSFERFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
#define MueLu_maxAll(rcpComm, in, out)
Exception indicating invalid cast attempted.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
bool IsAvailable(Level &level, const std::string &varName) const
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()
virtual const Teuchos::ParameterList & GetParameterList() const
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.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
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.