MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_RebalanceBlockRestrictionFactory_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_REBALANCEBLOCKRESTRICTIONFACTORY_DEF_HPP_
11#define MUELU_REBALANCEBLOCKRESTRICTIONFACTORY_DEF_HPP_
12
13#include <Teuchos_Tuple.hpp>
14
15#include "Xpetra_MultiVector.hpp"
16#include "Xpetra_MultiVectorFactory.hpp"
17#include "Xpetra_Vector.hpp"
18#include "Xpetra_VectorFactory.hpp"
19#include <Xpetra_Matrix.hpp>
20#include <Xpetra_BlockedCrsMatrix.hpp>
21#include <Xpetra_MapFactory.hpp>
22#include <Xpetra_MapExtractor.hpp>
23#include <Xpetra_MapExtractorFactory.hpp>
24#include <Xpetra_MatrixFactory.hpp>
25#include <Xpetra_Import.hpp>
26#include <Xpetra_ImportFactory.hpp>
27
29
30#include "MueLu_HierarchyUtils.hpp"
32#include "MueLu_Level.hpp"
33#include "MueLu_MasterList.hpp"
34#include "MueLu_Monitor.hpp"
35#include "MueLu_PerfUtils.hpp"
36
37namespace MueLu {
38
39template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
41 RCP<ParameterList> validParamList = rcp(new ParameterList());
42
43#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
44 SET_VALID_ENTRY("repartition: use subcommunicators");
45#undef SET_VALID_ENTRY
46
47 // validParamList->set< RCP<const FactoryBase> >("repartition: use subcommunicators", Teuchos::null, "test");
48
49 validParamList->set<RCP<const FactoryBase> >("R", Teuchos::null, "Factory of the restriction operator that need to be rebalanced (only used if type=Restriction)");
50 validParamList->set<RCP<const FactoryBase> >("Importer", Teuchos::null, "Generating factory of the matrix Importer for rebalancing");
51 validParamList->set<RCP<const FactoryBase> >("SubImporters", Teuchos::null, "Generating factory of the matrix sub-Importers for rebalancing");
52 validParamList->set<RCP<const FactoryBase> >("Nullspace", Teuchos::null, "Generating factory of the Nullspace operator");
53
54 return validParamList;
55}
56
57template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
59 FactManager_.push_back(FactManager);
60}
61
62template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64 Input(coarseLevel, "R");
65
66 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
67 for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
68 SetFactoryManager fineSFM(rcpFromRef(fineLevel), *it);
69 SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
70
71 if (!UseSingleSourceImporters_) coarseLevel.DeclareInput("Importer", (*it)->GetFactory("Importer").get(), this);
72 coarseLevel.DeclareInput("Nullspace", (*it)->GetFactory("Nullspace").get(), this);
73 }
74
75 // Use the non-manager path if the maps / importers are generated in one place
76 if (UseSingleSourceImporters_) {
77 Input(coarseLevel, "SubImporters");
78 }
79}
80
81template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 FactoryMonitor m(*this, "Build", coarseLevel);
84
85 RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
86
87 Teuchos::RCP<Matrix> originalTransferOp = Teuchos::null;
88 originalTransferOp = Get<RCP<Matrix> >(coarseLevel, "R");
89
90 RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > bOriginalTransferOp =
91 Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(originalTransferOp);
92 TEUCHOS_TEST_FOR_EXCEPTION(bOriginalTransferOp == Teuchos::null, Exceptions::BadCast, "MueLu::RebalanceBlockTransferFactory::Build: input matrix P or R is not of type BlockedCrsMatrix! error.");
93
94 RCP<const MapExtractor> rangeMapExtractor = bOriginalTransferOp->getRangeMapExtractor();
95 RCP<const MapExtractor> domainMapExtractor = bOriginalTransferOp->getDomainMapExtractor();
96
97 // restrict communicator?
98 bool bRestrictComm = false;
99 const ParameterList &pL = GetParameterList();
100 if (pL.get<bool>("repartition: use subcommunicators") == true)
101 bRestrictComm = true;
102
103 // check if GIDs for full maps have to be sorted:
104 // For the Thyra mode ordering they do not have to be sorted since the GIDs are
105 // numbered as 0...n1,0...,n2 (starting with zero for each subblock). The MapExtractor
106 // generates unique GIDs during the construction.
107 // For Xpetra style, the GIDs have to be reordered. Such that one obtains a ordered
108 // list of GIDs in an increasing ordering. In Xpetra, the GIDs are all unique through
109 // out all submaps.
110 bool bThyraRangeGIDs = rangeMapExtractor->getThyraMode();
111 bool bThyraDomainGIDs = domainMapExtractor->getThyraMode();
112
113 // rebuild rebalanced blocked P operator
114 std::vector<GO> fullRangeMapVector;
115 std::vector<GO> fullDomainMapVector;
116 std::vector<RCP<const Map> > subBlockRRangeMaps;
117 std::vector<RCP<const Map> > subBlockRDomainMaps;
118 subBlockRRangeMaps.reserve(bOriginalTransferOp->Rows()); // reserve size for block P operators
119 subBlockRDomainMaps.reserve(bOriginalTransferOp->Cols()); // reserve size for block P operators
120
121 std::vector<Teuchos::RCP<Matrix> > subBlockRebR;
122 subBlockRebR.reserve(bOriginalTransferOp->Cols());
123
124 std::vector<RCP<const Import> > importers = std::vector<RCP<const Import> >(bOriginalTransferOp->Rows(), Teuchos::null);
125 if (UseSingleSourceImporters_) {
126 importers = Get<std::vector<RCP<const Import> > >(coarseLevel, "SubImporters");
127 }
128
129 int curBlockId = 0;
130 Teuchos::RCP<const Import> rebalanceImporter;
131 std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
132 for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
133 // begin SubFactoryManager environment
134 SetFactoryManager fineSFM(rcpFromRef(fineLevel), *it);
135 SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
136
137 if (UseSingleSourceImporters_)
138 rebalanceImporter = importers[curBlockId];
139 else
140 rebalanceImporter = coarseLevel.Get<Teuchos::RCP<const Import> >("Importer", (*it)->GetFactory("Importer").get());
141
142 // extract matrix block
143 Teuchos::RCP<Matrix> Rii = bOriginalTransferOp->getMatrix(curBlockId, curBlockId);
144
145 // TODO run this only in the debug version
146 TEUCHOS_TEST_FOR_EXCEPTION(bThyraRangeGIDs == true && Rii->getRowMap()->getMinAllGlobalIndex() != 0,
148 "MueLu::RebalanceBlockRestrictionFactory::Build: inconsistent Thyra GIDs. Thyra global ids for block range " << curBlockId << " start with " << Rii->getRowMap()->getMinAllGlobalIndex() << " but should start with 0");
149 TEUCHOS_TEST_FOR_EXCEPTION(bThyraDomainGIDs == true && Rii->getColMap()->getMinAllGlobalIndex() != 0,
151 "MueLu::RebalanceBlockRestrictionFactory::Build: inconsistent Thyra GIDs. Thyra global ids for block domain " << curBlockId << " start with " << Rii->getColMap()->getMinAllGlobalIndex() << " but should start with 0");
152
153 Teuchos::RCP<Matrix> rebRii;
154 if (rebalanceImporter != Teuchos::null) {
155 std::stringstream ss;
156 ss << "Rebalancing restriction block R(" << curBlockId << "," << curBlockId << ")";
157 SubFactoryMonitor m1(*this, ss.str(), coarseLevel);
158 {
159 SubFactoryMonitor subM(*this, "Rebalancing restriction -- fusedImport", coarseLevel);
160 // Note: The 3rd argument says to use originalR's domain map.
161
162 RCP<Map> dummy;
163 rebRii = MatrixFactory::Build(Rii, *rebalanceImporter, dummy, rebalanceImporter->getTargetMap());
164 }
165
166 RCP<ParameterList> params = rcp(new ParameterList());
167 params->set("printLoadBalancingInfo", true);
168 std::stringstream ss2;
169 ss2 << "R(" << curBlockId << "," << curBlockId << ") rebalanced:";
170 GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*rebRii, ss2.str(), params);
171 } else {
172 rebRii = Rii;
173 RCP<ParameterList> params = rcp(new ParameterList());
174 params->set("printLoadBalancingInfo", true);
175 std::stringstream ss2;
176 ss2 << "R(" << curBlockId << "," << curBlockId << ") not rebalanced:";
177 GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*rebRii, ss2.str(), params);
178 }
179
180 // fix striding information for rebalanced diagonal block rebRii
181 Teuchos::RCP<const StridedMap> orig_stridedRgMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getMap(Teuchos::as<size_t>(curBlockId), rangeMapExtractor->getThyraMode()));
182 Teuchos::RCP<const Map> stridedRgMap = Teuchos::null;
183 if (orig_stridedRgMap != Teuchos::null) {
184 std::vector<size_t> stridingData = orig_stridedRgMap->getStridingData();
185 Teuchos::ArrayView<const GlobalOrdinal> nodeRangeMapii = rebRii->getRangeMap()->getLocalElementList();
186 stridedRgMap = StridedMapFactory::Build(
187 originalTransferOp->getRangeMap()->lib(),
188 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
189 nodeRangeMapii,
190 rebRii->getRangeMap()->getIndexBase(),
191 stridingData,
192 originalTransferOp->getRangeMap()->getComm(),
193 orig_stridedRgMap->getStridedBlockId(),
194 orig_stridedRgMap->getOffset());
195 } else
196 stridedRgMap = Rii->getRangeMap();
197
198 Teuchos::RCP<const StridedMap> orig_stridedDoMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getMap(Teuchos::as<size_t>(curBlockId), domainMapExtractor->getThyraMode()));
199 Teuchos::RCP<const Map> stridedDoMap = Teuchos::null;
200 if (orig_stridedDoMap != Teuchos::null) {
201 std::vector<size_t> stridingData = orig_stridedDoMap->getStridingData();
202 Teuchos::ArrayView<const GlobalOrdinal> nodeDomainMapii = rebRii->getDomainMap()->getLocalElementList();
203 stridedDoMap = StridedMapFactory::Build(
204 originalTransferOp->getDomainMap()->lib(),
205 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
206 nodeDomainMapii,
207 rebRii->getDomainMap()->getIndexBase(),
208 stridingData,
209 originalTransferOp->getDomainMap()->getComm(),
210 orig_stridedDoMap->getStridedBlockId(),
211 orig_stridedDoMap->getOffset());
212 } else
213 stridedDoMap = Rii->getDomainMap();
214
215 if (bRestrictComm) {
216 stridedRgMap->removeEmptyProcesses();
217 stridedDoMap->removeEmptyProcesses();
218 }
219
220 TEUCHOS_TEST_FOR_EXCEPTION(stridedRgMap == Teuchos::null, Exceptions::RuntimeError, "MueLu::RebalanceBlockRestrictionFactory::Build: failed to generate striding information. error.");
221 TEUCHOS_TEST_FOR_EXCEPTION(stridedDoMap == Teuchos::null, Exceptions::RuntimeError, "MueLu::RebalanceBlockRestrictionFactory::Build: failed to generate striding information. error.");
222
223 // replace stridedMaps view in diagonal sub block
224 if (rebRii->IsView("stridedMaps")) rebRii->RemoveView("stridedMaps");
225 rebRii->CreateView("stridedMaps", stridedRgMap, stridedDoMap);
226
227 // store rebalanced subblock
228 subBlockRebR.push_back(rebRii);
229
230 // append strided row map (= range map) to list of range maps.
231 Teuchos::RCP<const Map> rangeMapii = rebRii->getRowMap("stridedMaps");
232 subBlockRRangeMaps.push_back(rangeMapii);
233 Teuchos::ArrayView<const GlobalOrdinal> nodeRangeMapii = rebRii->getRangeMap()->getLocalElementList();
234 // append the GIDs in the end. Do not sort if we have Thyra style GIDs
235 fullRangeMapVector.insert(fullRangeMapVector.end(), nodeRangeMapii.begin(), nodeRangeMapii.end());
236 if (bThyraRangeGIDs == false)
237 sort(fullRangeMapVector.begin(), fullRangeMapVector.end());
238
239 // append strided col map (= domain map) to list of range maps.
240 Teuchos::RCP<const Map> domainMapii = rebRii->getColMap("stridedMaps");
241 subBlockRDomainMaps.push_back(domainMapii);
242 Teuchos::ArrayView<const GlobalOrdinal> nodeDomainMapii = rebRii->getDomainMap()->getLocalElementList();
243 // append the GIDs in the end. Do not sort if we have Thyra style GIDs
244 fullDomainMapVector.insert(fullDomainMapVector.end(), nodeDomainMapii.begin(), nodeDomainMapii.end());
245 if (bThyraDomainGIDs == false)
246 sort(fullDomainMapVector.begin(), fullDomainMapVector.end());
247
249
250 // rebalance null space
251 // This rebalances the null space partial vector associated with the current block (generated by the NullspaceFactory
252 // associated with the block)
253 if (rebalanceImporter != Teuchos::null) { // rebalance null space
254 std::stringstream ss2;
255 ss2 << "Rebalancing nullspace block(" << curBlockId << "," << curBlockId << ")";
256 SubFactoryMonitor subM(*this, ss2.str(), coarseLevel);
257
258 RCP<MultiVector> nullspace = coarseLevel.Get<RCP<MultiVector> >("Nullspace", (*it)->GetFactory("Nullspace").get());
259 RCP<MultiVector> permutedNullspace = MultiVectorFactory::Build(rebalanceImporter->getTargetMap(), nullspace->getNumVectors());
260 permutedNullspace->doImport(*nullspace, *rebalanceImporter, Xpetra::INSERT);
261
262 // TODO subcomm enabled everywhere or nowhere
263 if (bRestrictComm)
264 permutedNullspace->replaceMap(permutedNullspace->getMap()->removeEmptyProcesses());
265
266 coarseLevel.Set<RCP<MultiVector> >("Nullspace", permutedNullspace, (*it)->GetFactory("Nullspace").get());
267
268 } // end rebalance null space
269 else { // do nothing
270 RCP<MultiVector> nullspace = coarseLevel.Get<RCP<MultiVector> >("Nullspace", (*it)->GetFactory("Nullspace").get());
271 coarseLevel.Set<RCP<MultiVector> >("Nullspace", nullspace, (*it)->GetFactory("Nullspace").get());
272 }
273
275
276 curBlockId++;
277 } // end for loop
278
279 // extract map index base from maps of blocked P
280 GO rangeIndexBase = originalTransferOp->getRangeMap()->getIndexBase();
281 GO domainIndexBase = originalTransferOp->getDomainMap()->getIndexBase();
282
283 // check this
284 Teuchos::ArrayView<GO> fullRangeMapGIDs(fullRangeMapVector.size() ? &fullRangeMapVector[0] : 0, fullRangeMapVector.size());
285 Teuchos::RCP<const StridedMap> stridedRgFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getFullMap());
286 Teuchos::RCP<const Map> fullRangeMap = Teuchos::null;
287 if (stridedRgFullMap != Teuchos::null) {
288 std::vector<size_t> stridedData = stridedRgFullMap->getStridingData();
289 fullRangeMap =
290 StridedMapFactory::Build(
291 originalTransferOp->getRangeMap()->lib(),
292 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
293 fullRangeMapGIDs,
294 rangeIndexBase,
295 stridedData,
296 originalTransferOp->getRangeMap()->getComm(),
297 stridedRgFullMap->getStridedBlockId(),
298 stridedRgFullMap->getOffset());
299 } else {
300 fullRangeMap =
301 MapFactory::Build(
302 originalTransferOp->getRangeMap()->lib(),
303 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
304 fullRangeMapGIDs,
305 rangeIndexBase,
306 originalTransferOp->getRangeMap()->getComm());
307 }
308
309 Teuchos::ArrayView<GO> fullDomainMapGIDs(fullDomainMapVector.size() ? &fullDomainMapVector[0] : 0, fullDomainMapVector.size());
310 Teuchos::RCP<const StridedMap> stridedDoFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getFullMap());
311 Teuchos::RCP<const Map> fullDomainMap = Teuchos::null;
312 if (stridedDoFullMap != Teuchos::null) {
313 TEUCHOS_TEST_FOR_EXCEPTION(stridedDoFullMap == Teuchos::null, Exceptions::BadCast, "MueLu::BlockedPFactory::Build: full map in domain map extractor has no striding information! error.");
314 std::vector<size_t> stridedData2 = stridedDoFullMap->getStridingData();
315 fullDomainMap =
316 StridedMapFactory::Build(
317 originalTransferOp->getDomainMap()->lib(),
318 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
319 fullDomainMapGIDs,
320 domainIndexBase,
321 stridedData2,
322 originalTransferOp->getDomainMap()->getComm(),
323 stridedDoFullMap->getStridedBlockId(),
324 stridedDoFullMap->getOffset());
325 } else {
326 fullDomainMap =
327 MapFactory::Build(
328 originalTransferOp->getDomainMap()->lib(),
329 Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
330 fullDomainMapGIDs,
331 domainIndexBase,
332 originalTransferOp->getDomainMap()->getComm());
333 }
334
335 if (bRestrictComm) {
336 fullRangeMap->removeEmptyProcesses();
337 fullDomainMap->removeEmptyProcesses();
338 }
339
340 // build map extractors
341 Teuchos::RCP<const MapExtractor> rebrangeMapExtractor =
342 MapExtractorFactory::Build(fullRangeMap, subBlockRRangeMaps, bThyraRangeGIDs);
343 Teuchos::RCP<const MapExtractor> rebdomainMapExtractor =
344 MapExtractorFactory::Build(fullDomainMap, subBlockRDomainMaps, bThyraDomainGIDs);
345
346 Teuchos::RCP<BlockedCrsMatrix> bRebR = Teuchos::rcp(new BlockedCrsMatrix(rebrangeMapExtractor, rebdomainMapExtractor, 10));
347 for (size_t i = 0; i < subBlockRRangeMaps.size(); i++) {
348 Teuchos::RCP<CrsMatrixWrap> crsOpii = Teuchos::rcp_dynamic_cast<CrsMatrixWrap>(subBlockRebR[i]);
349 bRebR->setMatrix(i, i, crsOpii);
350 }
351
352 bRebR->fillComplete();
353
354 Set(coarseLevel, "R", Teuchos::rcp_dynamic_cast<Matrix>(bRebR)); // do nothing // TODO remove this!
355
356} // Build
357
358} // namespace MueLu
359
360#endif /* MUELU_REBALANCEBLOCKRESTRICTIONFACTORY_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 DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager)
Add a factory manager.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
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.
@ Statistics0
Print statistics that do not involve significant additional computation.