MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Hierarchy_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_HIERARCHY_DEF_HPP
11#define MUELU_HIERARCHY_DEF_HPP
12
13#include <time.h>
14
15#include <algorithm>
16#include <sstream>
17#include <list>
18
19#include <Xpetra_Matrix.hpp>
20#include <Xpetra_MultiVectorFactory.hpp>
21#include <Xpetra_Operator.hpp>
22#include <Xpetra_IO.hpp>
23
25
26#include "MueLu_FactoryManager.hpp"
27#include "MueLu_HierarchyUtils.hpp"
28#include "MueLu_TopRAPFactory.hpp"
29#include "MueLu_TopSmootherFactory.hpp"
30#include "MueLu_Level.hpp"
31#include "MueLu_Monitor.hpp"
32#include "MueLu_PerfUtils.hpp"
33#include "MueLu_PFactory.hpp"
34#include "MueLu_SmootherFactory.hpp"
36#include "MueLu_Behavior.hpp"
37
38#include "Teuchos_TimeMonitor.hpp"
39
40namespace MueLu {
41
42template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
44 : maxCoarseSize_(GetDefaultMaxCoarseSize())
45 , implicitTranspose_(GetDefaultImplicitTranspose())
46 , fuseProlongationAndUpdate_(GetDefaultFuseProlongationAndUpdate())
47 , doPRrebalance_(GetDefaultPRrebalance())
48 , doPRViaCopyrebalance_(false)
49 , isPreconditioner_(true)
50 , Cycle_(GetDefaultCycle())
51 , WCycleStartLevel_(0)
52 , scalingFactor_(Teuchos::ScalarTraits<double>::one())
53 , lib_(Xpetra::UseTpetra)
54 , isDumpingEnabled_(false)
55 , dumpLevel_(-2)
56 , rate_(-1)
57 , sizeOfAllocatedLevelMultiVectors_(0) {
58 AddLevel(rcp(new Level));
59}
60
61template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63 : Hierarchy() {
64 SetLabel(label);
65 setObjectLabel(label);
66 Levels_[0]->setObjectLabel(label);
67}
68
69template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71 : maxCoarseSize_(GetDefaultMaxCoarseSize())
72 , implicitTranspose_(GetDefaultImplicitTranspose())
73 , fuseProlongationAndUpdate_(GetDefaultFuseProlongationAndUpdate())
74 , doPRrebalance_(GetDefaultPRrebalance())
75 , doPRViaCopyrebalance_(false)
76 , isPreconditioner_(true)
77 , Cycle_(GetDefaultCycle())
78 , WCycleStartLevel_(0)
79 , scalingFactor_(Teuchos::ScalarTraits<double>::one())
80 , isDumpingEnabled_(false)
81 , dumpLevel_(-2)
82 , rate_(-1)
83 , sizeOfAllocatedLevelMultiVectors_(0) {
84 lib_ = A->getDomainMap()->lib();
85
86 RCP<Level> Finest = rcp(new Level);
87 AddLevel(Finest);
88
89 Finest->Set("A", A);
90}
91
92template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
93Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Hierarchy(const RCP<Matrix>& A, const std::string& label)
94 : Hierarchy(A) {
95 setObjectLabel(label);
96 Levels_[0]->setObjectLabel(label);
97}
98
99template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
101
102template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
104 int levelID = LastLevelID() + 1; // ID of the inserted level
105
106 if (level->GetLevelID() != -1 && (level->GetLevelID() != levelID))
107 GetOStream(Warnings1) << "Hierarchy::AddLevel(): Level with ID=" << level->GetLevelID() << " have been added at the end of the hierarchy\n but its ID have been redefined"
108 << " because last level ID of the hierarchy was " << LastLevelID() << "." << std::endl;
109
110 Levels_.push_back(level);
111 level->SetLevelID(levelID);
112 level->setlib(lib_);
113
114 level->SetPreviousLevel((levelID == 0) ? Teuchos::null : Levels_[LastLevelID() - 1]);
115 level->setObjectLabel(this->getObjectLabel());
116}
117
118template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
120 RCP<Level> newLevel = Levels_[LastLevelID()]->Build(); // new coarse level, using copy constructor
121 newLevel->setlib(lib_);
122 this->AddLevel(newLevel); // add to hierarchy
123}
124
125template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
127 TEUCHOS_TEST_FOR_EXCEPTION(levelID < 0 || levelID > LastLevelID(), Exceptions::RuntimeError,
128 "MueLu::Hierarchy::GetLevel(): invalid input parameter value: LevelID = " << levelID);
129 return Levels_[levelID];
130}
131
132template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134 return Levels_.size();
135}
136
137template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator>>("A");
140 RCP<const Teuchos::Comm<int>> comm = A->getDomainMap()->getComm();
141
142 int numLevels = GetNumLevels();
143 int numGlobalLevels;
144 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, numLevels, Teuchos::ptr(&numGlobalLevels));
145
146 return numGlobalLevels;
147}
148
149template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
151 double totalNnz = 0, lev0Nnz = 1;
152 for (int i = 0; i < GetNumLevels(); ++i) {
153 TEUCHOS_TEST_FOR_EXCEPTION(!(Levels_[i]->IsAvailable("A")), Exceptions::RuntimeError,
154 "Operator complexity cannot be calculated because A is unavailable on level " << i);
155 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator>>("A");
156 if (A.is_null())
157 break;
158
159 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
160 if (Am.is_null()) {
161 GetOStream(Warnings0) << "Some level operators are not matrices, operator complexity calculation aborted" << std::endl;
162 return 0.0;
163 }
164
165 totalNnz += as<double>(Am->getGlobalNumEntries());
166 if (i == 0)
167 lev0Nnz = totalNnz;
168 }
169 return totalNnz / lev0Nnz;
170}
171
172template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174 double node_sc = 0, global_sc = 0;
175 double a0_nnz = 0;
176 const size_t INVALID = Teuchos::OrdinalTraits<size_t>::invalid();
177 // Get cost of fine matvec
178 if (GetNumLevels() <= 0) return -1.0;
179 if (!Levels_[0]->IsAvailable("A")) return -1.0;
180
181 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator>>("A");
182 if (A.is_null()) return -1.0;
183 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
184 if (Am.is_null()) return -1.0;
185 a0_nnz = as<double>(Am->getGlobalNumEntries());
186
187 // Get smoother complexity at each level
188 for (int i = 0; i < GetNumLevels(); ++i) {
189 size_t level_sc = 0;
190 if (!Levels_[i]->IsAvailable("PreSmoother")) continue;
191 RCP<SmootherBase> S = Levels_[i]->template Get<RCP<SmootherBase>>("PreSmoother");
192 if (S.is_null()) continue;
193 level_sc = S->getNodeSmootherComplexity();
194 if (level_sc == INVALID) {
195 global_sc = -1.0;
196 break;
197 }
198
199 node_sc += as<double>(level_sc);
200 }
201
202 double min_sc = 0.0;
203 RCP<const Teuchos::Comm<int>> comm = A->getDomainMap()->getComm();
204 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, node_sc, Teuchos::ptr(&global_sc));
205 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, node_sc, Teuchos::ptr(&min_sc));
206
207 if (min_sc < 0.0)
208 return -1.0;
209 else
210 return global_sc / a0_nnz;
211}
212
213// Coherence checks todo in Setup() (using an helper function):
214template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
216 TEUCHOS_TEST_FOR_EXCEPTION(level.lib() != lib_, Exceptions::RuntimeError,
217 "MueLu::Hierarchy::CheckLevel(): wrong underlying linear algebra library.");
218 TEUCHOS_TEST_FOR_EXCEPTION(level.GetLevelID() != levelID, Exceptions::RuntimeError,
219 "MueLu::Hierarchy::CheckLevel(): wrong level ID");
220 TEUCHOS_TEST_FOR_EXCEPTION(levelID != 0 && level.GetPreviousLevel() != Levels_[levelID - 1], Exceptions::RuntimeError,
221 "MueLu::Hierarchy::Setup(): wrong level parent");
222}
223
224template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226 for (int i = 0; i < GetNumLevels(); ++i) {
227 RCP<Level> level = Levels_[i];
228 if (level->IsAvailable("A")) {
229 RCP<Operator> Aop = level->Get<RCP<Operator>>("A");
230 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Aop);
231 if (!A.is_null()) {
232 RCP<const Import> xpImporter = A->getCrsGraph()->getImporter();
233 if (!xpImporter.is_null())
234 xpImporter->setDistributorParameters(matvecParams);
235 RCP<const Export> xpExporter = A->getCrsGraph()->getExporter();
236 if (!xpExporter.is_null())
237 xpExporter->setDistributorParameters(matvecParams);
238 }
239 }
240 const std::list<std::string> matrices = {"P", "R", "D0", "NodeMatrix"};
241 for (auto it = matrices.begin(); it != matrices.end(); ++it) {
242 if (level->IsAvailable(*it)) {
243 RCP<Matrix> mat = level->Get<RCP<Matrix>>(*it);
244 if (!mat.is_null()) {
245 RCP<const Import> xpImporter = mat->getCrsGraph()->getImporter();
246 if (!xpImporter.is_null()) {
247 xpImporter->setDistributorParameters(matvecParams);
248 }
249 RCP<const Export> xpExporter = mat->getCrsGraph()->getExporter();
250 if (!xpExporter.is_null())
251 xpExporter->setDistributorParameters(matvecParams);
252 }
253 }
254 }
255 if (level->IsAvailable("Importer")) {
256 RCP<const Import> xpImporter = level->Get<RCP<const Import>>("Importer");
257 if (!xpImporter.is_null())
258 xpImporter->setDistributorParameters(matvecParams);
259 }
260 }
261}
262
263// The function uses three managers: fine, coarse and next coarse
264// We construct the data for the coarse level, and do requests for the next coarse
265template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
267 const RCP<const FactoryManagerBase> fineLevelManager,
268 const RCP<const FactoryManagerBase> coarseLevelManager,
269 const RCP<const FactoryManagerBase> nextLevelManager) {
270 // Use PrintMonitor/TimerMonitor instead of just a FactoryMonitor to print "Level 0" instead of Hierarchy(0)
271 // Print is done after the requests for next coarse level
272
273 TEUCHOS_TEST_FOR_EXCEPTION(LastLevelID() < coarseLevelID, Exceptions::RuntimeError,
274 "MueLu::Hierarchy:Setup(): level " << coarseLevelID << " (specified by coarseLevelID argument) "
275 "must be built before calling this function.");
276
277 Level& level = *Levels_[coarseLevelID];
278
279 bool useStackedTimer = !Teuchos::TimeMonitor::stackedTimerNameIsDefault();
280
281 std::string label = FormattingHelper::getColonLabel(level.getObjectLabel());
282 RCP<TimeMonitor> m1;
283 if (!useStackedTimer)
284 m1 = rcp(new TimeMonitor(*this, label + this->ShortClassName() + ": " + "Setup (total)"));
285 TimeMonitor m2(*this, label + this->ShortClassName() + ": " + "Setup" + " (total, level=" + Teuchos::toString(coarseLevelID) + ")");
286
287 // TODO: pass coarseLevelManager by reference
288 TEUCHOS_TEST_FOR_EXCEPTION(coarseLevelManager == Teuchos::null, Exceptions::RuntimeError,
289 "MueLu::Hierarchy::Setup(): argument coarseLevelManager cannot be null");
290
293
294 if (levelManagers_.size() < coarseLevelID + 1)
295 levelManagers_.resize(coarseLevelID + 1);
296 levelManagers_[coarseLevelID] = coarseLevelManager;
297
298 bool isFinestLevel = (fineLevelManager.is_null());
299 bool isLastLevel = (nextLevelManager.is_null());
300
301 int oldRank = -1;
302 if (isFinestLevel) {
303 RCP<Operator> A = level.Get<RCP<Operator>>("A");
304 RCP<const Map> domainMap = A->getDomainMap();
305 RCP<const Teuchos::Comm<int>> comm = domainMap->getComm();
306
307 // Initialize random seed for reproducibility
309
310 // Record the communicator on the level (used for timers sync)
311 level.SetComm(comm);
312 oldRank = SetProcRankVerbose(comm->getRank());
313
314 // Set the Hierarchy library to match that of the finest level matrix,
315 // even if it was already set
316 lib_ = domainMap->lib();
317 level.setlib(lib_);
318
319 } else {
320 // Permeate library to a coarser level
321 level.setlib(lib_);
322
323 Level& prevLevel = *Levels_[coarseLevelID - 1];
324 oldRank = SetProcRankVerbose(prevLevel.GetComm()->getRank());
325 }
326
327 CheckLevel(level, coarseLevelID);
328
329 // Attach FactoryManager to the fine level
330 RCP<SetFactoryManager> SFMFine;
331 if (!isFinestLevel)
332 SFMFine = rcp(new SetFactoryManager(Levels_[coarseLevelID - 1], fineLevelManager));
333
334 if (isFinestLevel && Levels_[coarseLevelID]->IsAvailable("Coordinates"))
335 ReplaceCoordinateMap(*Levels_[coarseLevelID]);
336
337 // Attach FactoryManager to the coarse level
338 SetFactoryManager SFMCoarse(Levels_[coarseLevelID], coarseLevelManager);
339
340 if (isDumpingEnabled_ && (dumpLevel_ == 0 || dumpLevel_ == -1) && coarseLevelID == 1)
341 DumpCurrentGraph(0);
342
343 RCP<TopSmootherFactory> coarseFact;
344 RCP<TopSmootherFactory> smootherFact = rcp(new TopSmootherFactory(coarseLevelManager, "Smoother"));
345
346 int nextLevelID = coarseLevelID + 1;
347
348 RCP<SetFactoryManager> SFMNext;
349 if (isLastLevel == false) {
350 // We are not at the coarsest level, so there is going to be another level ("next coarse") after this one ("coarse")
351 if (nextLevelID > LastLevelID())
352 AddNewLevel();
353 CheckLevel(*Levels_[nextLevelID], nextLevelID);
354
355 // Attach FactoryManager to the next level (level after coarse)
356 SFMNext = rcp(new SetFactoryManager(Levels_[nextLevelID], nextLevelManager));
357 Levels_[nextLevelID]->Request(TopRAPFactory(coarseLevelManager, nextLevelManager));
358
359 // Do smoother requests here. We don't know whether this is going to be
360 // the coarsest level or not, but we need to DeclareInput before we call
361 // coarseRAPFactory.Build(), otherwise some stuff may be erased after
362 // level releases
363 level.Request(*smootherFact);
364
365 } else {
366 // Similar to smoother above, do the coarse solver request here. We don't
367 // know whether this is going to be the coarsest level or not, but we
368 // need to DeclareInput before we call coarseRAPFactory.Build(),
369 // otherwise some stuff may be erased after level releases. This is
370 // actually evident on ProjectorSmoother. It requires both "A" and
371 // "Nullspace". However, "Nullspace" is erased after all releases, so if
372 // we call the coarse factory request after RAP build we would not have
373 // any data, and cannot get it as we don't have previous managers. The
374 // typical trace looks like this:
375 //
376 // MueLu::Level(0)::GetFactory(Aggregates, 0): No FactoryManager
377 // during request for data " Aggregates" on level 0 by factory TentativePFactory
378 // during request for data " P" on level 1 by factory EminPFactory
379 // during request for data " P" on level 1 by factory TransPFactory
380 // during request for data " R" on level 1 by factory RAPFactory
381 // during request for data " A" on level 1 by factory TentativePFactory
382 // during request for data " Nullspace" on level 2 by factory NullspaceFactory
383 // during request for data " Nullspace" on level 2 by factory NullspacePresmoothFactory
384 // during request for data " Nullspace" on level 2 by factory ProjectorSmoother
385 // during request for data " PreSmoother" on level 2 by factory NoFactory
386 if (coarseFact.is_null())
387 coarseFact = rcp(new TopSmootherFactory(coarseLevelManager, "CoarseSolver"));
388 level.Request(*coarseFact);
389 }
390
391 GetOStream(Runtime0) << std::endl;
392 PrintMonitor m0(*this, "Level " + Teuchos::toString(coarseLevelID), static_cast<MsgType>(Runtime0 | Test));
393
394 // Build coarse level hierarchy
395 RCP<Operator> Ac = Teuchos::null;
396 TopRAPFactory coarseRAPFactory(fineLevelManager, coarseLevelManager);
397
398 if (level.IsAvailable("A")) {
399 Ac = level.Get<RCP<Operator>>("A");
400 } else if (!isFinestLevel) {
401 // We only build here, the release is done later
402 coarseRAPFactory.Build(*level.GetPreviousLevel(), level);
403 }
404
405 bool setLastLevelviaMaxCoarseSize = false;
406 if (level.IsAvailable("A"))
407 Ac = level.Get<RCP<Operator>>("A");
408 RCP<Matrix> Acm = rcp_dynamic_cast<Matrix>(Ac);
409
410 // Record the communicator on the level
411 if (!Ac.is_null())
412 level.SetComm(Ac->getDomainMap()->getComm());
413
414 // Test if we reach the end of the hierarchy
415 bool isOrigLastLevel = isLastLevel;
416 if (isLastLevel) {
417 // Last level as we have achieved the max limit
418 isLastLevel = true;
419
420 } else if (Ac.is_null()) {
421 // Last level for this processor, as it does not belong to the next
422 // subcommunicator. Other processors may continue working on the
423 // hierarchy
424 isLastLevel = true;
425
426 } else {
427 if (!Acm.is_null() && Acm->getGlobalNumRows() <= maxCoarseSize_) {
428 // Last level as the size of the coarse matrix became too small
429 GetOStream(Runtime0) << "Max coarse size (<= " << maxCoarseSize_ << ") achieved" << std::endl;
430 isLastLevel = true;
431 if (Acm->getGlobalNumRows() != 0) setLastLevelviaMaxCoarseSize = true;
432 }
433 }
434
435 if (!Ac.is_null() && !isFinestLevel) {
436 RCP<Operator> A = Levels_[coarseLevelID - 1]->template Get<RCP<Operator>>("A");
437 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
438
439 const double maxCoarse2FineRatio = 0.8;
440 if (!Acm.is_null() && !Am.is_null() && Acm->getGlobalNumRows() > maxCoarse2FineRatio * Am->getGlobalNumRows()) {
441 // We could abort here, but for now we simply notify user.
442 // Couple of additional points:
443 // - if repartitioning is delayed until level K, but the aggregation
444 // procedure stagnates between levels K-1 and K. In this case,
445 // repartitioning could enable faster coarsening once again, but the
446 // hierarchy construction will abort due to the stagnation check.
447 // - if the matrix is small enough, we could move it to one processor.
448 GetOStream(Warnings0) << "Aggregation stagnated. Please check your matrix and/or adjust your configuration file."
449 << "Possible fixes:\n"
450 << " - reduce the maximum number of levels\n"
451 << " - enable repartitioning\n"
452 << " - increase the minimum coarse size." << std::endl;
453 }
454 }
455
456 if (isLastLevel) {
457 if (!isOrigLastLevel) {
458 // We did not expect to finish this early so we did request a smoother.
459 // We need a coarse solver instead. Do the magic.
460 level.Release(*smootherFact);
461 if (coarseFact.is_null())
462 coarseFact = rcp(new TopSmootherFactory(coarseLevelManager, "CoarseSolver"));
463 level.Request(*coarseFact);
464 }
465
466 // Do the actual build, if we have any data.
467 // NOTE: this is not a great check, we may want to call Build() regardless.
468 if (!Ac.is_null())
469 coarseFact->Build(level);
470
471 // Once the dirty deed is done, release stuff. The smoother has already
472 // been released.
473 level.Release(*coarseFact);
474
475 } else {
476 // isLastLevel = false => isOrigLastLevel = false, meaning that we have
477 // requested the smoother. Now we need to build it and to release it.
478 // We don't need to worry about the coarse solver, as we didn't request it.
479 if (!Ac.is_null())
480 smootherFact->Build(level);
481
482 level.Release(*smootherFact);
483 }
484
485 if (isLastLevel == true) {
486 int actualNumLevels = nextLevelID;
487 if (isOrigLastLevel == false) {
488 // Earlier in the function, we constructed the next coarse level, and requested data for the that level,
489 // assuming that we are not at the coarsest level. Now, we changed our mind, so we have to release those.
490 Levels_[nextLevelID]->Release(TopRAPFactory(coarseLevelManager, nextLevelManager));
491
492 // We truncate/resize the hierarchy and possibly remove the last created level if there is
493 // something wrong with it as indicated by its P not being valid. This might happen
494 // if the global number of aggregates turns out to be zero
495
496 if (!setLastLevelviaMaxCoarseSize) {
497 if (Levels_[nextLevelID - 1]->IsAvailable("P")) {
498 if (Levels_[nextLevelID - 1]->template Get<RCP<Matrix>>("P") == Teuchos::null) actualNumLevels = nextLevelID - 1;
499 } else
500 actualNumLevels = nextLevelID - 1;
501 }
502 }
503 if (actualNumLevels == nextLevelID - 1) {
504 // Didn't expect to finish early so we requested smoother but need coarse solver instead.
505 Levels_[nextLevelID - 2]->Release(*smootherFact);
506
507 if (Levels_[nextLevelID - 2]->IsAvailable("PreSmoother")) Levels_[nextLevelID - 2]->RemoveKeepFlag("PreSmoother", NoFactory::get());
508 if (Levels_[nextLevelID - 2]->IsAvailable("PostSmoother")) Levels_[nextLevelID - 2]->RemoveKeepFlag("PostSmoother", NoFactory::get());
509 if (coarseFact.is_null())
510 coarseFact = rcp(new TopSmootherFactory(coarseLevelManager, "CoarseSolver"));
511 Levels_[nextLevelID - 2]->Request(*coarseFact);
512 if (!(Levels_[nextLevelID - 2]->template Get<RCP<Matrix>>("A").is_null()))
513 coarseFact->Build(*(Levels_[nextLevelID - 2]));
514 Levels_[nextLevelID - 2]->Release(*coarseFact);
515 }
516 Levels_.resize(actualNumLevels);
517 levelManagers_.resize(actualNumLevels);
518 }
519
520 // I think this is the proper place for graph so that it shows every dependence
521 if (isDumpingEnabled_ && ((dumpLevel_ > 0 && coarseLevelID == dumpLevel_) || dumpLevel_ == -1))
522 DumpCurrentGraph(coarseLevelID);
523
524 if (!isFinestLevel) {
525 // Release the hierarchy data
526 // We release so late to help blocked solvers, as the smoothers for them need A blocks
527 // which we construct in RAPFactory
528 level.Release(coarseRAPFactory);
529 }
530
531 if (oldRank != -1)
532 SetProcRankVerbose(oldRank);
533
534 return isLastLevel;
535}
536
537template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
539 int numLevels = Levels_.size();
540 TEUCHOS_TEST_FOR_EXCEPTION(levelManagers_.size() != numLevels, Exceptions::RuntimeError,
541 "Hierarchy::SetupRe: " << Levels_.size() << " levels, but " << levelManagers_.size() << " level factory managers");
542
543 const int startLevel = 0;
544 Clear(startLevel);
545
546#ifdef HAVE_MUELU_DEBUG
547 // Reset factories' data used for debugging
548 for (int i = 0; i < numLevels; i++)
549 levelManagers_[i]->ResetDebugData();
550
551#endif
552
553 int levelID;
554 for (levelID = startLevel; levelID < numLevels;) {
555 bool r = Setup(levelID,
556 (levelID != 0 ? levelManagers_[levelID - 1] : Teuchos::null),
557 levelManagers_[levelID],
558 (levelID + 1 != numLevels ? levelManagers_[levelID + 1] : Teuchos::null));
559 levelID++;
560 if (r) break;
561 }
562 // We may construct fewer levels for some reason, make sure we continue
563 // doing that in the future
564 Levels_.resize(levelID);
565 levelManagers_.resize(levelID);
566
567 int sizeOfVecs = sizeOfAllocatedLevelMultiVectors_;
568
569 AllocateLevelMultiVectors(sizeOfVecs, true);
570
571 // since the # of levels, etc. may have changed, force re-determination of description during next call to description()
572 ResetDescription();
573
574 describe(GetOStream(Statistics0), GetVerbLevel());
575
576 CheckForEmptySmoothersAndCoarseSolve();
577}
578
579template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
580void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Setup(const FactoryManagerBase& manager, int startLevel, int numDesiredLevels) {
581 // Use MueLu::BaseClass::description() to avoid printing "{numLevels = 1}" (numLevels is increasing...)
582 PrintMonitor m0(*this, "Setup (" + this->MueLu::BaseClass::description() + ")", Runtime0);
583
584 Clear(startLevel);
585
586 // Check Levels_[startLevel] exists.
587 TEUCHOS_TEST_FOR_EXCEPTION(Levels_.size() <= startLevel, Exceptions::RuntimeError,
588 "MueLu::Hierarchy::Setup(): fine level (" << startLevel << ") does not exist");
589
590 TEUCHOS_TEST_FOR_EXCEPTION(numDesiredLevels <= 0, Exceptions::RuntimeError,
591 "Constructing non-positive (" << numDesiredLevels << ") number of levels does not make sense.");
592
593 // Check for fine level matrix A
594 TEUCHOS_TEST_FOR_EXCEPTION(!Levels_[startLevel]->IsAvailable("A"), Exceptions::RuntimeError,
595 "MueLu::Hierarchy::Setup(): fine level (" << startLevel << ") has no matrix A! "
596 "Set fine level matrix A using Level.Set()");
597
598 RCP<Operator> A = Levels_[startLevel]->template Get<RCP<Operator>>("A");
599 lib_ = A->getDomainMap()->lib();
600
601 if (IsPrint(Statistics2)) {
602 RCP<Matrix> Amat = rcp_dynamic_cast<Matrix>(A);
603
604 if (!Amat.is_null()) {
605 RCP<ParameterList> params = rcp(new ParameterList());
606 params->set("printLoadBalancingInfo", true);
607 params->set("printCommInfo", true);
608
609 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Amat, "A0", params);
610 } else {
611 GetOStream(Warnings1) << "Fine level operator is not a matrix, statistics are not available" << std::endl;
612 }
613 }
614
615 RCP<const FactoryManagerBase> rcpmanager = rcpFromRef(manager);
616
617 const int lastLevel = startLevel + numDesiredLevels - 1;
618 GetOStream(Runtime0) << "Setup loop: startLevel = " << startLevel << ", lastLevel = " << lastLevel
619 << " (stop if numLevels = " << numDesiredLevels << " or Ac.size() < " << maxCoarseSize_ << ")" << std::endl;
620
621 // Setup multigrid levels
622 int iLevel = 0;
623 if (numDesiredLevels == 1) {
624 iLevel = 0;
625 Setup(startLevel, Teuchos::null, rcpmanager, Teuchos::null); // setup finest==coarsest level (first and last managers are Teuchos::null)
626
627 } else {
628 bool bIsLastLevel = Setup(startLevel, Teuchos::null, rcpmanager, rcpmanager); // setup finest level (level 0) (first manager is Teuchos::null)
629 if (bIsLastLevel == false) {
630 for (iLevel = startLevel + 1; iLevel < lastLevel; iLevel++) {
631 bIsLastLevel = Setup(iLevel, rcpmanager, rcpmanager, rcpmanager); // setup intermediate levels
632 if (bIsLastLevel == true)
633 break;
634 }
635 if (bIsLastLevel == false)
636 Setup(lastLevel, rcpmanager, rcpmanager, Teuchos::null); // setup coarsest level (last manager is Teuchos::null)
637 }
638 }
639
640 // TODO: some check like this should be done at the beginning of the routine
641 TEUCHOS_TEST_FOR_EXCEPTION(iLevel != Levels_.size() - 1, Exceptions::RuntimeError,
642 "MueLu::Hierarchy::Setup(): number of level");
643
644 // TODO: this is not exception safe: manager will still hold default
645 // factories if you exit this function with an exception
646 manager.Clean();
647
648 describe(GetOStream(Statistics0), GetVerbLevel());
649
650 CheckForEmptySmoothersAndCoarseSolve();
651}
652
653template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
655 for (LO levelNo = 0; levelNo < GetNumLevels(); ++levelNo) {
656 auto level = Levels_[levelNo];
657 if ((level->IsAvailable("A") && !level->template Get<RCP<Operator>>("A").is_null()) && (!level->IsAvailable("PreSmoother")) && (!level->IsAvailable("PostSmoother"))) {
658 GetOStream(Warnings1) << "No " << (levelNo == as<LO>(Levels_.size()) - 1 ? "coarse grid solver" : "smoother") << " on level " << level->GetLevelID() << std::endl;
659 }
660 }
661}
662
663template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
665 if (startLevel < GetNumLevels())
666 GetOStream(Runtime0) << "Clearing old data (if any)" << std::endl;
667
668 for (int iLevel = startLevel; iLevel < GetNumLevels(); iLevel++)
669 Levels_[iLevel]->Clear();
670}
671
672template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
674 GetOStream(Runtime0) << "Clearing old data (expert)" << std::endl;
675 for (int iLevel = 0; iLevel < GetNumLevels(); iLevel++)
676 Levels_[iLevel]->ExpertClear();
677}
678
679#if defined(HAVE_MUELU_EXPERIMENTAL) && defined(HAVE_MUELU_ADDITIVE_VARIANT)
680template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
681ConvergenceStatus Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Iterate(const MultiVector& B, MultiVector& X, ConvData conv,
682 bool InitialGuessIsZero, LO startLevel) {
683 LO nIts = conv.maxIts_;
684 MagnitudeType tol = conv.tol_;
685
686 std::string prefix = this->ShortClassName() + ": ";
687 std::string levelSuffix = " (level=" + toString(startLevel) + ")";
688 std::string levelSuffix1 = " (level=" + toString(startLevel + 1) + ")";
689
690 using namespace Teuchos;
691 RCP<Time> CompTime = Teuchos::TimeMonitor::getNewCounter(prefix + "Computational Time (total)");
692 RCP<Time> Concurrent = Teuchos::TimeMonitor::getNewCounter(prefix + "Concurrent portion");
693 RCP<Time> ApplyR = Teuchos::TimeMonitor::getNewCounter(prefix + "R: Computational Time");
694 RCP<Time> ApplyPbar = Teuchos::TimeMonitor::getNewCounter(prefix + "Pbar: Computational Time");
695 RCP<Time> CompFine = Teuchos::TimeMonitor::getNewCounter(prefix + "Fine: Computational Time");
696 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix + "Coarse: Computational Time");
697 RCP<Time> ApplySum = Teuchos::TimeMonitor::getNewCounter(prefix + "Sum: Computational Time");
698 RCP<Time> Synchronize_beginning = Teuchos::TimeMonitor::getNewCounter(prefix + "Synchronize_beginning");
699 RCP<Time> Synchronize_center = Teuchos::TimeMonitor::getNewCounter(prefix + "Synchronize_center");
700 RCP<Time> Synchronize_end = Teuchos::TimeMonitor::getNewCounter(prefix + "Synchronize_end");
701
702 RCP<Level> Fine = Levels_[0];
703 RCP<Level> Coarse;
704
705 RCP<Operator> A = Fine->Get<RCP<Operator>>("A");
706 Teuchos::RCP<const Teuchos::Comm<int>> communicator = A->getDomainMap()->getComm();
707
708 // Synchronize_beginning->start();
709 // communicator->barrier();
710 // Synchronize_beginning->stop();
711
712 CompTime->start();
713
714 SC one = STS::one(), zero = STS::zero();
715
716 bool zeroGuess = InitialGuessIsZero;
717
718 // ======= UPFRONT DEFINITION OF COARSE VARIABLES ===========
719
720 // RCP<const Map> origMap;
721 RCP<Operator> P;
722 RCP<Operator> Pbar;
723 RCP<Operator> R;
724 RCP<MultiVector> coarseRhs, coarseX;
725 RCP<Operator> Ac;
726 RCP<SmootherBase> preSmoo_coarse, postSmoo_coarse;
727 bool emptyCoarseSolve = true;
728 RCP<MultiVector> coarseX_prolonged = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(), true);
729
730 RCP<const Import> importer;
731
732 if (Levels_.size() > 1) {
733 Coarse = Levels_[1];
734 if (Coarse->IsAvailable("Importer"))
735 importer = Coarse->Get<RCP<const Import>>("Importer");
736
737 R = Coarse->Get<RCP<Operator>>("R");
738 P = Coarse->Get<RCP<Operator>>("P");
739
740 // if(Coarse->IsAvailable("Pbar"))
741 Pbar = Coarse->Get<RCP<Operator>>("Pbar");
742
743 coarseRhs = MultiVectorFactory::Build(R->getRangeMap(), B.getNumVectors(), true);
744
745 Ac = Coarse->Get<RCP<Operator>>("A");
746
747 ApplyR->start();
748 R->apply(B, *coarseRhs, Teuchos::NO_TRANS, one, zero);
749 // P->apply(B, *coarseRhs, Teuchos::TRANS, one, zero);
750 ApplyR->stop();
751
752 if (doPRrebalance_ || importer.is_null()) {
753 coarseX = MultiVectorFactory::Build(coarseRhs->getMap(), X.getNumVectors(), true);
754
755 } else {
756 RCP<TimeMonitor> ITime = rcp(new TimeMonitor(*this, prefix + "Solve : import (total)", Timings0));
757 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : import" + levelSuffix1, Timings0));
758
759 // Import: range map of R --> domain map of rebalanced Ac (before subcomm replacement)
760 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getTargetMap(), coarseRhs->getNumVectors());
761 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
762 coarseRhs.swap(coarseTmp);
763
764 coarseX = MultiVectorFactory::Build(importer->getTargetMap(), X.getNumVectors(), true);
765 }
766
767 if (Coarse->IsAvailable("PreSmoother"))
768 preSmoo_coarse = Coarse->Get<RCP<SmootherBase>>("PreSmoother");
769 if (Coarse->IsAvailable("PostSmoother"))
770 postSmoo_coarse = Coarse->Get<RCP<SmootherBase>>("PostSmoother");
771 }
772
773 // ==========================================================
774
775 MagnitudeType prevNorm = STS::magnitude(STS::one()), curNorm = STS::magnitude(STS::one());
776 rate_ = 1.0;
777
778 if (Behavior::debug()) {
779 if (A->getDomainMap()->isCompatible(*(X.getMap())) == false) {
780 std::ostringstream ss;
781 ss << "Level " << startLevel << ": level A's domain map is not compatible with X";
782 throw Exceptions::Incompatible(ss.str());
783 }
784
785 if (A->getRangeMap()->isCompatible(*(B.getMap())) == false) {
786 std::ostringstream ss;
787 ss << "Level " << startLevel << ": level A's range map is not compatible with B";
788 throw Exceptions::Incompatible(ss.str());
789 }
790 }
791
792 bool emptyFineSolve = true;
793
794 RCP<MultiVector> fineX;
795 fineX = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(), true);
796
797 // Synchronize_center->start();
798 // communicator->barrier();
799 // Synchronize_center->stop();
800
801 Concurrent->start();
802
803 // NOTE: we need to check using IsAvailable before Get here to avoid building default smoother
804 if (Fine->IsAvailable("PreSmoother")) {
805 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>("PreSmoother");
806 CompFine->start();
807 preSmoo->Apply(*fineX, B, zeroGuess);
808 CompFine->stop();
809 emptyFineSolve = false;
810 }
811 if (Fine->IsAvailable("PostSmoother")) {
812 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>("PostSmoother");
813 CompFine->start();
814 postSmoo->Apply(*fineX, B, zeroGuess);
815 CompFine->stop();
816
817 emptyFineSolve = false;
818 }
819 if (emptyFineSolve == true) {
820 // Fine grid smoother is identity
821 fineX->update(one, B, zero);
822 }
823
824 if (Levels_.size() > 1) {
825 // NOTE: we need to check using IsAvailable before Get here to avoid building default smoother
826 if (Coarse->IsAvailable("PreSmoother")) {
827 CompCoarse->start();
828 preSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
829 CompCoarse->stop();
830 emptyCoarseSolve = false;
831 }
832 if (Coarse->IsAvailable("PostSmoother")) {
833 CompCoarse->start();
834 postSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
835 CompCoarse->stop();
836 emptyCoarseSolve = false;
837 }
838 if (emptyCoarseSolve == true) {
839 // Coarse operator is identity
840 coarseX->update(one, *coarseRhs, zero);
841 }
842 Concurrent->stop();
843 // Synchronize_end->start();
844 // communicator->barrier();
845 // Synchronize_end->stop();
846
847 if (!doPRrebalance_ && !importer.is_null()) {
848 RCP<TimeMonitor> ITime = rcp(new TimeMonitor(*this, prefix + "Solve : export (total)", Timings0));
849 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : export" + levelSuffix1, Timings0));
850
851 // Import: range map of rebalanced Ac (before subcomm replacement) --> domain map of P
852 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getSourceMap(), coarseX->getNumVectors());
853 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
854 coarseX.swap(coarseTmp);
855 }
856
857 ApplyPbar->start();
858 Pbar->apply(*coarseX, *coarseX_prolonged, Teuchos::NO_TRANS, one, zero);
859 ApplyPbar->stop();
860 }
861
862 ApplySum->start();
863 X.update(1.0, *fineX, 1.0, *coarseX_prolonged, 0.0);
864 ApplySum->stop();
865
866 CompTime->stop();
867
868 // communicator->barrier();
869
871}
872#else
873// ---------------------------------------- Iterate -------------------------------------------------------
874template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
876 bool InitialGuessIsZero, LO startLevel) {
877 LO nIts = conv.maxIts_;
878 MagnitudeType tol = conv.tol_;
879
880 // These timers work as follows. "iterateTime" records total time spent in
881 // iterate. "levelTime" records time on a per level basis. The label is
882 // crafted to mimic the per-level messages used in Monitors. Note that a
883 // given level is timed with a TimeMonitor instead of a Monitor or
884 // SubMonitor. This is mainly because I want to time each level
885 // separately, and Monitors/SubMonitors print "(total) xx yy zz" ,
886 // "(sub,total) xx yy zz", respectively, which is subject to
887 // misinterpretation. The per-level TimeMonitors are stopped/started
888 // manually before/after a recursive call to Iterate. A side artifact to
889 // this approach is that the counts for intermediate level timers are twice
890 // the counts for the finest and coarsest levels.
891
892 RCP<Level> Fine = Levels_[startLevel];
893
894 std::string label = FormattingHelper::getColonLabel(Fine->getObjectLabel());
895 std::string prefix = label + this->ShortClassName() + ": ";
896 std::string levelSuffix = " (level=" + toString(startLevel) + ")";
897 std::string levelSuffix1 = " (level=" + toString(startLevel + 1) + ")";
898
899 bool useStackedTimer = !Teuchos::TimeMonitor::stackedTimerNameIsDefault();
900
901 RCP<Monitor> iterateTime;
902 RCP<TimeMonitor> iterateTime1;
903 if (startLevel == 0)
904 iterateTime = rcp(new Monitor(*this, "Solve", label, (nIts == 1) ? None : Runtime0, Timings0));
905 else if (!useStackedTimer)
906 iterateTime1 = rcp(new TimeMonitor(*this, prefix + "Solve (total, level=" + toString(startLevel) + ")", Timings0));
907
908 std::string iterateLevelTimeLabel = prefix + "Solve" + levelSuffix;
909 RCP<TimeMonitor> iterateLevelTime = rcp(new TimeMonitor(*this, iterateLevelTimeLabel, Timings0));
910
911 bool zeroGuess = InitialGuessIsZero;
912
913 RCP<Operator> A = Fine->Get<RCP<Operator>>("A");
914 using namespace Teuchos;
915 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix + "Coarse: Computational Time");
916
917 if (A.is_null()) {
918 // This processor does not have any data for this process on coarser
919 // levels. This can only happen when there are multiple processors and
920 // we use repartitioning.
922 }
923
924 // If we switched the number of vectors, we'd need to reallocate here.
925 // If the number of vectors is unchanged, this is a noop.
926 // NOTE: We need to check against B because the tests in AllocateLevelMultiVectors
927 // will fail on Stokhos Scalar types (due to the so-called 'hidden dimension')
928 const BlockedMultiVector* Bblocked = dynamic_cast<const BlockedMultiVector*>(&B);
929 if (residual_.size() > startLevel &&
930 ((Bblocked && !Bblocked->isSameSize(*residual_[startLevel])) ||
931 (!Bblocked && !residual_[startLevel]->isSameSize(B))))
932 DeleteLevelMultiVectors();
933 AllocateLevelMultiVectors(X.getNumVectors());
934
935 // Print residual information before iterating
936 typedef Teuchos::ScalarTraits<typename STS::magnitudeType> STM;
937 MagnitudeType prevNorm = STM::one();
938 rate_ = 1.0;
939 if (IsCalculationOfResidualRequired(startLevel, conv)) {
940 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, Teuchos::ScalarTraits<LO>::zero(), startLevel, conv, prevNorm);
941 if (convergenceStatus == MueLu::ConvergenceStatus::Converged)
942 return convergenceStatus;
943 }
944
945 SC one = STS::one(), zero = STS::zero();
946 for (LO iteration = 1; iteration <= nIts; iteration++) {
947#ifdef HAVE_MUELU_DEBUG
948#if 0 // TODO fix me
949 if (A->getDomainMap()->isCompatible(*(X.getMap())) == false) {
950 std::ostringstream ss;
951 ss << "Level " << startLevel << ": level A's domain map is not compatible with X";
952 throw Exceptions::Incompatible(ss.str());
953 }
954
955 if (A->getRangeMap()->isCompatible(*(B.getMap())) == false) {
956 std::ostringstream ss;
957 ss << "Level " << startLevel << ": level A's range map is not compatible with B";
958 throw Exceptions::Incompatible(ss.str());
959 }
960#endif
961#endif
962
963 if (startLevel == as<LO>(Levels_.size()) - 1) {
964 // On the coarsest level, we do either smoothing (if defined) or a direct solve.
965 RCP<TimeMonitor> CLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : coarse" + levelSuffix, Timings0));
966
967 bool emptySolve = true;
968
969 // NOTE: we need to check using IsAvailable before Get here to avoid building default smoother
970 if (Fine->IsAvailable("PreSmoother")) {
971 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>("PreSmoother");
972 CompCoarse->start();
973 preSmoo->Apply(X, B, zeroGuess);
974 CompCoarse->stop();
975 zeroGuess = false;
976 emptySolve = false;
977 }
978 if (Fine->IsAvailable("PostSmoother")) {
979 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>("PostSmoother");
980 CompCoarse->start();
981 postSmoo->Apply(X, B, zeroGuess);
982 CompCoarse->stop();
983 emptySolve = false;
984 zeroGuess = false;
985 }
986 if (emptySolve == true) {
987 // Coarse operator is identity
988 X.update(one, B, zero);
989 }
990
991 } else {
992 // On intermediate levels, we do cycles
993 RCP<Level> Coarse = Levels_[startLevel + 1];
994 {
995 // ============== PRESMOOTHING ==============
996 RCP<TimeMonitor> STime;
997 if (!useStackedTimer)
998 STime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing (total)", Timings0));
999 RCP<TimeMonitor> SLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing" + levelSuffix, Timings0));
1000
1001 if (Fine->IsAvailable("PreSmoother")) {
1002 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>("PreSmoother");
1003 preSmoo->Apply(X, B, zeroGuess);
1004 zeroGuess = false;
1005 }
1006 }
1007
1008 RCP<MultiVector> residual;
1009 {
1010 RCP<TimeMonitor> ATime;
1011 if (!useStackedTimer)
1012 ATime = rcp(new TimeMonitor(*this, prefix + "Solve : residual calculation (total)", Timings0));
1013 RCP<TimeMonitor> ALevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : residual calculation" + levelSuffix, Timings0));
1014 if (zeroGuess) {
1015 // If there's a pre-smoother, then zeroGuess is false. If there isn't and people still have zeroGuess set,
1016 // then then X still has who knows what, so we need to zero it before we go to the coarse grid.
1017 X.putScalar(zero);
1018 }
1019
1020 Utilities::Residual(*A, X, B, *residual_[startLevel]);
1021 residual = residual_[startLevel];
1022 }
1023
1024 RCP<Operator> P = Coarse->Get<RCP<Operator>>("P");
1025 if (Coarse->IsAvailable("Pbar"))
1026 P = Coarse->Get<RCP<Operator>>("Pbar");
1027
1028 RCP<MultiVector> coarseRhs, coarseX;
1029 // const bool initializeWithZeros = true;
1030 {
1031 // ============== RESTRICTION ==============
1032 RCP<TimeMonitor> RTime;
1033 if (!useStackedTimer)
1034 RTime = rcp(new TimeMonitor(*this, prefix + "Solve : restriction (total)", Timings0));
1035 RCP<TimeMonitor> RLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : restriction" + levelSuffix, Timings0));
1036 coarseRhs = coarseRhs_[startLevel];
1037
1038 if (implicitTranspose_) {
1039 P->apply(*residual, *coarseRhs, Teuchos::TRANS, one, zero);
1040
1041 } else {
1042 RCP<Operator> R = Coarse->Get<RCP<Operator>>("R");
1043 R->apply(*residual, *coarseRhs, Teuchos::NO_TRANS, one, zero);
1044 }
1045 }
1046
1047 RCP<const Import> importer;
1048 if (Coarse->IsAvailable("Importer"))
1049 importer = Coarse->Get<RCP<const Import>>("Importer");
1050
1051 coarseX = coarseX_[startLevel];
1052 if (!doPRrebalance_ && !importer.is_null()) {
1053 RCP<TimeMonitor> ITime;
1054 if (!useStackedTimer)
1055 ITime = rcp(new TimeMonitor(*this, prefix + "Solve : import (total)", Timings0));
1056 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : import" + levelSuffix1, Timings0));
1057
1058 // Import: range map of R --> domain map of rebalanced Ac (before subcomm replacement)
1059 RCP<MultiVector> coarseTmp = coarseImport_[startLevel];
1060 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
1061 coarseRhs.swap(coarseTmp);
1062 }
1063
1064 RCP<Operator> Ac = Coarse->Get<RCP<Operator>>("A");
1065 if (!Ac.is_null()) {
1066 RCP<const Map> origXMap = coarseX->getMap();
1067 RCP<const Map> origRhsMap = coarseRhs->getMap();
1068
1069 // Replace maps with maps with a subcommunicator
1070 coarseRhs->replaceMap(Ac->getRangeMap());
1071 coarseX->replaceMap(Ac->getDomainMap());
1072
1073 {
1074 iterateLevelTime = Teuchos::null; // stop timing this level
1075
1076 Iterate(*coarseRhs, *coarseX, 1, true, startLevel + 1);
1077 // ^^ zero initial guess
1078 if (Cycle_ == WCYCLE && WCycleStartLevel_ <= startLevel)
1079 Iterate(*coarseRhs, *coarseX, 1, false, startLevel + 1);
1080 // ^^ nonzero initial guess
1081
1082 iterateLevelTime = rcp(new TimeMonitor(*this, iterateLevelTimeLabel)); // restart timing this level
1083 }
1084 coarseX->replaceMap(origXMap);
1085 coarseRhs->replaceMap(origRhsMap);
1086 }
1087
1088 if (!doPRrebalance_ && !importer.is_null()) {
1089 RCP<TimeMonitor> ITime;
1090 if (!useStackedTimer)
1091 ITime = rcp(new TimeMonitor(*this, prefix + "Solve : export (total)", Timings0));
1092 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : export" + levelSuffix1, Timings0));
1093
1094 // Import: range map of rebalanced Ac (before subcomm replacement) --> domain map of P
1095 RCP<MultiVector> coarseTmp = coarseExport_[startLevel];
1096 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
1097 coarseX.swap(coarseTmp);
1098 }
1099
1100 {
1101 // ============== PROLONGATION ==============
1102 RCP<TimeMonitor> PTime;
1103 if (!useStackedTimer)
1104 PTime = rcp(new TimeMonitor(*this, prefix + "Solve : prolongation (total)", Timings0));
1105 RCP<TimeMonitor> PLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : prolongation" + levelSuffix, Timings0));
1106 // Update X += P * coarseX
1107 // Note that due to what may be round-off error accumulation, use of the fused kernel
1108 // P->apply(*coarseX, X, Teuchos::NO_TRANS, one, one);
1109 // can in some cases result in slightly higher iteration counts.
1110 if (fuseProlongationAndUpdate_) {
1111 P->apply(*coarseX, X, Teuchos::NO_TRANS, scalingFactor_, one);
1112 } else {
1113 RCP<MultiVector> correction = correction_[startLevel];
1114 P->apply(*coarseX, *correction, Teuchos::NO_TRANS, one, zero);
1115 X.update(scalingFactor_, *correction, one);
1116 }
1117 }
1118
1119 {
1120 // ============== POSTSMOOTHING ==============
1121 RCP<TimeMonitor> STime;
1122 if (!useStackedTimer)
1123 STime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing (total)", Timings0));
1124 RCP<TimeMonitor> SLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing" + levelSuffix, Timings0));
1125
1126 if (Fine->IsAvailable("PostSmoother")) {
1127 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>("PostSmoother");
1128 postSmoo->Apply(X, B, false);
1129 }
1130 }
1131 }
1132 zeroGuess = false;
1133
1134 if (IsCalculationOfResidualRequired(startLevel, conv)) {
1135 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, iteration, startLevel, conv, prevNorm);
1136 if (convergenceStatus == MueLu::ConvergenceStatus::Converged)
1137 return convergenceStatus;
1138 }
1139 }
1141}
1142#endif
1143
1144template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1145void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(const LO& start, const LO& end, const std::string& suffix) {
1146 LO startLevel = (start != -1 ? start : 0);
1147 LO endLevel = (end != -1 ? end : Levels_.size() - 1);
1148
1149 TEUCHOS_TEST_FOR_EXCEPTION(startLevel > endLevel, Exceptions::RuntimeError,
1150 "MueLu::Hierarchy::Write : startLevel must be <= endLevel");
1151
1152 TEUCHOS_TEST_FOR_EXCEPTION(startLevel < 0 || endLevel >= Levels_.size(), Exceptions::RuntimeError,
1153 "MueLu::Hierarchy::Write bad start or end level");
1154
1155 for (LO i = startLevel; i < endLevel + 1; i++) {
1156 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Levels_[i]->template Get<RCP<Operator>>("A")), P, R;
1157 if (i > 0) {
1158 P = rcp_dynamic_cast<Matrix>(Levels_[i]->template Get<RCP<Operator>>("P"));
1159 if (!implicitTranspose_)
1160 R = rcp_dynamic_cast<Matrix>(Levels_[i]->template Get<RCP<Operator>>("R"));
1161 }
1162
1163 if (!A.is_null()) Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write("A_" + toString(i) + suffix + ".m", *A);
1164 if (!P.is_null()) {
1165 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write("P_" + toString(i) + suffix + ".m", *P);
1166 }
1167 if (!R.is_null()) {
1168 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write("R_" + toString(i) + suffix + ".m", *R);
1169 }
1170 }
1171}
1172
1173template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1174void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Keep(const std::string& ename, const FactoryBase* factory) {
1175 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1176 (*it)->Keep(ename, factory);
1177}
1178
1179template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1181 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1182 (*it)->Delete(ename, factory);
1183}
1184
1185template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1187 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1188 (*it)->AddKeepFlag(ename, factory, keep);
1189}
1190
1191template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1193 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1194 (*it)->RemoveKeepFlag(ename, factory, keep);
1195}
1196
1197template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1199 if (description_.empty()) {
1200 std::ostringstream out;
1201 out << BaseClass::description();
1202 out << "{#levels = " << GetGlobalNumLevels() << ", complexity = " << GetOperatorComplexity() << "}";
1203 description_ = out.str();
1204 }
1205 return description_;
1206}
1207
1208template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1209void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel tVerbLevel) const {
1210 describe(out, toMueLuVerbLevel(tVerbLevel));
1211}
1212
1213template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1214void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream& out, const VerbLevel verbLevel) const {
1215 RCP<Operator> A0 = Levels_[0]->template Get<RCP<Operator>>("A");
1216 RCP<const Teuchos::Comm<int>> comm = A0->getDomainMap()->getComm();
1217
1218 int numLevels = GetNumLevels();
1219 RCP<Operator> Ac = Levels_[numLevels - 1]->template Get<RCP<Operator>>("A");
1220 if (Ac.is_null()) {
1221 // It may happen that we do repartition on the last level, but the matrix
1222 // is small enough to satisfy "max coarse size" requirement. Then, even
1223 // though we have the level, the matrix would be null on all but one processors
1224 numLevels--;
1225 }
1226 int root = comm->getRank();
1227
1228#ifdef HAVE_MPI
1229 int smartData = numLevels * comm->getSize() + comm->getRank(), maxSmartData;
1230 reduceAll(*comm, Teuchos::REDUCE_MAX, smartData, Teuchos::ptr(&maxSmartData));
1231 root = maxSmartData % comm->getSize();
1232#endif
1233
1234 // Compute smoother complexity, if needed
1235 double smoother_comp = -1.0;
1236 if (verbLevel & (Statistics0 | Test))
1237 smoother_comp = GetSmootherComplexity();
1238
1239 std::string outstr;
1240 if (comm->getRank() == root && verbLevel & (Statistics0 | Test)) {
1241 std::vector<Xpetra::global_size_t> nnzPerLevel;
1242 std::vector<Xpetra::global_size_t> rowsPerLevel;
1243 std::vector<int> numProcsPerLevel;
1244 bool someOpsNotMatrices = false;
1245 const Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1246 for (int i = 0; i < numLevels; i++) {
1247 TEUCHOS_TEST_FOR_EXCEPTION(!(Levels_[i]->IsAvailable("A")), Exceptions::RuntimeError,
1248 "Operator A is not available on level " << i);
1249
1250 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator>>("A");
1251 TEUCHOS_TEST_FOR_EXCEPTION(A.is_null(), Exceptions::RuntimeError,
1252 "Operator A on level " << i << " is null.");
1253
1254 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
1255 if (Am.is_null()) {
1256 someOpsNotMatrices = true;
1257 nnzPerLevel.push_back(INVALID);
1258 rowsPerLevel.push_back(A->getDomainMap()->getGlobalNumElements());
1259 numProcsPerLevel.push_back(A->getDomainMap()->getComm()->getSize());
1260 } else {
1261 LO storageblocksize = Am->GetStorageBlockSize();
1262 Xpetra::global_size_t nnz = Am->getGlobalNumEntries() * storageblocksize * storageblocksize;
1263 nnzPerLevel.push_back(nnz);
1264 rowsPerLevel.push_back(Am->getGlobalNumRows() * storageblocksize);
1265 numProcsPerLevel.push_back(Am->getRowMap()->getComm()->getSize());
1266 }
1267 }
1268 if (someOpsNotMatrices)
1269 GetOStream(Warnings0) << "Some level operators are not matrices, statistics calculation are incomplete" << std::endl;
1270
1271 {
1272 std::string label = Levels_[0]->getObjectLabel();
1273 std::ostringstream oss;
1274 oss << std::setfill(' ');
1275 oss << "\n--------------------------------------------------------------------------------\n";
1276 oss << "--- Multigrid Summary " << std::setw(32) << "---\n";
1277 oss << "--------------------------------------------------------------------------------" << std::endl;
1278 if (hierarchyLabel_ != "") oss << "Label = " << hierarchyLabel_ << std::endl;
1279 if (verbLevel & Parameters1)
1280 oss << "Scalar = " << Teuchos::ScalarTraits<Scalar>::name() << std::endl;
1281 oss << "Number of levels = " << numLevels << std::endl;
1282 oss << "Operator complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1283 if (!someOpsNotMatrices)
1284 oss << GetOperatorComplexity() << std::endl;
1285 else
1286 oss << "not available (Some operators in hierarchy are not matrices.)" << std::endl;
1287
1288 if (smoother_comp != -1.0) {
1289 oss << "Smoother complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed)
1290 << smoother_comp << std::endl;
1291 }
1292
1293 switch (Cycle_) {
1294 case VCYCLE:
1295 oss << "Cycle type = V" << std::endl;
1296 break;
1297 case WCYCLE:
1298 oss << "Cycle type = W" << std::endl;
1299 if (WCycleStartLevel_ > 0)
1300 oss << "Cycle start level = " << WCycleStartLevel_ << std::endl;
1301 break;
1302 default:
1303 break;
1304 };
1305 oss << std::endl;
1306
1307 Xpetra::global_size_t tt = rowsPerLevel[0];
1308 int rowspacer = 2;
1309 while (tt != 0) {
1310 tt /= 10;
1311 rowspacer++;
1312 }
1313 for (size_t i = 0; i < nnzPerLevel.size(); ++i) {
1314 tt = nnzPerLevel[i];
1315 if (tt != INVALID)
1316 break;
1317 tt = 100; // This will get used if all levels are operators.
1318 }
1319 int nnzspacer = 2;
1320 while (tt != 0) {
1321 tt /= 10;
1322 nnzspacer++;
1323 }
1324 tt = numProcsPerLevel[0];
1325 int npspacer = 2;
1326 while (tt != 0) {
1327 tt /= 10;
1328 npspacer++;
1329 }
1330 oss << "level " << std::setw(rowspacer) << " rows " << std::setw(nnzspacer) << " nnz "
1331 << " nnz/row" << std::setw(npspacer) << " c ratio"
1332 << " procs" << std::endl;
1333 for (size_t i = 0; i < nnzPerLevel.size(); ++i) {
1334 oss << " " << i << " ";
1335 oss << std::setw(rowspacer) << rowsPerLevel[i];
1336 if (nnzPerLevel[i] != INVALID) {
1337 oss << std::setw(nnzspacer) << nnzPerLevel[i];
1338 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1339 oss << std::setw(9) << as<double>(nnzPerLevel[i]) / rowsPerLevel[i];
1340 } else {
1341 oss << std::setw(nnzspacer) << "Operator";
1342 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1343 oss << std::setw(9) << " ";
1344 }
1345 if (i)
1346 oss << std::setw(9) << as<double>(rowsPerLevel[i - 1]) / rowsPerLevel[i];
1347 else
1348 oss << std::setw(9) << " ";
1349 oss << " " << std::setw(npspacer) << numProcsPerLevel[i] << std::endl;
1350 }
1351 oss << std::endl;
1352 for (int i = 0; i < GetNumLevels(); ++i) {
1353 RCP<SmootherBase> preSmoo, postSmoo;
1354 if (Levels_[i]->IsAvailable("PreSmoother"))
1355 preSmoo = Levels_[i]->template Get<RCP<SmootherBase>>("PreSmoother");
1356 if (Levels_[i]->IsAvailable("PostSmoother"))
1357 postSmoo = Levels_[i]->template Get<RCP<SmootherBase>>("PostSmoother");
1358
1359 if (preSmoo != null && preSmoo == postSmoo)
1360 oss << "Smoother (level " << i << ") both : " << preSmoo->description() << std::endl;
1361 else {
1362 oss << "Smoother (level " << i << ") pre : "
1363 << (preSmoo != null ? preSmoo->description() : "no smoother") << std::endl;
1364 oss << "Smoother (level " << i << ") post : "
1365 << (postSmoo != null ? postSmoo->description() : "no smoother") << std::endl;
1366 }
1367
1368 oss << std::endl;
1369 }
1370
1371 outstr = oss.str();
1372 }
1373 }
1374
1375#ifdef HAVE_MPI
1376 RCP<const Teuchos::MpiComm<int>> mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm);
1377 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
1378
1379 int strLength = outstr.size();
1380 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
1381 if (comm->getRank() != root)
1382 outstr.resize(strLength);
1383 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
1384#endif
1385
1386 out << outstr;
1387}
1388
1389// NOTE: at some point this should be replaced by a friend operator <<
1390template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1391void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(std::ostream& out, const VerbLevel verbLevel) const {
1392 Teuchos::OSTab tab2(out);
1393 for (int i = 0; i < GetNumLevels(); ++i)
1394 Levels_[i]->print(out, verbLevel);
1395}
1396
1397template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1399 isPreconditioner_ = flag;
1400}
1401
1402template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1404 if (GetProcRankVerbose() != 0)
1405 return;
1406#if defined(HAVE_MUELU_BOOST) && defined(HAVE_MUELU_BOOST_FOR_REAL) && defined(BOOST_VERSION) && (BOOST_VERSION >= 104400)
1407
1408 BoostGraph graph;
1409
1410 BoostProperties dp;
1411 dp.property("label", boost::get(boost::vertex_name, graph));
1412 dp.property("id", boost::get(boost::vertex_index, graph));
1413 dp.property("label", boost::get(boost::edge_name, graph));
1414 dp.property("color", boost::get(boost::edge_color, graph));
1415
1416 // create local maps
1417 std::map<const FactoryBase*, BoostVertex> vindices;
1418 typedef std::map<std::pair<BoostVertex, BoostVertex>, std::string> emap;
1419 emap edges;
1420
1421 static int call_id = 0;
1422
1423 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator>>("A");
1424 int rank = A->getDomainMap()->getComm()->getRank();
1425
1426 // printf("[%d] CMS: ----------------------\n",rank);
1427 for (int i = currLevel; i <= currLevel + 1 && i < GetNumLevels(); i++) {
1428 edges.clear();
1429 Levels_[i]->UpdateGraph(vindices, edges, dp, graph);
1430
1431 for (emap::const_iterator eit = edges.begin(); eit != edges.end(); eit++) {
1432 std::pair<BoostEdge, bool> boost_edge = boost::add_edge(eit->first.first, eit->first.second, graph);
1433 // printf("[%d] CMS: Hierarchy, adding edge (%d->%d) %d\n",rank,(int)eit->first.first,(int)eit->first.second,(int)boost_edge.second);
1434 // Because xdot.py views 'Graph' as a keyword
1435 if (eit->second == std::string("Graph"))
1436 boost::put("label", dp, boost_edge.first, std::string("Graph_"));
1437 else
1438 boost::put("label", dp, boost_edge.first, eit->second);
1439 if (i == currLevel)
1440 boost::put("color", dp, boost_edge.first, std::string("red"));
1441 else
1442 boost::put("color", dp, boost_edge.first, std::string("blue"));
1443 }
1444 }
1445
1446 std::ofstream out(dumpFile_.c_str() + std::string("_") + std::to_string(currLevel) + std::string("_") + std::to_string(call_id) + std::string("_") + std::to_string(rank) + std::string(".dot"));
1447 boost::write_graphviz_dp(out, graph, dp, std::string("id"));
1448 out.close();
1449 call_id++;
1450#else
1451 GetOStream(Errors) << "Dependency graph output requires boost and MueLu_ENABLE_Boost_for_real" << std::endl;
1452#endif
1453}
1454
1455// Enforce that coordinate vector's map is consistent with that of A
1456template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1458 RCP<Operator> Ao = level.Get<RCP<Operator>>("A");
1459 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Ao);
1460 if (A.is_null()) {
1461 GetOStream(Runtime1) << "Hierarchy::ReplaceCoordinateMap: operator is not a matrix, skipping..." << std::endl;
1462 return;
1463 }
1464 if (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A) != Teuchos::null) {
1465 GetOStream(Runtime1) << "Hierarchy::ReplaceCoordinateMap: operator is a BlockedCrsMatrix, skipping..." << std::endl;
1466 return;
1467 }
1468
1469 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LO, GO, NO> xdMV;
1470
1471 RCP<xdMV> coords = level.Get<RCP<xdMV>>("Coordinates");
1472
1473 if (A->getRowMap()->isSameAs(*(coords->getMap()))) {
1474 GetOStream(Runtime1) << "Hierarchy::ReplaceCoordinateMap: matrix and coordinates maps are same, skipping..." << std::endl;
1475 return;
1476 }
1477
1478 if (A->IsView("stridedMaps") && rcp_dynamic_cast<const StridedMap>(A->getRowMap("stridedMaps")) != Teuchos::null) {
1479 RCP<const StridedMap> stridedRowMap = rcp_dynamic_cast<const StridedMap>(A->getRowMap("stridedMaps"));
1480
1481 // It is better to through an exceptions if maps may be inconsistent, than to ignore it and experience unfathomable breakdowns
1482 TEUCHOS_TEST_FOR_EXCEPTION(stridedRowMap->getStridedBlockId() != -1 || stridedRowMap->getOffset() != 0,
1483 Exceptions::RuntimeError, "Hierarchy::ReplaceCoordinateMap: nontrivial maps (block id = " << stridedRowMap->getStridedBlockId() << ", offset = " << stridedRowMap->getOffset() << ")");
1484 }
1485
1486 GetOStream(Runtime1) << "Replacing coordinate map" << std::endl;
1487 TEUCHOS_TEST_FOR_EXCEPTION(A->GetFixedBlockSize() % A->GetStorageBlockSize() != 0, Exceptions::RuntimeError, "Hierarchy::ReplaceCoordinateMap: Storage block size does not evenly divide fixed block size");
1488
1489 size_t blkSize = A->GetFixedBlockSize() / A->GetStorageBlockSize();
1490
1491 RCP<const Map> nodeMap = A->getRowMap();
1492 if (blkSize > 1) {
1493 // Create a nodal map, as coordinates have not been expanded to a DOF map yet.
1494 RCP<const Map> dofMap = A->getRowMap();
1495 GO indexBase = dofMap->getIndexBase();
1496 size_t numLocalDOFs = dofMap->getLocalNumElements();
1497 TEUCHOS_TEST_FOR_EXCEPTION(numLocalDOFs % blkSize, Exceptions::RuntimeError,
1498 "Hierarchy::ReplaceCoordinateMap: block size (" << blkSize << ") is incompatible with the number of local dofs in a row map (" << numLocalDOFs);
1499 ArrayView<const GO> GIDs = dofMap->getLocalElementList();
1500
1501 Array<GO> nodeGIDs(numLocalDOFs / blkSize);
1502 for (size_t i = 0; i < numLocalDOFs; i += blkSize)
1503 nodeGIDs[i / blkSize] = (GIDs[i] - indexBase) / blkSize + indexBase;
1504
1505 Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1506 nodeMap = MapFactory::Build(dofMap->lib(), INVALID, nodeGIDs(), indexBase, dofMap->getComm());
1507 } else {
1508 // blkSize == 1
1509 // Check whether the length of vectors fits to the size of A
1510 // If yes, make sure that the maps are matching
1511 // If no, throw a warning but do not touch the Coordinates
1512 if (coords->getLocalLength() != A->getRowMap()->getLocalNumElements()) {
1513 GetOStream(Warnings) << "Coordinate vector does not match row map of matrix A!" << std::endl;
1514 return;
1515 }
1516 }
1517
1518 Array<ArrayView<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>> coordDataView;
1519 std::vector<ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>> coordData;
1520 for (size_t i = 0; i < coords->getNumVectors(); i++) {
1521 coordData.push_back(coords->getData(i));
1522 coordDataView.push_back(coordData[i]());
1523 }
1524
1525 RCP<xdMV> newCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LO, GO, NO>::Build(nodeMap, coordDataView(), coords->getNumVectors());
1526 level.Set("Coordinates", newCoords);
1527}
1528
1529template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1531 int N = Levels_.size();
1532 if (((sizeOfAllocatedLevelMultiVectors_ == numvecs && residual_.size() == N) || numvecs <= 0) && !forceMapCheck) return;
1533
1534 // If, somehow, we changed the number of levels, delete everything first
1535 if (residual_.size() != N) {
1536 DeleteLevelMultiVectors();
1537
1538 residual_.resize(N);
1539 coarseRhs_.resize(N);
1540 coarseX_.resize(N);
1541 coarseImport_.resize(N);
1542 coarseExport_.resize(N);
1543 correction_.resize(N);
1544 }
1545
1546 for (int i = 0; i < N; i++) {
1547 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator>>("A");
1548 if (!A.is_null()) {
1549 // This dance is because we allow A to have a BlockedMap and X/B to have (compatible) non-blocked map
1550 RCP<const BlockedCrsMatrix> A_as_blocked = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(A);
1551 RCP<const Map> Arm = A->getRangeMap();
1552 RCP<const Map> Adm = A->getDomainMap();
1553 if (!A_as_blocked.is_null()) {
1554 Adm = A_as_blocked->getFullDomainMap();
1555 }
1556
1557 if (residual_[i].is_null() || !residual_[i]->getMap()->isSameAs(*Arm))
1558 // This is zero'd by default since it is filled via an operator apply
1559 residual_[i] = MultiVectorFactory::Build(Arm, numvecs, true);
1560 if (correction_[i].is_null() || !correction_[i]->getMap()->isSameAs(*Adm))
1561 correction_[i] = MultiVectorFactory::Build(Adm, numvecs, false);
1562 }
1563
1564 if (i + 1 < N) {
1565 // This is zero'd by default since it is filled via an operator apply
1566 if (implicitTranspose_) {
1567 RCP<Operator> P = Levels_[i + 1]->template Get<RCP<Operator>>("P");
1568 if (!P.is_null()) {
1569 RCP<const Map> map = P->getDomainMap();
1570 if (coarseRhs_[i].is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1571 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs, true);
1572 }
1573 } else {
1574 RCP<Operator> R = Levels_[i + 1]->template Get<RCP<Operator>>("R");
1575 if (!R.is_null()) {
1576 RCP<const Map> map = R->getRangeMap();
1577 if (coarseRhs_[i].is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1578 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs, true);
1579 }
1580 }
1581
1582 RCP<const Import> importer;
1583 if (Levels_[i + 1]->IsAvailable("Importer"))
1584 importer = Levels_[i + 1]->template Get<RCP<const Import>>("Importer");
1585 if (doPRrebalance_ || importer.is_null()) {
1586 RCP<const Map> map = coarseRhs_[i]->getMap();
1587 if (coarseX_[i].is_null() || !coarseX_[i]->getMap()->isSameAs(*map))
1588 coarseX_[i] = MultiVectorFactory::Build(map, numvecs, true);
1589 } else {
1590 RCP<const Map> map;
1591 map = importer->getTargetMap();
1592 if (coarseImport_[i].is_null() || !coarseImport_[i]->getMap()->isSameAs(*map)) {
1593 coarseImport_[i] = MultiVectorFactory::Build(map, numvecs, false);
1594 coarseX_[i] = MultiVectorFactory::Build(map, numvecs, false);
1595 }
1596 map = importer->getSourceMap();
1597 if (coarseExport_[i].is_null() || !coarseExport_[i]->getMap()->isSameAs(*map))
1598 coarseExport_[i] = MultiVectorFactory::Build(map, numvecs, false);
1599 }
1600 }
1601 }
1602 sizeOfAllocatedLevelMultiVectors_ = numvecs;
1603}
1604
1605template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1607 if (sizeOfAllocatedLevelMultiVectors_ == 0) return;
1608 residual_.resize(0);
1609 coarseRhs_.resize(0);
1610 coarseX_.resize(0);
1611 coarseImport_.resize(0);
1612 coarseExport_.resize(0);
1613 correction_.resize(0);
1614 sizeOfAllocatedLevelMultiVectors_ = 0;
1615}
1616
1617template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1619 const LO startLevel, const ConvData& conv) const {
1620 return (startLevel == 0 && !isPreconditioner_ && (IsPrint(Statistics1) || conv.tol_ > 0));
1621}
1622
1623template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1625 const Teuchos::Array<MagnitudeType>& residualNorm, const MagnitudeType convergenceTolerance) const {
1627
1628 if (convergenceTolerance > Teuchos::ScalarTraits<MagnitudeType>::zero()) {
1629 bool passed = true;
1630 for (LO k = 0; k < residualNorm.size(); k++)
1631 if (residualNorm[k] >= convergenceTolerance)
1632 passed = false;
1633
1634 if (passed)
1635 convergenceStatus = ConvergenceStatus::Converged;
1636 else
1637 convergenceStatus = ConvergenceStatus::Unconverged;
1638 }
1639
1640 return convergenceStatus;
1641}
1642
1643template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1645 const LO iteration, const Teuchos::Array<MagnitudeType>& residualNorm) const {
1646 GetOStream(Statistics1) << "iter: "
1647 << std::setiosflags(std::ios::left)
1648 << std::setprecision(3) << std::setw(4) << iteration
1649 << " residual = "
1650 << std::setprecision(10) << residualNorm
1651 << std::endl;
1652}
1653
1654template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1656 const Operator& A, const MultiVector& X, const MultiVector& B, const LO iteration,
1657 const LO startLevel, const ConvData& conv, MagnitudeType& previousResidualNorm) {
1658 Teuchos::Array<MagnitudeType> residualNorm;
1659 residualNorm = Utilities::ResidualNorm(A, X, B, *residual_[startLevel]);
1660
1661 const MagnitudeType currentResidualNorm = residualNorm[0];
1662 rate_ = currentResidualNorm / previousResidualNorm;
1663 previousResidualNorm = currentResidualNorm;
1664
1665 if (IsPrint(Statistics1))
1666 PrintResidualHistory(iteration, residualNorm);
1667
1668 return IsConverged(residualNorm, conv.tol_);
1669}
1670
1671} // namespace MueLu
1672
1673#endif // MUELU_HIERARCHY_DEF_HPP
static bool debug()
Whether MueLu is in debug mode.
virtual std::string description() const
Return a simple one-line description of this object.
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
Base class for factories (e.g., R, P, and A_coarse).
Class that provides default factories within Needs class.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
void AddLevel(const RCP< Level > &level)
Add a level at the end of the hierarchy.
double GetSmootherComplexity() const
void Write(const LO &start=-1, const LO &end=-1, const std::string &suffix="")
Print matrices in the multigrid hierarchy to file.
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
virtual ~Hierarchy()
Destructor.
void CheckLevel(Level &level, int levelID)
Helper function.
std::string description() const
Return a simple one-line description of this object.
void CheckForEmptySmoothersAndCoarseSolve()
void IsPreconditioner(const bool flag)
Array< RCP< Level > > Levels_
Container for Level objects.
bool Setup(int coarseLevelID, const RCP< const FactoryManagerBase > fineLevelManager, const RCP< const FactoryManagerBase > coarseLevelManager, const RCP< const FactoryManagerBase > nextLevelManager=Teuchos::null)
Multi-level setup phase: build a new level of the hierarchy.
STS::magnitudeType MagnitudeType
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
ConvergenceStatus IsConverged(const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const
Decide if the multigrid iteration is converged.
ConvergenceStatus Iterate(const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0)
Apply the multigrid preconditioner.
void DumpCurrentGraph(int level) const
void SetMatvecParams(RCP< ParameterList > matvecParams)
Xpetra::UnderlyingLib lib_
Epetra/Tpetra mode.
void Clear(int startLevel=0)
Clear impermanent data from previous setup.
bool IsCalculationOfResidualRequired(const LO startLevel, const ConvData &conv) const
Decide if the residual needs to be computed.
ConvergenceStatus ComputeResidualAndPrintHistory(const Operator &A, const MultiVector &X, const MultiVector &B, const LO iteration, const LO startLevel, const ConvData &conv, MagnitudeType &previousResidualNorm)
Compute the residual norm and print it depending on the verbosity level.
double GetOperatorComplexity() const
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck=false)
void print(std::ostream &out=std::cout, const VerbLevel verbLevel=(MueLu::Parameters|MueLu::Statistics0)) const
Hierarchy::print is local hierarchy function, thus the statistics can be different from global ones.
void Delete(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Delete(ename, factory) for each level of the Hierarchy.
void Keep(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Keep(ename, factory) for each level of the Hierarchy.
void SetLabel(const std::string &hierarchyLabel)
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
void AddNewLevel()
Add a new level at the end of the hierarchy.
void RemoveKeepFlag(const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All)
Call Level::RemoveKeepFlag for each level of the Hierarchy.
void ReplaceCoordinateMap(Level &level)
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void SetComm(RCP< const Teuchos::Comm< int > > const &comm)
void setlib(Xpetra::UnderlyingLib lib2)
RCP< Level > & GetPreviousLevel()
Previous level.
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
RCP< const Teuchos::Comm< int > > GetComm() const
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())
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.
Xpetra::UnderlyingLib lib()
Timer to be used in non-factories.
static const NoFactory * get()
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
An exception safe way to call the method 'Level::SetFactoryManager()'.
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Statistics2
Print even more statistics.
@ Warnings
Print all warning messages.
@ Statistics1
Print more statistics.
@ Timings0
High level timing information (use Teuchos::TimeMonitor::summarize() to print)
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
@ Warnings1
Additional warnings.
@ Statistics0
Print statistics that do not involve significant additional computation.
@ Parameters1
Print class parameters (more parameters, more verbose)
short KeepType
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
VerbLevel toMueLuVerbLevel(const Teuchos::EVerbosityLevel verbLevel)
Translate Teuchos verbosity level to MueLu verbosity level.
static std::string getColonLabel(const std::string &label)
Helper function for object label.
Data struct for defining stopping criteria of multigrid iteration.