47#ifndef MUELU_HIERARCHY_DEF_HPP
48#define MUELU_HIERARCHY_DEF_HPP
55#include <Xpetra_Matrix.hpp>
56#include <Xpetra_MultiVectorFactory.hpp>
57#include <Xpetra_Operator.hpp>
58#include <Xpetra_IO.hpp>
62#include "MueLu_FactoryManager.hpp"
63#include "MueLu_HierarchyUtils.hpp"
64#include "MueLu_TopRAPFactory.hpp"
65#include "MueLu_TopSmootherFactory.hpp"
68#include "MueLu_PerfUtils.hpp"
69#include "MueLu_PFactory.hpp"
70#include "MueLu_SmootherFactory.hpp"
73#include "Teuchos_TimeMonitor.hpp"
79 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
90 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
94 setObjectLabel(label);
95 Levels_[0]->setObjectLabel(label);
98 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
106 lib_ = A->getDomainMap()->
lib();
108 RCP<Level> Finest = rcp(
new Level);
114 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 setObjectLabel(label);
119 Levels_[0]->setObjectLabel(label);
122 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
126 if (level->GetLevelID() != -1 && (level->GetLevelID() != levelID))
128 " have been added at the end of the hierarchy\n but its ID have been redefined" <<
129 " because last level ID of the hierarchy was " <<
LastLevelID() <<
"." << std::endl;
132 level->SetLevelID(levelID);
135 level->SetPreviousLevel( (levelID == 0) ? Teuchos::null :
Levels_[
LastLevelID() - 1] );
136 level->setObjectLabel(this->getObjectLabel());
139 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
142 newLevel->setlib(
lib_);
146 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
149 "MueLu::Hierarchy::GetLevel(): invalid input parameter value: LevelID = " << levelID);
153 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
158 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
160 RCP<Operator> A =
Levels_[0]->template Get<RCP<Operator> >(
"A");
161 RCP<const Teuchos::Comm<int> > comm = A->getDomainMap()->getComm();
165 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, numLevels, Teuchos::ptr(&numGlobalLevels));
167 return numGlobalLevels;
170 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
172 double totalNnz = 0, lev0Nnz = 1;
175 "Operator complexity cannot be calculated because A is unavailable on level " << i);
176 RCP<Operator> A =
Levels_[i]->template Get<RCP<Operator> >(
"A");
180 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
182 GetOStream(
Warnings0) <<
"Some level operators are not matrices, operator complexity calculation aborted" << std::endl;
186 totalNnz += as<double>(Am->getGlobalNumEntries());
190 return totalNnz / lev0Nnz;
193 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
195 double node_sc = 0, global_sc=0;
197 const size_t INVALID = Teuchos::OrdinalTraits<size_t>::invalid();
200 if (!
Levels_[0]->IsAvailable(
"A"))
return -1.0;
202 RCP<Operator> A =
Levels_[0]->template Get<RCP<Operator> >(
"A");
203 if (A.is_null())
return -1.0;
204 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
205 if(Am.is_null())
return -1.0;
206 a0_nnz = as<double>(Am->getGlobalNumEntries());
211 if(!
Levels_[i]->IsAvailable(
"PreSmoother"))
continue;
212 RCP<SmootherBase> S =
Levels_[i]->template Get<RCP<SmootherBase> >(
"PreSmoother");
213 if (S.is_null())
continue;
214 level_sc = S->getNodeSmootherComplexity();
215 if(level_sc == INVALID) {global_sc=-1.0;
break;}
217 node_sc += as<double>(level_sc);
221 RCP<const Teuchos::Comm<int> > comm =A->getDomainMap()->getComm();
222 Teuchos::reduceAll(*comm,Teuchos::REDUCE_SUM,node_sc,Teuchos::ptr(&global_sc));
223 Teuchos::reduceAll(*comm,Teuchos::REDUCE_MIN,node_sc,Teuchos::ptr(&min_sc));
225 if(min_sc < 0.0)
return -1.0;
226 else return global_sc / a0_nnz;
233 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236 "MueLu::Hierarchy::CheckLevel(): wrong underlying linear algebra library.");
238 "MueLu::Hierarchy::CheckLevel(): wrong level ID");
240 "MueLu::Hierarchy::Setup(): wrong level parent");
243 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
247 if (level->IsAvailable(
"A")) {
248 RCP<Operator> Aop = level->Get<RCP<Operator> >(
"A");
249 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Aop);
251 RCP<const Import> xpImporter = A->getCrsGraph()->getImporter();
252 if (!xpImporter.is_null())
253 xpImporter->setDistributorParameters(matvecParams);
254 RCP<const Export> xpExporter = A->getCrsGraph()->getExporter();
255 if (!xpExporter.is_null())
256 xpExporter->setDistributorParameters(matvecParams);
259 if (level->IsAvailable(
"P")) {
260 RCP<Matrix> P = level->Get<RCP<Matrix> >(
"P");
261 RCP<const Import> xpImporter = P->getCrsGraph()->getImporter();
262 if (!xpImporter.is_null())
263 xpImporter->setDistributorParameters(matvecParams);
264 RCP<const Export> xpExporter = P->getCrsGraph()->getExporter();
265 if (!xpExporter.is_null())
266 xpExporter->setDistributorParameters(matvecParams);
268 if (level->IsAvailable(
"R")) {
269 RCP<Matrix> R = level->Get<RCP<Matrix> >(
"R");
270 RCP<const Import> xpImporter = R->getCrsGraph()->getImporter();
271 if (!xpImporter.is_null())
272 xpImporter->setDistributorParameters(matvecParams);
273 RCP<const Export> xpExporter = R->getCrsGraph()->getExporter();
274 if (!xpExporter.is_null())
275 xpExporter->setDistributorParameters(matvecParams);
277 if (level->IsAvailable(
"Importer")) {
278 RCP<const Import> xpImporter = level->Get< RCP<const Import> >(
"Importer");
279 if (!xpImporter.is_null())
280 xpImporter->setDistributorParameters(matvecParams);
287 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
289 const RCP<const FactoryManagerBase> fineLevelManager,
290 const RCP<const FactoryManagerBase> coarseLevelManager,
291 const RCP<const FactoryManagerBase> nextLevelManager) {
296 "MueLu::Hierarchy:Setup(): level " << coarseLevelID <<
" (specified by coarseLevelID argument) "
297 "must be built before calling this function.");
303 TimeMonitor m2(*
this, label + this->
ShortClassName() +
": " +
"Setup" +
" (total, level=" + Teuchos::toString(coarseLevelID) +
")");
307 "MueLu::Hierarchy::Setup(): argument coarseLevelManager cannot be null");
316 bool isFinestLevel = (fineLevelManager.is_null());
317 bool isLastLevel = (nextLevelManager.is_null());
321 RCP<Operator> A = level.
Get< RCP<Operator> >(
"A");
322 RCP<const Map> domainMap = A->getDomainMap();
323 RCP<const Teuchos::Comm<int> > comm = domainMap->getComm();
334 lib_ = domainMap->lib();
348 RCP<SetFactoryManager> SFMFine;
352 if (isFinestLevel &&
Levels_[coarseLevelID]->IsAvailable(
"Coordinates"))
361 RCP<TopSmootherFactory> coarseFact;
362 RCP<TopSmootherFactory> smootherFact = rcp(
new TopSmootherFactory(coarseLevelManager,
"Smoother"));
364 int nextLevelID = coarseLevelID + 1;
366 RCP<SetFactoryManager> SFMNext;
367 if (isLastLevel ==
false) {
404 if (coarseFact.is_null())
413 RCP<Operator> Ac = Teuchos::null;
414 TopRAPFactory coarseRAPFactory(fineLevelManager, coarseLevelManager);
417 Ac = level.
Get<RCP<Operator> >(
"A");
418 }
else if (!isFinestLevel) {
423 bool setLastLevelviaMaxCoarseSize =
false;
425 Ac = level.
Get<RCP<Operator> >(
"A");
426 RCP<Matrix> Acm = rcp_dynamic_cast<Matrix>(Ac);
430 level.
SetComm(Ac->getDomainMap()->getComm());
433 bool isOrigLastLevel = isLastLevel;
438 }
else if (Ac.is_null()) {
445 if (!Acm.is_null() && Acm->getGlobalNumRows() <=
maxCoarseSize_) {
449 if (Acm->getGlobalNumRows() != 0) setLastLevelviaMaxCoarseSize =
true;
453 if (!Ac.is_null() && !isFinestLevel) {
454 RCP<Operator> A =
Levels_[coarseLevelID-1]->template Get< RCP<Operator> >(
"A");
455 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
457 const double maxCoarse2FineRatio = 0.8;
458 if (!Acm.is_null() && !Am.is_null() && Acm->getGlobalNumRows() > maxCoarse2FineRatio * Am->getGlobalNumRows()) {
466 GetOStream(
Warnings0) <<
"Aggregation stagnated. Please check your matrix and/or adjust your configuration file."
467 <<
"Possible fixes:\n"
468 <<
" - reduce the maximum number of levels\n"
469 <<
" - enable repartitioning\n"
470 <<
" - increase the minimum coarse size." << std::endl;
476 if (!isOrigLastLevel) {
480 if (coarseFact.is_null())
488 coarseFact->Build(level);
499 smootherFact->Build(level);
504 if (isLastLevel ==
true) {
505 int actualNumLevels = nextLevelID;
506 if (isOrigLastLevel ==
false) {
516 if (!setLastLevelviaMaxCoarseSize) {
517 if (
Levels_[nextLevelID-1]->IsAvailable(
"P")) {
518 if (
Levels_[nextLevelID-1]->
template Get<RCP<Matrix> >(
"P") == Teuchos::null) actualNumLevels = nextLevelID-1;
520 else actualNumLevels = nextLevelID-1;
523 if (actualNumLevels == nextLevelID-1) {
525 Levels_[nextLevelID-2]->Release(*smootherFact);
529 if (coarseFact.is_null())
531 Levels_[nextLevelID-2]->Request(*coarseFact);
532 if ( !(
Levels_[nextLevelID-2]->
template Get<RCP<Matrix> >(
"A").is_null() ))
533 coarseFact->Build( *(
Levels_[nextLevelID-2]));
534 Levels_[nextLevelID-2]->Release(*coarseFact);
536 Levels_.resize(actualNumLevels);
543 if (!isFinestLevel) {
547 level.
Release(coarseRAPFactory);
556 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
558 int numLevels =
Levels_.size();
560 "Hierarchy::SetupRe: " <<
Levels_.size() <<
" levels, but " <<
levelManagers_.size() <<
" level factory managers");
562 const int startLevel = 0;
565#ifdef HAVE_MUELU_DEBUG
567 for (
int i = 0; i < numLevels; i++)
573 for (levelID = startLevel; levelID < numLevels;) {
574 bool r =
Setup(levelID,
577 (levelID+1 != numLevels ?
levelManagers_[levelID+1] : Teuchos::null));
596 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
605 "MueLu::Hierarchy::Setup(): fine level (" << startLevel <<
") does not exist");
608 "Constructing non-positive (" << numDesiredLevels <<
") number of levels does not make sense.");
612 "MueLu::Hierarchy::Setup(): fine level (" << startLevel <<
") has no matrix A! "
613 "Set fine level matrix A using Level.Set()");
615 RCP<Operator> A =
Levels_[startLevel]->template Get<RCP<Operator> >(
"A");
616 lib_ = A->getDomainMap()->lib();
619 RCP<Matrix> Amat = rcp_dynamic_cast<Matrix>(A);
621 if (!Amat.is_null()) {
622 RCP<ParameterList> params = rcp(
new ParameterList());
623 params->set(
"printLoadBalancingInfo",
true);
624 params->set(
"printCommInfo",
true);
628 GetOStream(
Warnings1) <<
"Fine level operator is not a matrix, statistics are not available" << std::endl;
632 RCP<const FactoryManagerBase> rcpmanager = rcpFromRef(manager);
634 const int lastLevel = startLevel + numDesiredLevels - 1;
635 GetOStream(
Runtime0) <<
"Setup loop: startLevel = " << startLevel <<
", lastLevel = " << lastLevel
636 <<
" (stop if numLevels = " << numDesiredLevels <<
" or Ac.size() < " <<
maxCoarseSize_ <<
")" << std::endl;
640 if (numDesiredLevels == 1) {
642 Setup(startLevel, Teuchos::null, rcpmanager, Teuchos::null);
645 bool bIsLastLevel =
Setup(startLevel, Teuchos::null, rcpmanager, rcpmanager);
646 if (bIsLastLevel ==
false) {
647 for (iLevel = startLevel + 1; iLevel < lastLevel; iLevel++) {
648 bIsLastLevel =
Setup(iLevel, rcpmanager, rcpmanager, rcpmanager);
649 if (bIsLastLevel ==
true)
652 if (bIsLastLevel ==
false)
653 Setup(lastLevel, rcpmanager, rcpmanager, Teuchos::null);
659 "MueLu::Hierarchy::Setup(): number of level");
668 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
673 for (
int iLevel = startLevel; iLevel <
GetNumLevels(); iLevel++)
677 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
684#if defined(HAVE_MUELU_EXPERIMENTAL) && defined(HAVE_MUELU_ADDITIVE_VARIANT)
685 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
687 bool InitialGuessIsZero, LO startLevel) {
688 LO nIts = conv.maxIts_;
689 MagnitudeType tol = conv.tol_;
691 std::string prefix = this->ShortClassName() +
": ";
692 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
693 std::string levelSuffix1 =
" (level=" +
toString(startLevel+1) +
")";
696 RCP<Time> CompTime = Teuchos::TimeMonitor::getNewCounter(prefix +
"Computational Time (total)");
697 RCP<Time> Concurrent = Teuchos::TimeMonitor::getNewCounter(prefix +
"Concurrent portion");
698 RCP<Time> ApplyR = Teuchos::TimeMonitor::getNewCounter(prefix +
"R: Computational Time");
699 RCP<Time> ApplyPbar = Teuchos::TimeMonitor::getNewCounter(prefix +
"Pbar: Computational Time");
700 RCP<Time> CompFine = Teuchos::TimeMonitor::getNewCounter(prefix +
"Fine: Computational Time");
701 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix +
"Coarse: Computational Time");
702 RCP<Time> ApplySum = Teuchos::TimeMonitor::getNewCounter(prefix +
"Sum: Computational Time");
703 RCP<Time> Synchronize_beginning = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_beginning");
704 RCP<Time> Synchronize_center = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_center");
705 RCP<Time> Synchronize_end = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_end");
707 RCP<Level> Fine = Levels_[0];
710 RCP<Operator> A = Fine->Get< RCP<Operator> >(
"A");
711 Teuchos::RCP< const Teuchos::Comm< int > > communicator = A->getDomainMap()->getComm();
720 SC one = STS::one(), zero = STS::zero();
722 bool zeroGuess = InitialGuessIsZero;
728 RCP< Operator > Pbar;
730 RCP< MultiVector > coarseRhs, coarseX;
732 RCP<SmootherBase> preSmoo_coarse, postSmoo_coarse;
733 bool emptyCoarseSolve =
true;
734 RCP<MultiVector> coarseX_prolonged = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
736 RCP<const Import> importer;
738 if (Levels_.size() > 1) {
740 if (Coarse->IsAvailable(
"Importer"))
741 importer = Coarse->Get< RCP<const Import> >(
"Importer");
743 R = Coarse->Get< RCP<Operator> >(
"R");
744 P = Coarse->Get< RCP<Operator> >(
"P");
748 Pbar = Coarse->Get< RCP<Operator> >(
"Pbar");
750 coarseRhs = MultiVectorFactory::Build(R->getRangeMap(), B.getNumVectors(),
true);
752 Ac = Coarse->Get< RCP< Operator > >(
"A");
755 R->apply(B, *coarseRhs, Teuchos::NO_TRANS, one, zero);
759 if (doPRrebalance_ || importer.is_null()) {
760 coarseX = MultiVectorFactory::Build(coarseRhs->getMap(), X.getNumVectors(),
true);
764 RCP<TimeMonitor> ITime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import (total)" ,
Timings0));
765 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import" + levelSuffix1,
Timings0));
768 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getTargetMap(), coarseRhs->getNumVectors());
769 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
770 coarseRhs.swap(coarseTmp);
772 coarseX = MultiVectorFactory::Build(importer->getTargetMap(), X.getNumVectors(),
true);
775 if (Coarse->IsAvailable(
"PreSmoother"))
776 preSmoo_coarse = Coarse->Get< RCP<SmootherBase> >(
"PreSmoother");
777 if (Coarse->IsAvailable(
"PostSmoother"))
778 postSmoo_coarse = Coarse->Get< RCP<SmootherBase> >(
"PostSmoother");
784 MagnitudeType prevNorm = STS::magnitude(STS::one()), curNorm = STS::magnitude(STS::one());
787 for (LO i = 1; i <= nIts; i++) {
788#ifdef HAVE_MUELU_DEBUG
789 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
790 std::ostringstream ss;
791 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
795 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
796 std::ostringstream ss;
797 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
803 bool emptyFineSolve =
true;
805 RCP< MultiVector > fineX;
806 fineX = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
815 if (Fine->IsAvailable(
"PreSmoother")) {
816 RCP<SmootherBase> preSmoo = Fine->Get< RCP<SmootherBase> >(
"PreSmoother");
818 preSmoo->Apply(*fineX, B, zeroGuess);
820 emptyFineSolve =
false;
822 if (Fine->IsAvailable(
"PostSmoother")) {
823 RCP<SmootherBase> postSmoo = Fine->Get< RCP<SmootherBase> >(
"PostSmoother");
825 postSmoo->Apply(*fineX, B, zeroGuess);
828 emptyFineSolve =
false;
830 if (emptyFineSolve ==
true) {
831 GetOStream(
Warnings1) <<
"No fine grid smoother" << std::endl;
833 fineX->update(one, B, zero);
836 if (Levels_.size() > 1) {
838 if (Coarse->IsAvailable(
"PreSmoother")) {
840 preSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
842 emptyCoarseSolve =
false;
845 if (Coarse->IsAvailable(
"PostSmoother")) {
847 postSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
849 emptyCoarseSolve =
false;
852 if (emptyCoarseSolve ==
true) {
853 GetOStream(
Warnings1) <<
"No coarse grid solver" << std::endl;
855 coarseX->update(one, *coarseRhs, zero);
862 if (!doPRrebalance_ && !importer.is_null()) {
863 RCP<TimeMonitor> ITime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export (total)" ,
Timings0));
864 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export" + levelSuffix1,
Timings0));
867 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getSourceMap(), coarseX->getNumVectors());
868 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
869 coarseX.swap(coarseTmp);
873 Pbar->apply(*coarseX, *coarseX_prolonged, Teuchos::NO_TRANS, one, zero);
878 X.update(1.0, *fineX, 1.0, *coarseX_prolonged, 0.0);
889 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
891 bool InitialGuessIsZero, LO startLevel) {
907 RCP<Level> Fine =
Levels_[startLevel];
911 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
912 std::string levelSuffix1 =
" (level=" +
toString(startLevel+1) +
")";
914 bool useStackedTimer = !Teuchos::TimeMonitor::getStackedTimer().is_null();
916 RCP<Monitor> iterateTime;
917 RCP<TimeMonitor> iterateTime1;
920 else if (!useStackedTimer)
923 std::string iterateLevelTimeLabel = prefix +
"Solve" + levelSuffix;
924 RCP<TimeMonitor> iterateLevelTime = rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel,
Timings0));
926 bool zeroGuess = InitialGuessIsZero;
928 RCP<Operator> A = Fine->Get< RCP<Operator> >(
"A");
930 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix +
"Coarse: Computational Time");
943 const BlockedMultiVector * Bblocked =
dynamic_cast<const BlockedMultiVector*
>(&B);
945 ( ( Bblocked && !Bblocked->isSameSize(*
residual_[startLevel])) ||
946 (!Bblocked && !
residual_[startLevel]->isSameSize(B))))
951 typedef Teuchos::ScalarTraits<typename STS::magnitudeType> STM;
957 return convergenceStatus;
960 SC one = STS::one(), zero = STS::zero();
961 for (LO iteration = 1; iteration <= nIts; iteration++) {
962#ifdef HAVE_MUELU_DEBUG
964 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
965 std::ostringstream ss;
966 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
970 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
971 std::ostringstream ss;
972 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
978 if (startLevel == as<LO>(
Levels_.size())-1) {
980 RCP<TimeMonitor> CLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : coarse" + levelSuffix,
Timings0));
982 bool emptySolve =
true;
985 if (Fine->IsAvailable(
"PreSmoother")) {
986 RCP<SmootherBase> preSmoo = Fine->Get< RCP<SmootherBase> >(
"PreSmoother");
988 preSmoo->Apply(X, B, zeroGuess);
993 if (Fine->IsAvailable(
"PostSmoother")) {
994 RCP<SmootherBase> postSmoo = Fine->Get< RCP<SmootherBase> >(
"PostSmoother");
996 postSmoo->Apply(X, B, zeroGuess);
1001 if (emptySolve ==
true) {
1004 X.update(one, B, zero);
1009 RCP<Level> Coarse =
Levels_[startLevel+1];
1012 RCP<TimeMonitor> STime;
1013 if (!useStackedTimer)
1015 RCP<TimeMonitor> SLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : smoothing" + levelSuffix,
Timings0));
1017 if (Fine->IsAvailable(
"PreSmoother")) {
1018 RCP<SmootherBase> preSmoo = Fine->Get< RCP<SmootherBase> >(
"PreSmoother");
1019 preSmoo->Apply(X, B, zeroGuess);
1024 RCP<MultiVector> residual;
1026 RCP<TimeMonitor> ATime;
1027 if (!useStackedTimer)
1028 ATime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : residual calculation (total)" ,
Timings0));
1029 RCP<TimeMonitor> ALevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : residual calculation" + levelSuffix,
Timings0));
1040 RCP<Operator> P = Coarse->Get< RCP<Operator> >(
"P");
1041 if (Coarse->IsAvailable(
"Pbar"))
1042 P = Coarse->Get< RCP<Operator> >(
"Pbar");
1044 RCP<MultiVector> coarseRhs, coarseX;
1048 RCP<TimeMonitor> RTime;
1049 if (!useStackedTimer)
1051 RCP<TimeMonitor> RLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : restriction" + levelSuffix,
Timings0));
1055 P->apply(*residual, *coarseRhs, Teuchos::TRANS, one, zero);
1058 RCP<Operator> R = Coarse->Get< RCP<Operator> >(
"R");
1059 R->apply(*residual, *coarseRhs, Teuchos::NO_TRANS, one, zero);
1063 RCP<const Import> importer;
1064 if (Coarse->IsAvailable(
"Importer"))
1065 importer = Coarse->Get< RCP<const Import> >(
"Importer");
1069 RCP<TimeMonitor> ITime;
1070 if (!useStackedTimer)
1072 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import" + levelSuffix1,
Timings0));
1076 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
1077 coarseRhs.swap(coarseTmp);
1080 RCP<Operator> Ac = Coarse->Get< RCP<Operator> >(
"A");
1081 if (!Ac.is_null()) {
1082 RCP<const Map> origXMap = coarseX->getMap();
1083 RCP<const Map> origRhsMap = coarseRhs->getMap();
1086 coarseRhs->replaceMap(Ac->getRangeMap());
1087 coarseX ->replaceMap(Ac->getDomainMap());
1090 iterateLevelTime = Teuchos::null;
1092 Iterate(*coarseRhs, *coarseX, 1,
true, startLevel+1);
1095 Iterate(*coarseRhs, *coarseX, 1,
false, startLevel+1);
1098 iterateLevelTime = rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel));
1100 coarseX->replaceMap(origXMap);
1101 coarseRhs->replaceMap(origRhsMap);
1105 RCP<TimeMonitor> ITime;
1106 if (!useStackedTimer)
1108 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export" + levelSuffix1,
Timings0));
1112 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
1113 coarseX.swap(coarseTmp);
1118 RCP<TimeMonitor> PTime;
1119 if (!useStackedTimer)
1121 RCP<TimeMonitor> PLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : prolongation" + levelSuffix,
Timings0));
1129 RCP<MultiVector> correction =
correction_[startLevel];
1130 P->apply(*coarseX, *correction, Teuchos::NO_TRANS, one, zero);
1137 RCP<TimeMonitor> STime;
1138 if (!useStackedTimer)
1140 RCP<TimeMonitor> SLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : smoothing" + levelSuffix,
Timings0));
1142 if (Fine->IsAvailable(
"PostSmoother")) {
1143 RCP<SmootherBase> postSmoo = Fine->Get< RCP<SmootherBase> >(
"PostSmoother");
1144 postSmoo->Apply(X, B,
false);
1154 return convergenceStatus;
1162 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1164 LO startLevel = (start != -1 ? start : 0);
1165 LO endLevel = (end != -1 ? end :
Levels_.size()-1);
1168 "MueLu::Hierarchy::Write : startLevel must be <= endLevel");
1171 "MueLu::Hierarchy::Write bad start or end level");
1173 for (LO i = startLevel; i < endLevel + 1; i++) {
1174 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(
Levels_[i]->
template Get< RCP< Operator> >(
"A")), P, R;
1176 P = rcp_dynamic_cast<Matrix>(
Levels_[i]->
template Get< RCP< Operator> >(
"P"));
1178 R = rcp_dynamic_cast<Matrix>(
Levels_[i]->
template Get< RCP< Operator> >(
"R"));
1181 if (!A.is_null()) Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"A_" +
toString(i) + suffix +
".m", *A);
1183 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"P_" +
toString(i) + suffix +
".m", *P);
1186 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"R_" +
toString(i) + suffix +
".m", *R);
1191 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1193 for (Array<RCP<Level> >::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1194 (*it)->Keep(ename, factory);
1197 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1199 for (Array<RCP<Level> >::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1200 (*it)->Delete(ename, factory);
1203 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1205 for (Array<RCP<Level> >::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1206 (*it)->AddKeepFlag(ename, factory, keep);
1209 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1211 for (Array<RCP<Level> >::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1212 (*it)->RemoveKeepFlag(ename, factory, keep);
1215 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1219 std::ostringstream out;
1227 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1232 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1234 RCP<Operator> A0 =
Levels_[0]->template Get<RCP<Operator> >(
"A");
1235 RCP<const Teuchos::Comm<int> > comm = A0->getDomainMap()->getComm();
1238 RCP<Operator> Ac =
Levels_[numLevels-1]->template Get<RCP<Operator> >(
"A");
1245 int root = comm->getRank();
1248 int smartData = numLevels*comm->getSize() + comm->getRank(), maxSmartData;
1249 reduceAll(*comm, Teuchos::REDUCE_MAX, smartData, Teuchos::ptr(&maxSmartData));
1250 root = maxSmartData % comm->getSize();
1254 double smoother_comp =-1.0;
1260 std::vector<Xpetra::global_size_t> nnzPerLevel;
1261 std::vector<Xpetra::global_size_t> rowsPerLevel;
1262 std::vector<int> numProcsPerLevel;
1263 bool someOpsNotMatrices =
false;
1264 const Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1265 for (
int i = 0; i < numLevels; i++) {
1267 "Operator A is not available on level " << i);
1269 RCP<Operator> A =
Levels_[i]->template Get<RCP<Operator> >(
"A");
1271 "Operator A on level " << i <<
" is null.");
1273 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
1275 someOpsNotMatrices =
true;
1276 nnzPerLevel .push_back(INVALID);
1277 rowsPerLevel .push_back(A->getDomainMap()->getGlobalNumElements());
1278 numProcsPerLevel.push_back(A->getDomainMap()->getComm()->getSize());
1280 LO storageblocksize=Am->GetStorageBlockSize();
1281 Xpetra::global_size_t nnz = Am->getGlobalNumEntries()*storageblocksize*storageblocksize;
1282 nnzPerLevel .push_back(nnz);
1283 rowsPerLevel .push_back(Am->getGlobalNumRows()*storageblocksize);
1284 numProcsPerLevel.push_back(Am->getRowMap()->getComm()->getSize());
1287 if (someOpsNotMatrices)
1288 GetOStream(
Warnings0) <<
"Some level operators are not matrices, statistics calculation are incomplete" << std::endl;
1291 std::string label =
Levels_[0]->getObjectLabel();
1292 std::ostringstream oss;
1293 oss << std::setfill(
' ');
1294 oss <<
"\n--------------------------------------------------------------------------------\n";
1295 oss <<
"--- Multigrid Summary " << std::setw(28) << std::left << label <<
"---\n";
1296 oss <<
"--------------------------------------------------------------------------------" << std::endl;
1298 oss <<
"Scalar = " << Teuchos::ScalarTraits<Scalar>::name() << std::endl;
1299 oss <<
"Number of levels = " << numLevels << std::endl;
1300 oss <<
"Operator complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1301 if (!someOpsNotMatrices)
1304 oss <<
"not available (Some operators in hierarchy are not matrices.)" << std::endl;
1306 if(smoother_comp!=-1.0) {
1307 oss <<
"Smoother complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed)
1308 << smoother_comp << std::endl;
1313 oss <<
"Cycle type = V" << std::endl;
1316 oss <<
"Cycle type = W" << std::endl;
1325 Xpetra::global_size_t tt = rowsPerLevel[0];
1326 int rowspacer = 2;
while (tt != 0) { tt /= 10; rowspacer++; }
1327 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1328 tt = nnzPerLevel[i];
1333 int nnzspacer = 2;
while (tt != 0) { tt /= 10; nnzspacer++; }
1334 tt = numProcsPerLevel[0];
1335 int npspacer = 2;
while (tt != 0) { tt /= 10; npspacer++; }
1336 oss <<
"level " << std::setw(rowspacer) <<
" rows " << std::setw(nnzspacer) <<
" nnz " <<
" nnz/row" << std::setw(npspacer) <<
" c ratio" <<
" procs" << std::endl;
1337 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1338 oss <<
" " << i <<
" ";
1339 oss << std::setw(rowspacer) << rowsPerLevel[i];
1340 if (nnzPerLevel[i] != INVALID) {
1341 oss << std::setw(nnzspacer) << nnzPerLevel[i];
1342 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1343 oss << std::setw(9) << as<double>(nnzPerLevel[i]) / rowsPerLevel[i];
1345 oss << std::setw(nnzspacer) <<
"Operator";
1346 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1347 oss << std::setw(9) <<
" ";
1349 if (i) oss << std::setw(9) << as<double>(rowsPerLevel[i-1])/rowsPerLevel[i];
1350 else oss << std::setw(9) <<
" ";
1351 oss <<
" " << std::setw(npspacer) << numProcsPerLevel[i] << std::endl;
1355 RCP<SmootherBase> preSmoo, postSmoo;
1356 if (
Levels_[i]->IsAvailable(
"PreSmoother"))
1357 preSmoo =
Levels_[i]->template Get< RCP<SmootherBase> >(
"PreSmoother");
1358 if (
Levels_[i]->IsAvailable(
"PostSmoother"))
1359 postSmoo =
Levels_[i]->template Get< RCP<SmootherBase> >(
"PostSmoother");
1361 if (preSmoo != null && preSmoo == postSmoo)
1362 oss <<
"Smoother (level " << i <<
") both : " << preSmoo->description() << std::endl;
1364 oss <<
"Smoother (level " << i <<
") pre : "
1365 << (preSmoo != null ? preSmoo->description() :
"no smoother") << std::endl;
1366 oss <<
"Smoother (level " << i <<
") post : "
1367 << (postSmoo != null ? postSmoo->description() :
"no smoother") << std::endl;
1378 RCP<const Teuchos::MpiComm<int> > mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
1379 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
1381 int strLength = outstr.size();
1382 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
1383 if (comm->getRank() != root)
1384 outstr.resize(strLength);
1385 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
1392 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1394 Teuchos::OSTab tab2(out);
1399 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1404 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1408#if defined(HAVE_MUELU_BOOST) && defined(HAVE_MUELU_BOOST_FOR_REAL) && defined(BOOST_VERSION) && (BOOST_VERSION >= 104400)
1413 dp.property(
"label", boost::get(boost::vertex_name, graph));
1414 dp.property(
"id", boost::get(boost::vertex_index, graph));
1415 dp.property(
"label", boost::get(boost::edge_name, graph));
1416 dp.property(
"color", boost::get(boost::edge_color, graph));
1419 std::map<const FactoryBase*, BoostVertex> vindices;
1420 typedef std::map<std::pair<BoostVertex,BoostVertex>, std::string> emap; emap edges;
1422 static int call_id=0;
1424 RCP<Operator> A =
Levels_[0]->template Get<RCP<Operator> >(
"A");
1425 int rank = A->getDomainMap()->getComm()->getRank();
1428 for (
int i = currLevel; i <= currLevel+1 && i <
GetNumLevels(); i++) {
1430 Levels_[i]->UpdateGraph(vindices, edges, dp, graph);
1432 for (emap::const_iterator eit = edges.begin(); eit != edges.end(); eit++) {
1433 std::pair<BoostEdge, bool> boost_edge = boost::add_edge(eit->first.first, eit->first.second, graph);
1436 if(eit->second==std::string(
"Graph")) boost::put(
"label", dp, boost_edge.first, std::string(
"Graph_"));
1437 else boost::put(
"label", dp, boost_edge.first, eit->second);
1439 boost::put(
"color", dp, boost_edge.first, std::string(
"red"));
1441 boost::put(
"color", dp, boost_edge.first, std::string(
"blue"));
1445 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"));
1446 boost::write_graphviz_dp(out, graph, dp, std::string(
"id"));
1450 GetOStream(
Errors) <<
"Dependency graph output requires boost and MueLu_ENABLE_Boost_for_real" << std::endl;
1455 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1457 RCP<Operator> Ao = level.
Get<RCP<Operator> >(
"A");
1458 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Ao);
1460 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is not a matrix, skipping..." << std::endl;
1463 if(Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A) != Teuchos::null) {
1464 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is a BlockedCrsMatrix, skipping..." << std::endl;
1468 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO> xdMV;
1470 RCP<xdMV> coords = level.
Get<RCP<xdMV> >(
"Coordinates");
1472 if (A->getRowMap()->isSameAs(*(coords->getMap()))) {
1473 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: matrix and coordinates maps are same, skipping..." << std::endl;
1477 if (A->IsView(
"stridedMaps") && rcp_dynamic_cast<const StridedMap>(A->getRowMap(
"stridedMaps")) != Teuchos::null) {
1478 RCP<const StridedMap> stridedRowMap = rcp_dynamic_cast<const StridedMap>(A->getRowMap(
"stridedMaps"));
1481 TEUCHOS_TEST_FOR_EXCEPTION(stridedRowMap->getStridedBlockId() != -1 || stridedRowMap->getOffset() != 0,
1482 Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: nontrivial maps (block id = " << stridedRowMap->getStridedBlockId()
1483 <<
", offset = " << stridedRowMap->getOffset() <<
")");
1487 TEUCHOS_TEST_FOR_EXCEPTION(A->GetFixedBlockSize() % A->GetStorageBlockSize() != 0,
Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: Storage block size does not evenly divide fixed block size");
1489 size_t blkSize = A->GetFixedBlockSize() / A->GetStorageBlockSize();
1491 RCP<const Map> nodeMap = A->getRowMap();
1494 RCP<const Map> dofMap = A->getRowMap();
1495 GO indexBase = dofMap->getIndexBase();
1496 size_t numLocalDOFs = dofMap->getLocalNumElements();
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();
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;
1505 Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1506 nodeMap = MapFactory::Build(dofMap->lib(), INVALID, nodeGIDs(), indexBase, dofMap->getComm());
1512 if(coords->getLocalLength() != A->getRowMap()->getLocalNumElements()) {
1513 GetOStream(
Warnings) <<
"Coordinate vector does not match row map of matrix A!" << std::endl;
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]());
1525 RCP<xdMV> newCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO>::Build(nodeMap, coordDataView(), coords->getNumVectors());
1526 level.
Set(
"Coordinates", newCoords);
1529 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1546 for(
int i=0; i<N; i++) {
1547 RCP<Operator> A =
Levels_[i]->template Get< RCP<Operator> >(
"A");
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();
1559 residual_[i] = MultiVectorFactory::Build(Arm, numvecs,
true);
1561 correction_[i] = MultiVectorFactory::Build(Adm, numvecs,
false);
1567 RCP<Operator> P =
Levels_[i+1]->template Get< RCP<Operator> >(
"P");
1569 RCP<const Map> map = P->getDomainMap();
1571 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1574 RCP<Operator> R =
Levels_[i+1]->template Get< RCP<Operator> >(
"R");
1576 RCP<const Map> map = R->getRangeMap();
1578 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1583 RCP<const Import> importer;
1584 if(
Levels_[i+1]->IsAvailable(
"Importer"))
1585 importer =
Levels_[i+1]->template Get< RCP<const Import> >(
"Importer");
1587 RCP<const Map> map =
coarseRhs_[i]->getMap();
1589 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1592 map = importer->getTargetMap();
1594 coarseImport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1595 coarseX_[i] = MultiVectorFactory::Build(map,numvecs,
false);
1597 map = importer->getSourceMap();
1599 coarseExport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1607template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1620template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1622 const LO startLevel,
const ConvData& conv)
const
1628template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1630 const Teuchos::Array<MagnitudeType>& residualNorm,
const MagnitudeType convergenceTolerance)
const
1634 if (convergenceTolerance > Teuchos::ScalarTraits<MagnitudeType>::zero())
1637 for (LO k = 0; k < residualNorm.size(); k++)
1638 if (residualNorm[k] >= convergenceTolerance)
1647 return convergenceStatus;
1651template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1653 const LO iteration,
const Teuchos::Array<MagnitudeType>& residualNorm)
const
1656 << std::setiosflags(std::ios::left)
1657 << std::setprecision(3) << std::setw(4) << iteration
1659 << std::setprecision(10) << residualNorm
1663template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1665 const Operator& A,
const MultiVector& X,
const MultiVector& B,
const LO iteration,
1668 Teuchos::Array<MagnitudeType> residualNorm;
1672 rate_ = currentResidualNorm / previousResidualNorm;
1673 previousResidualNorm = currentResidualNorm;
virtual std::string ShortClassName() const
Return the class name of the object, without template parameters and without namespace.
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.
Class that provides default factories within Needs class.
virtual void Clean() const
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.
Array< RCP< MultiVector > > residual_
void CheckLevel(Level &level, int levelID)
Helper function.
int WCycleStartLevel_
Level at which to start W-cycle.
std::string description() const
Return a simple one-line description of this object.
std::string description_
cache description to avoid recreating in each call to description() - use ResetDescription() to force...
void IsPreconditioner(const bool flag)
static CycleType GetDefaultCycle()
Array< RCP< Level > > Levels_
Container for Level objects.
static bool GetDefaultPRrebalance()
Array< RCP< MultiVector > > coarseRhs_
Array< RCP< MultiVector > > coarseExport_
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.
void ResetDescription()
force recreation of cached description_ next time description() is called:
STS::magnitudeType MagnitudeType
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
bool isPreconditioner_
Hierarchy may be used in a standalone mode, or as a preconditioner.
ConvergenceStatus IsConverged(const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const
Decide if the multigrid iteration is converged.
void DeleteLevelMultiVectors()
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
int sizeOfAllocatedLevelMultiVectors_
Caching (Multi)Vectors used in Hierarchy::Iterate().
static bool GetDefaultImplicitTranspose()
void SetMatvecParams(RCP< ParameterList > matvecParams)
static Xpetra::global_size_t GetDefaultMaxCoarseSize()
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.
CycleType Cycle_
V- or W-cycle.
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
MagnitudeType rate_
Convergece rate.
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
Array< RCP< MultiVector > > coarseX_
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.
double scalingFactor_
Scaling factor to be applied to coarse grid correction.
void Delete(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Delete(ename, factory) for each level of the Hierarchy.
int GetGlobalNumLevels() const
void Keep(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Keep(ename, factory) for each level of the Hierarchy.
Xpetra::UnderlyingLib lib()
static bool GetDefaultFuseProlongationAndUpdate()
bool doPRViaCopyrebalance_
Array< RCP< MultiVector > > coarseImport_
bool isDumpingEnabled_
Graph dumping.
Xpetra::global_size_t maxCoarseSize_
MueLu::FactoryBase FactoryBase
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
Array< RCP< MultiVector > > correction_
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.
bool fuseProlongationAndUpdate_
Array< RCP< const FactoryManagerBase > > levelManagers_
Level managers used during the Setup.
void ReplaceCoordinateMap(Level &level)
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, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
VerbLevel GetVerbLevel() const
Get the verbosity level.
int SetProcRankVerbose(int procRank) const
Set proc rank used for printing.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
int GetProcRankVerbose() const
Get proc rank used for printing. Do not use this information for any other purpose....
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).
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.
Data struct for defining stopping criteria of multigrid iteration.