MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RebalanceBlockAcFactory_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_REBALANCEBLOCKACFACTORY_DEF_HPP_
11#define MUELU_REBALANCEBLOCKACFACTORY_DEF_HPP_
12
13#include <Xpetra_Matrix.hpp>
14#include <Xpetra_CrsMatrix.hpp>
15#include <Xpetra_CrsMatrixWrap.hpp>
16#include <Xpetra_MatrixFactory.hpp>
17#include <Xpetra_MapExtractor.hpp>
18#include <Xpetra_MapExtractorFactory.hpp>
19#include <Xpetra_StridedMap.hpp>
20#include <Xpetra_StridedMapFactory.hpp>
21#include <Xpetra_BlockedCrsMatrix.hpp>
22
23#include <Xpetra_VectorFactory.hpp>
24
26
28#include "MueLu_HierarchyUtils.hpp"
29#include "MueLu_MasterList.hpp"
30#include "MueLu_Monitor.hpp"
31#include "MueLu_PerfUtils.hpp"
32
33namespace MueLu {
34
35template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37 RCP<ParameterList> validParamList = rcp(new ParameterList());
38
39#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
40 SET_VALID_ENTRY("repartition: use subcommunicators");
41#undef SET_VALID_ENTRY
42
43 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A for rebalancing");
44 validParamList->set<RCP<const FactoryBase> >("Importer", Teuchos::null, "Generating factory of the matrix Importer for rebalancing");
45 validParamList->set<RCP<const FactoryBase> >("SubImporters", Teuchos::null, "Generating factory of the matrix sub-Importers for rebalancing");
46
47 return validParamList;
48}
49
50template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
52 FactManager_.push_back(FactManager);
53}
54
55template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57 Input(coarseLevel, "A");
58
59 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
60 for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
61 SetFactoryManager fineSFM(rcpFromRef(fineLevel), *it);
62 SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
63
64 coarseLevel.DeclareInput("Importer", (*it)->GetFactory("Importer").get(), this);
65 }
66
67 // Use the non-manager path if the maps / importers are generated in one place
68 if (FactManager_.size() == 0) {
69 Input(coarseLevel, "SubImporters");
70 }
71}
72
73template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75 FactoryMonitor m(*this, "Computing blocked Ac", coarseLevel);
76
77 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
78
79 RCP<Matrix> originalAc = Get<RCP<Matrix> >(coarseLevel, "A");
80
81 RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > bA = Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(originalAc);
82 TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::RebalanceBlockAcFactory::Build: input matrix A is not of type BlockedCrsMatrix! error.");
83 TEUCHOS_TEST_FOR_EXCEPTION(bA->Rows() != bA->Cols(), Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: Blocked operator has " << bA->Rows() << " and " << bA->Cols() << ". We only support square matrices (with same number of blocks and columns).");
84
85 // Variables to set up map extractors for blocked operators
86 std::vector<GO> fullRangeMapVector;
87 std::vector<GO> fullDomainMapVector;
88 std::vector<RCP<const Map> > subBlockARangeMaps;
89 std::vector<RCP<const Map> > subBlockADomainMaps;
90 subBlockARangeMaps.reserve(bA->Rows());
91 subBlockADomainMaps.reserve(bA->Cols());
92
93 // store map extractors
94 RCP<const MapExtractor> rangeMapExtractor = bA->getRangeMapExtractor();
95 RCP<const MapExtractor> domainMapExtractor = bA->getDomainMapExtractor();
96
97 // check if GIDs for full maps have to be sorted:
98 // For the Thyra mode ordering they do not have to be sorted since the GIDs are
99 // numbered as 0...n1,0...,n2 (starting with zero for each subblock). The MapExtractor
100 // generates unique GIDs during the construction.
101 // For Xpetra style, the GIDs have to be reordered. Such that one obtains a ordered
102 // list of GIDs in an increasing ordering. In Xpetra, the GIDs are all unique through
103 // out all submaps.
104 bool bThyraRangeGIDs = rangeMapExtractor->getThyraMode();
105 bool bThyraDomainGIDs = domainMapExtractor->getThyraMode();
106
107 // vector containing rebalanced blocks (for final output)
108 std::vector<RCP<Matrix> > subBlockRebA =
109 std::vector<RCP<Matrix> >(bA->Cols() * bA->Rows(), Teuchos::null);
110
111 // vector with Import objects from the different
112 // RepartitionFactory instances
113 std::vector<RCP<const Import> > importers = std::vector<RCP<const Import> >(bA->Rows(), Teuchos::null);
114 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
115 size_t idx = 0;
116 for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
117 SetFactoryManager fineSFM(rcpFromRef(fineLevel), *it);
118 SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
119
120 RCP<const Import> rebalanceImporter = coarseLevel.Get<RCP<const Import> >("Importer", (*it)->GetFactory("Importer").get());
121 importers[idx] = rebalanceImporter;
122 idx++;
123 }
124 if (FactManager_.size() == 0) {
125 importers = Get<std::vector<RCP<const Import> > >(coarseLevel, "SubImporters");
126 }
127
128 // restrict communicator?
129 bool bRestrictComm = false;
130 const ParameterList &pL = GetParameterList();
131 if (pL.get<bool>("repartition: use subcommunicators") == true)
132 bRestrictComm = true;
133
134 RCP<ParameterList> XpetraList = Teuchos::rcp(new ParameterList());
135 if (bRestrictComm)
136 XpetraList->set("Restrict Communicator", true);
137 else
138 XpetraList->set("Restrict Communicator", false);
139
140 // communicator for final (rebalanced) operator.
141 // If communicator is not restricted it should be identical to the communicator in bA
142 RCP<const Teuchos::Comm<int> > rebalancedComm = Teuchos::null;
143
144 // loop through all blocks and rebalance blocks
145 // Note: so far we do not support rebalancing of nested operators
146 // TODO add a check for this
147 for (size_t i = 0; i < bA->Rows(); i++) {
148 for (size_t j = 0; j < bA->Cols(); j++) {
149 // extract matrix block
150 RCP<Matrix> Aij = bA->getMatrix(i, j);
151
152 std::stringstream ss;
153 ss << "Rebalancing matrix block A(" << i << "," << j << ")";
154 SubFactoryMonitor subM(*this, ss.str(), coarseLevel);
155
156 RCP<Matrix> rebAij = Teuchos::null;
157 // General rebalancing
158 if (importers[i] != Teuchos::null &&
159 importers[j] != Teuchos::null &&
160 Aij != Teuchos::null) {
161 RCP<const Map> targetRangeMap = importers[i]->getTargetMap();
162 RCP<const Map> targetDomainMap = importers[j]->getTargetMap();
163
164 // Copy the block Aij
165 // TAW: Do we need a copy or can we do in-place rebalancing?
166 // If we do in-place rebalancing the original distribution is lost
167 // We don't really need it any more, though.
168 // RCP<Matrix> cAij = MatrixFactory::BuildCopy(Aij);
169 RCP<Matrix> cAij = Aij; // do not copy the matrix data (just an rcp pointer)
170
171 // create a new importer for column map needed for rebalanced Aij
172 Teuchos::RCP<const Import> rebAijImport = ImportFactory::Build(importers[j]->getTargetMap(), cAij->getColMap());
173 TEUCHOS_TEST_FOR_EXCEPTION(rebAijImport.is_null() == true, Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: Importer associated with block " << j << " is null.");
174
175 Teuchos::RCP<const CrsMatrixWrap> cAwij = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(cAij);
176 TEUCHOS_TEST_FOR_EXCEPTION(cAwij.is_null() == true, Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: Block (" << i << "," << j << ") is not of type CrsMatrix. We cannot rebalanced (nested) operators.");
177 Teuchos::RCP<CrsMatrix> cAmij = cAwij->getCrsMatrix();
178
179 // change domain map to rebalanced domain map (in-place). Update the importer to represent the column map
180 // cAmij->replaceDomainMapAndImporter(importers[j]->getTargetMap(),rebAijImport);
181
182 // rebalance rows of matrix block. Don't change the domain map (-> Teuchos::null)
183 // NOTE: If the communicator is restricted away, Build returns Teuchos::null.
184 rebAij = MatrixFactory::Build(cAij, *(importers[i]), *(importers[j]), targetDomainMap, targetRangeMap, XpetraList);
185 } // rebalance matrix block A(i,i)
186 else {
187 rebAij = Aij; // no rebalancing or empty block!
188 }
189
190 // store new block in output
191 subBlockRebA[i * bA->Cols() + j] = rebAij;
192
193 if (!rebAij.is_null()) {
194 // store communicator
195 if (rebalancedComm.is_null()) rebalancedComm = rebAij->getRowMap()->getComm();
196
197 // printout rebalancing information
198 RCP<ParameterList> params = rcp(new ParameterList());
199 params->set("printLoadBalancingInfo", true);
200 std::stringstream ss2;
201 ss2 << "A(" << i << "," << j << ") rebalanced:";
202 GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*rebAij, ss2.str(), params);
203 }
204 } // loop over columns j
205
206 // fix striding information of diagonal blocks
207 // Note: we do not care about the off-diagonal blocks. We just make sure, that the
208 // diagonal blocks have the corresponding striding information from the map extractors
209 // Note: the diagonal block never should be zero.
210 // TODO what if a diagonal block is Teuchos::null?
211 if (subBlockRebA[i * bA->Cols() + i].is_null() == false) {
212 RCP<Matrix> rebAii = subBlockRebA[i * bA->Cols() + i];
213 Teuchos::RCP<const StridedMap> orig_stridedRgMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getMap(i, rangeMapExtractor->getThyraMode()));
214 Teuchos::RCP<const Map> stridedRgMap = Teuchos::null;
215 if (orig_stridedRgMap != Teuchos::null) {
216 std::vector<size_t> stridingData = orig_stridedRgMap->getStridingData();
217 Teuchos::ArrayView<const GlobalOrdinal> nodeRangeMapii = rebAii->getRangeMap()->getLocalElementList();
218 stridedRgMap = StridedMapFactory::Build(
219 bA->getRangeMap()->lib(),
220 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
221 nodeRangeMapii,
222 rebAii->getRangeMap()->getIndexBase(),
223 stridingData,
224 rebalancedComm,
225 orig_stridedRgMap->getStridedBlockId(),
226 orig_stridedRgMap->getOffset());
227 }
228 Teuchos::RCP<const StridedMap> orig_stridedDoMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getMap(i, domainMapExtractor->getThyraMode()));
229 Teuchos::RCP<const Map> stridedDoMap = Teuchos::null;
230 if (orig_stridedDoMap != Teuchos::null) {
231 std::vector<size_t> stridingData = orig_stridedDoMap->getStridingData();
232 Teuchos::ArrayView<const GlobalOrdinal> nodeDomainMapii = rebAii->getDomainMap()->getLocalElementList();
233 stridedDoMap = StridedMapFactory::Build(
234 bA->getDomainMap()->lib(),
235 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
236 nodeDomainMapii,
237 rebAii->getDomainMap()->getIndexBase(),
238 stridingData,
239 rebalancedComm,
240 orig_stridedDoMap->getStridedBlockId(),
241 orig_stridedDoMap->getOffset());
242 }
243
244 if (bRestrictComm) {
245 stridedRgMap->removeEmptyProcesses();
246 stridedDoMap->removeEmptyProcesses();
247 }
248
249 TEUCHOS_TEST_FOR_EXCEPTION(stridedRgMap == Teuchos::null, Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: failed to generate striding information. error.");
250 TEUCHOS_TEST_FOR_EXCEPTION(stridedDoMap == Teuchos::null, Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: failed to generate striding information. error.");
251
252 // replace stridedMaps view in diagonal sub block
253 if (rebAii->IsView("stridedMaps")) rebAii->RemoveView("stridedMaps");
254 rebAii->CreateView("stridedMaps", stridedRgMap, stridedDoMap);
255 // collect Xpetra-based global row ids for map extractors
256 subBlockARangeMaps.push_back(rebAii->getRowMap("stridedMaps"));
257 Teuchos::ArrayView<const GlobalOrdinal> nodeRangeMap = rebAii->getRangeMap()->getLocalElementList();
258 // append the GIDs in the end. Do not sort if we have Thyra style GIDs
259 fullRangeMapVector.insert(fullRangeMapVector.end(), nodeRangeMap.begin(), nodeRangeMap.end());
260 if (bThyraRangeGIDs == false)
261 sort(fullRangeMapVector.begin(), fullRangeMapVector.end());
262
263 subBlockADomainMaps.push_back(rebAii->getColMap("stridedMaps"));
264 Teuchos::ArrayView<const GlobalOrdinal> nodeDomainMap = rebAii->getDomainMap()->getLocalElementList();
265 // append the GIDs in the end. Do not sort if we have Thyra style GIDs
266 fullDomainMapVector.insert(fullDomainMapVector.end(), nodeDomainMap.begin(), nodeDomainMap.end());
267 if (bThyraDomainGIDs == false)
268 sort(fullDomainMapVector.begin(), fullDomainMapVector.end());
269 } // end if rebAii != Teuchos::null
270 } // loop over rows i
271
272 // all sub blocks are rebalanced (if available)
273
274 // Short cut if this processor is not in the list of active processors
275 if (rebalancedComm == Teuchos::null) {
276 GetOStream(Debug, -1) << "RebalanceBlockedAc: deactivate proc " << originalAc->getRowMap()->getComm()->getRank() << std::endl;
277 // TAW: it is important that we create a dummy object of type BlockedCrsMatrix (even if we set it to Teuchos::null)
278 Teuchos::RCP<BlockedCrsMatrix> reb_bA = Teuchos::null;
279 coarseLevel.Set("A", Teuchos::rcp_dynamic_cast<Matrix>(reb_bA), this);
280 return;
281 }
282
283 // now, subBlockRebA contains all rebalanced matrix blocks
284 // extract map index base from maps of blocked A
285 GO rangeIndexBase = bA->getRangeMap()->getIndexBase();
286 GO domainIndexBase = bA->getDomainMap()->getIndexBase();
287
288 Teuchos::ArrayView<GO> fullRangeMapGIDs(fullRangeMapVector.size() ? &fullRangeMapVector[0] : 0, fullRangeMapVector.size());
289 Teuchos::RCP<const StridedMap> stridedRgFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getFullMap());
290 Teuchos::RCP<const Map> fullRangeMap = Teuchos::null;
291 if (stridedRgFullMap != Teuchos::null) {
292 std::vector<size_t> stridedData = stridedRgFullMap->getStridingData();
293 fullRangeMap =
294 StridedMapFactory::Build(
295 bA->getRangeMap()->lib(),
296 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
297 fullRangeMapGIDs,
298 rangeIndexBase,
299 stridedData,
300 rebalancedComm,
301 stridedRgFullMap->getStridedBlockId(),
302 stridedRgFullMap->getOffset());
303 } else {
304 fullRangeMap =
305 MapFactory::Build(
306 bA->getRangeMap()->lib(),
307 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
308 fullRangeMapGIDs,
309 rangeIndexBase,
310 rebalancedComm);
311 }
312 Teuchos::ArrayView<GO> fullDomainMapGIDs(fullDomainMapVector.size() ? &fullDomainMapVector[0] : 0, fullDomainMapVector.size());
313
314 Teuchos::RCP<const StridedMap> stridedDoFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getFullMap());
315 Teuchos::RCP<const Map> fullDomainMap = Teuchos::null;
316 if (stridedDoFullMap != Teuchos::null) {
317 TEUCHOS_TEST_FOR_EXCEPTION(stridedDoFullMap == Teuchos::null, Exceptions::BadCast, "MueLu::RebalanceBlockedAc::Build: full map in domain map extractor has no striding information! error.");
318 std::vector<size_t> stridedData2 = stridedDoFullMap->getStridingData();
319 fullDomainMap =
320 StridedMapFactory::Build(
321 bA->getDomainMap()->lib(),
322 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
323 fullDomainMapGIDs,
324 domainIndexBase,
325 stridedData2,
326 rebalancedComm,
327 stridedDoFullMap->getStridedBlockId(),
328 stridedDoFullMap->getOffset());
329 } else {
330 fullDomainMap =
331 MapFactory::Build(
332 bA->getDomainMap()->lib(),
333 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
334 fullDomainMapGIDs,
335 domainIndexBase,
336 rebalancedComm);
337 }
338
339 if (bRestrictComm) {
340 fullRangeMap->removeEmptyProcesses();
341 fullDomainMap->removeEmptyProcesses();
342 }
343
344 // build map extractors
345 RCP<const MapExtractor> rebRangeMapExtractor = MapExtractorFactory::Build(fullRangeMap, subBlockARangeMaps, bThyraRangeGIDs);
346 RCP<const MapExtractor> rebDomainMapExtractor = MapExtractorFactory::Build(fullDomainMap, subBlockADomainMaps, bThyraDomainGIDs);
347
348 TEUCHOS_TEST_FOR_EXCEPTION(rangeMapExtractor->NumMaps() != rebRangeMapExtractor->NumMaps(), Exceptions::RuntimeError,
349 "MueLu::RebalanceBlockedAc::Build: Rebalanced RangeMapExtractor has " << rebRangeMapExtractor->NumMaps()
350 << " sub maps. Original RangeMapExtractor has " << rangeMapExtractor->NumMaps() << ". They must match!");
351 TEUCHOS_TEST_FOR_EXCEPTION(domainMapExtractor->NumMaps() != rebDomainMapExtractor->NumMaps(), Exceptions::RuntimeError,
352 "MueLu::RebalanceBlockedAc::Build: Rebalanced DomainMapExtractor has " << rebDomainMapExtractor->NumMaps()
353 << " sub maps. Original DomainMapExtractor has " << domainMapExtractor->NumMaps() << ". They must match!");
354
355 Teuchos::RCP<BlockedCrsMatrix> reb_bA = Teuchos::rcp(new BlockedCrsMatrix(rebRangeMapExtractor, rebDomainMapExtractor, 10));
356 for (size_t i = 0; i < bA->Rows(); i++) {
357 for (size_t j = 0; j < bA->Cols(); j++) {
358 reb_bA->setMatrix(i, j, subBlockRebA[i * bA->Cols() + j]);
359 }
360 }
361
362 reb_bA->fillComplete();
363
364 coarseLevel.Set("A", Teuchos::rcp_dynamic_cast<Matrix>(reb_bA), this);
365 // rebalance additional data:
366 // be aware, that we just call the rebalance factories without switching to local
367 // factory managers, i.e. the rebalance factories have to be defined with the appropriate
368 // factories by the user!
369 if (rebalanceFacts_.begin() != rebalanceFacts_.end()) {
370 SubFactoryMonitor m2(*this, "Rebalance additional data", coarseLevel);
371
372 // call Build of all user-given transfer factories
373 for (std::vector<RCP<const FactoryBase> >::const_iterator it2 = rebalanceFacts_.begin(); it2 != rebalanceFacts_.end(); ++it2) {
374 GetOStream(Runtime0) << "RebalanceBlockedAc: call rebalance factory " << (*it2).get() << ": " << (*it2)->description() << std::endl;
375 (*it2)->CallBuild(coarseLevel);
376 }
377 }
378} // Build()
379
380template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
382 rebalanceFacts_.push_back(factory);
383} // AddRebalanceFactory()
384
385} // namespace MueLu
386
387#endif /* MUELU_REBALANCEBLOCKACFACTORY_DEF_HPP_ */
#define SET_VALID_ENTRY(name)
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
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 std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager)
Add a factory manager.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void AddRebalanceFactory(const RCP< const FactoryBase > &factory)
Add rebalancing factory in the end of list of rebalancing factories in RebalanceAcFactory.
An exception safe way to call the method 'Level::SetFactoryManager()'.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Debug
Print additional debugging information.
@ Runtime0
One-liner description of what is happening.
@ Statistics0
Print statistics that do not involve significant additional computation.