46 #ifndef MUELU_IFPACK2SMOOTHER_DEF_HPP
47 #define MUELU_IFPACK2SMOOTHER_DEF_HPP
51 #if defined(HAVE_MUELU_TPETRA) && defined(HAVE_MUELU_IFPACK2)
53 #include <Teuchos_ParameterList.hpp>
55 #include <Tpetra_RowMatrix.hpp>
57 #include <Ifpack2_Chebyshev.hpp>
58 #include <Ifpack2_RILUK.hpp>
59 #include <Ifpack2_Relaxation.hpp>
60 #include <Ifpack2_ILUT.hpp>
61 #include <Ifpack2_BlockRelaxation.hpp>
62 #include <Ifpack2_Factory.hpp>
63 #include <Ifpack2_Parameters.hpp>
65 #include <Xpetra_BlockedCrsMatrix.hpp>
66 #include <Xpetra_CrsMatrix.hpp>
67 #include <Xpetra_CrsMatrixWrap.hpp>
68 #include <Xpetra_Matrix.hpp>
69 #include <Xpetra_MultiVectorFactory.hpp>
70 #include <Xpetra_TpetraMultiVector.hpp>
75 #include "MueLu_Utilities.hpp"
78 #ifdef HAVE_MUELU_INTREPID2
81 #include "Intrepid2_Basis.hpp"
89 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
91 : type_(type), overlap_(overlap)
93 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
94 bool isSupported = Ifpack2::Factory::isSupported<tRowMatrix>(
type_) || (
type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
95 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
96 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
97 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
98 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
99 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
100 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
101 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
102 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
103 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
104 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
105 type_ ==
"LINESMOOTHING_BLOCKRELAXATION" ||
106 type_ ==
"TOPOLOGICAL");
112 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
123 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
126 paramList.setParameters(list);
128 RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
130 if(!precList.is_null() && precList->isParameter(
"partitioner: type") && precList->get<std::string>(
"partitioner: type") ==
"linear" &&
131 !precList->isParameter(
"partitioner: local parts")) {
132 precList->set(
"partitioner: local parts", (
int)A_->getNodeNumRows() / A_->GetFixedBlockSize());
135 prec_->setParameters(*precList);
137 paramList.setParameters(*precList);
140 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
142 this->Input(currentLevel,
"A");
144 if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
145 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
146 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
147 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
148 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
149 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
150 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
151 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
152 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
153 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
154 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
155 type_ ==
"LINESMOOTHING_BLOCKRELAXATION") {
156 this->Input(currentLevel,
"CoarseNumZLayers");
157 this->Input(currentLevel,
"LineDetection_VertLineIds");
159 else if (type_ ==
"BLOCK RELAXATION" ||
160 type_ ==
"BLOCK_RELAXATION" ||
161 type_ ==
"BLOCKRELAXATION" ||
163 type_ ==
"BANDED_RELAXATION" ||
164 type_ ==
"BANDED RELAXATION" ||
165 type_ ==
"BANDEDRELAXATION" ||
167 type_ ==
"TRIDI_RELAXATION" ||
168 type_ ==
"TRIDI RELAXATION" ||
169 type_ ==
"TRIDIRELAXATION" ||
170 type_ ==
"TRIDIAGONAL_RELAXATION" ||
171 type_ ==
"TRIDIAGONAL RELAXATION" ||
172 type_ ==
"TRIDIAGONALRELAXATION")
175 ParameterList precList = this->GetParameterList();
176 if(precList.isParameter(
"partitioner: type") &&
177 precList.get<std::string>(
"partitioner: type") ==
"line") {
178 this->Input(currentLevel,
"Coordinates");
181 else if (type_ ==
"TOPOLOGICAL")
184 this->Input(currentLevel,
"pcoarsen: element to node map");
188 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
191 A_ = Factory::Get< RCP<Matrix> >(currentLevel,
"A");
193 if (type_ ==
"SCHWARZ")
194 SetupSchwarz(currentLevel);
196 else if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
197 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
198 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
199 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
200 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
201 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
202 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
203 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
204 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
205 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
206 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
207 type_ ==
"LINESMOOTHING_BLOCKRELAXATION")
208 SetupLineSmoothing(currentLevel);
210 else if (type_ ==
"BLOCK_RELAXATION" ||
211 type_ ==
"BLOCK RELAXATION" ||
212 type_ ==
"BLOCKRELAXATION" ||
214 type_ ==
"BANDED_RELAXATION" ||
215 type_ ==
"BANDED RELAXATION" ||
216 type_ ==
"BANDEDRELAXATION" ||
218 type_ ==
"TRIDI_RELAXATION" ||
219 type_ ==
"TRIDI RELAXATION" ||
220 type_ ==
"TRIDIRELAXATION" ||
221 type_ ==
"TRIDIAGONAL_RELAXATION" ||
222 type_ ==
"TRIDIAGONAL RELAXATION" ||
223 type_ ==
"TRIDIAGONALRELAXATION")
224 SetupBlockRelaxation(currentLevel);
226 else if (type_ ==
"CHEBYSHEV")
227 SetupChebyshev(currentLevel);
229 else if (type_ ==
"TOPOLOGICAL")
231 #ifdef HAVE_MUELU_INTREPID2
232 SetupTopological(currentLevel);
234 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"'TOPOLOGICAL' smoother choice requires Intrepid2");
239 SetupGeneric(currentLevel);
244 this->GetOStream(
Statistics1) << description() << std::endl;
247 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
249 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
251 bool reusePreconditioner =
false;
252 if (this->IsSetup() ==
true) {
254 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
256 bool isTRowMatrix =
true;
257 RCP<const tRowMatrix> tA;
261 isTRowMatrix =
false;
264 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
265 if (!prec.is_null() && isTRowMatrix) {
266 #ifdef IFPACK2_HAS_PROPER_REUSE
267 prec->resetMatrix(tA);
268 reusePreconditioner =
true;
270 this->GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
274 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
275 "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
279 if (!reusePreconditioner) {
280 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
282 bool isBlockedMatrix =
false;
283 RCP<Matrix> merged2Mat;
285 std::string sublistName =
"subdomain solver parameters";
286 if (paramList.isSublist(sublistName)) {
295 ParameterList& subList = paramList.sublist(sublistName);
297 std::string partName =
"partitioner: type";
298 if (subList.isParameter(partName) && subList.get<std::string>(partName) ==
"user") {
299 isBlockedMatrix =
true;
301 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
303 "Matrix A must be of type BlockedCrsMatrix.");
305 size_t numVels = bA->getMatrix(0,0)->getNodeNumRows();
306 size_t numPres = bA->getMatrix(1,0)->getNodeNumRows();
307 size_t numRows = A_->getNodeNumRows();
309 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
311 size_t numBlocks = 0;
312 for (
size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
313 blockSeeds[rowOfB] = numBlocks++;
315 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
317 "Matrix A must be of type BlockedCrsMatrix.");
319 merged2Mat = bA2->Merge();
323 bool haveBoundary =
false;
324 for (LO i = 0; i < boundaryNodes.size(); i++)
325 if (boundaryNodes[i]) {
329 blockSeeds[i] = numBlocks;
335 subList.set(
"partitioner: map", blockSeeds);
336 subList.set(
"partitioner: local parts", as<int>(numBlocks));
339 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
341 isBlockedMatrix =
true;
342 merged2Mat = bA->Merge();
347 RCP<const tRowMatrix> tA;
360 #ifdef HAVE_MUELU_INTREPID2
361 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
371 if (this->IsSetup() ==
true) {
372 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
373 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
376 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
378 typedef typename Node::device_type::execution_space ES;
380 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO;
382 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
386 const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO> >(currentLevel,
"pcoarsen: element to node map");
388 string basisString = paramList.get<
string>(
"pcoarsen: hi basis");
394 auto basis = MueLuIntrepid::BasisFactory<double,ES>(basisString, degree);
396 string topologyTypeString = paramList.get<
string>(
"smoother: neighborhood type");
398 if (topologyTypeString ==
"node")
400 else if (topologyTypeString ==
"edge")
402 else if (topologyTypeString ==
"face")
404 else if (topologyTypeString ==
"cell")
405 dimension = basis->getBaseCellTopology().getDimension();
407 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
408 vector<vector<LocalOrdinal>> seeds;
413 int myNodeCount = A_->getRowMap()->getNodeNumElements();
414 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
415 int localPartitionNumber = 0;
418 nodeSeeds[seed] = localPartitionNumber++;
421 paramList.remove(
"smoother: neighborhood type");
422 paramList.remove(
"pcoarsen: hi basis");
424 paramList.set(
"partitioner: map", nodeSeeds);
425 paramList.set(
"partitioner: type",
"user");
426 paramList.set(
"partitioner: overlap", 1);
427 paramList.set(
"partitioner: local parts",
int(seeds[dimension].size()));
431 type_ =
"BLOCKRELAXATION";
439 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
441 if (this->IsSetup() ==
true) {
442 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
443 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
446 ParameterList& myparamList =
const_cast<ParameterList&
>(this->GetParameterList());
448 LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,
"CoarseNumZLayers");
449 if (CoarseNumZLayers > 0) {
450 Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get< Teuchos::ArrayRCP<LO> >(currentLevel,
"LineDetection_VertLineIds");
454 for(
size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
455 if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
457 size_t numLocalRows = A_->getNodeNumRows();
460 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
465 int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
466 LO matrixBlockSize = A_->GetFixedBlockSize();
467 if(matrixBlockSize > 1 && actualDofsPerNode > 1)
470 "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
472 else if(matrixBlockSize > 1)
474 actualDofsPerNode = A_->GetFixedBlockSize();
476 myparamList.set(
"partitioner: PDE equations", actualDofsPerNode);
478 if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
479 myparamList.set(
"partitioner: type",
"user");
480 myparamList.set(
"partitioner: map",TVertLineIdSmoo);
481 myparamList.set(
"partitioner: local parts",maxPart+1);
484 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
487 Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
488 for (
size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
489 for (
size_t dof = 0; dof < numDofsPerNode; dof++)
490 partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
491 myparamList.set(
"partitioner: type",
"user");
492 myparamList.set(
"partitioner: map",partitionerMap);
493 myparamList.set(
"partitioner: local parts",maxPart + 1);
496 if (type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
497 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
498 type_ ==
"LINESMOOTHING_BANDEDRELAXATION")
499 type_ =
"BANDEDRELAXATION";
500 else if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
501 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
502 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
503 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
504 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
505 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION")
506 type_ =
"TRIDIAGONALRELAXATION";
508 type_ =
"BLOCKRELAXATION";
511 this->GetOStream(
Runtime0) <<
"Line detection failed: fall back to point-wise relaxation" << std::endl;
512 myparamList.remove(
"partitioner: type",
false);
513 myparamList.remove(
"partitioner: map",
false);
514 myparamList.remove(
"partitioner: local parts",
false);
515 type_ =
"RELAXATION";
526 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
528 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
530 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
536 bool reusePreconditioner =
false;
537 if (this->IsSetup() ==
true) {
539 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
541 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
542 if (!prec.is_null()) {
543 #ifdef IFPACK2_HAS_PROPER_REUSE
544 prec->resetMatrix(tA);
545 reusePreconditioner =
true;
547 this->GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
551 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
552 "reverting to full construction" << std::endl;
556 if (!reusePreconditioner) {
557 ParameterList& myparamList =
const_cast<ParameterList&
>(this->GetParameterList());
559 if(myparamList.isParameter(
"partitioner: type") &&
560 myparamList.get<std::string>(
"partitioner: type") ==
"line") {
561 Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > xCoordinates =
562 Factory::Get<Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > >(currentLevel,
"Coordinates");
563 Teuchos::RCP<Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > coordinates = Teuchos::rcpFromRef(Xpetra::toTpetra<
typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>(*xCoordinates));
565 size_t numDofsPerNode = A_->getNodeNumRows() / xCoordinates->getMap()->getNodeNumElements();
566 myparamList.set(
"partitioner: coordinates", coordinates);
567 myparamList.set(
"partitioner: PDE equations", (
int) numDofsPerNode);
578 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
580 if (this->IsSetup() ==
true) {
581 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): SetupChebyshev() has already been called" << std::endl;
582 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available, reverting to full construction" << std::endl;
585 typedef Teuchos::ScalarTraits<SC> STS;
586 SC negone = -STS::one();
588 SC lambdaMax = negone;
590 std::string maxEigString =
"chebyshev: max eigenvalue";
591 std::string eigRatioString =
"chebyshev: ratio eigenvalue";
593 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
596 if (paramList.isParameter(maxEigString)) {
597 if (paramList.isType<
double>(maxEigString))
598 lambdaMax = paramList.get<
double>(maxEigString);
600 lambdaMax = paramList.get<SC>(maxEigString);
601 this->GetOStream(
Statistics1) << maxEigString <<
" (cached with smoother parameter list) = " << lambdaMax << std::endl;
604 lambdaMax = A_->GetMaxEigenvalueEstimate();
605 if (lambdaMax != negone) {
606 this->GetOStream(
Statistics1) << maxEigString <<
" (cached with matrix) = " << lambdaMax << std::endl;
607 paramList.set(maxEigString, lambdaMax);
612 const SC defaultEigRatio = 20;
614 SC ratio = defaultEigRatio;
615 if (paramList.isParameter(eigRatioString)) {
616 if (paramList.isType<
double>(eigRatioString))
617 ratio = paramList.get<
double>(eigRatioString);
619 ratio = paramList.get<SC>(eigRatioString);
626 RCP<const Matrix> fineA = currentLevel.
GetPreviousLevel()->Get<RCP<Matrix> >(
"A");
627 size_t nRowsFine = fineA->getGlobalNumRows();
628 size_t nRowsCoarse = A_->getGlobalNumRows();
630 SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
631 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
635 this->GetOStream(
Statistics1) << eigRatioString <<
" (computed) = " << ratio << std::endl;
636 paramList.set(eigRatioString, ratio);
652 if (lambdaMax == negone) {
653 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
655 Teuchos::RCP<Ifpack2::Chebyshev<MatrixType> > chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
656 if (chebyPrec != Teuchos::null) {
657 lambdaMax = chebyPrec->getLambdaMaxForApply();
658 A_->SetMaxEigenvalueEstimate(lambdaMax);
659 this->GetOStream(
Statistics1) <<
"chebyshev: max eigenvalue (calculated by Ifpack2)" <<
" = " << lambdaMax << std::endl;
661 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone,
Exceptions::RuntimeError,
"MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
665 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
667 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
668 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
674 bool reusePreconditioner =
false;
675 if (this->IsSetup() ==
true) {
677 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
679 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
680 if (!prec.is_null()) {
681 #ifdef IFPACK2_HAS_PROPER_REUSE
682 prec->resetMatrix(tA);
683 reusePreconditioner =
true;
685 this->GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
689 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
690 "reverting to full construction" << std::endl;
694 if (!reusePreconditioner) {
703 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
716 Teuchos::ParameterList paramList;
717 bool supportInitialGuess =
false;
718 if (type_ ==
"CHEBYSHEV") {
719 paramList.set(
"chebyshev: zero starting solution", InitialGuessIsZero);
720 SetPrecParameters(paramList);
721 supportInitialGuess =
true;
723 }
else if (type_ ==
"RELAXATION") {
724 paramList.set(
"relaxation: zero starting solution", InitialGuessIsZero);
725 SetPrecParameters(paramList);
726 supportInitialGuess =
true;
728 }
else if (type_ ==
"KRYLOV") {
729 paramList.set(
"krylov: zero starting solution", InitialGuessIsZero);
730 SetPrecParameters(paramList);
731 supportInitialGuess =
true;
733 }
else if (type_ ==
"SCHWARZ") {
734 paramList.set(
"schwarz: zero starting solution", InitialGuessIsZero);
739 prec_->setParameters(paramList);
740 supportInitialGuess =
true;
750 if (InitialGuessIsZero || supportInitialGuess) {
753 prec_->apply(tpB, tpX);
755 typedef Teuchos::ScalarTraits<Scalar> TST;
757 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
762 prec_->apply(tpB, tpX);
764 X.update(TST::one(), *Correction, TST::one());
768 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
771 smoother->SetParameterList(this->GetParameterList());
775 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
777 std::ostringstream out;
779 out << prec_->description();
782 out <<
"{type = " << type_ <<
"}";
787 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
792 out0 <<
"Prec. type: " << type_ << std::endl;
795 out0 <<
"Parameter list: " << std::endl;
796 Teuchos::OSTab tab2(out);
797 out << this->GetParameterList();
798 out0 <<
"Overlap: " << overlap_ << std::endl;
802 if (prec_ != Teuchos::null) {
803 Teuchos::OSTab tab2(out);
804 out << *prec_ << std::endl;
807 if (verbLevel &
Debug) {
810 <<
"RCP<prec_>: " << prec_ << std::endl;
814 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
816 typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
818 RCP<Ifpack2::Relaxation<MatrixType> > pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType> >(prec_);
819 if(!pr.is_null())
return pr->getNodeSmootherComplexity();
821 RCP<Ifpack2::Chebyshev<MatrixType> > pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
822 if(!pc.is_null())
return pc->getNodeSmootherComplexity();
824 RCP<Ifpack2::BlockRelaxation<MatrixType> > pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType> >(prec_);
825 if(!pb.is_null())
return pb->getNodeSmootherComplexity();
827 RCP<Ifpack2::ILUT<MatrixType> > pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType> >(prec_);
828 if(!pi.is_null())
return pi->getNodeSmootherComplexity();
830 RCP<Ifpack2::RILUK<MatrixType> > pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType> >(prec_);
831 if(!pk.is_null())
return pk->getNodeSmootherComplexity();
834 return Teuchos::OrdinalTraits<size_t>::invalid();
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
MueLu::DefaultLocalOrdinal LocalOrdinal
static Teuchos::RCP< Preconditioner< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > create(const std::string &precType, const Teuchos::RCP< const MatrixType > &matrix)
virtual std::string description() const
Return a simple one-line description of this object.
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that encapsulates Ifpack2 smoothers.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void SetupSchwarz(Level ¤tLevel)
void Setup(Level ¤tLevel)
Set up the smoother.
void SetupBlockRelaxation(Level ¤tLevel)
void SetupChebyshev(Level ¤tLevel)
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
RCP< SmootherPrototype > Copy() const
std::string description() const
Return a simple one-line description of this object.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
void SetupTopological(Level ¤tLevel)
std::string type_
ifpack2-specific key phrase that denote smoother type
friend class Ifpack2Smoother
Constructor.
void DeclareInput(Level ¤tLevel) const
Input.
void SetupLineSmoothing(Level ¤tLevel)
void SetupGeneric(Level ¤tLevel)
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
RCP< Level > & GetPreviousLevel()
Previous level.
virtual void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
void declareConstructionOutcome(bool fail, std::string msg)
bool IsSetup() const
Get the state of a smoother prototype.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::magnitude(0.), const bool count_twos_as_dirichlet=false)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Debug
Print additional debugging information.
@ Statistics1
Print more statistics.
@ External
Print external lib objects.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
@ Parameters0
Print class parameters.
@ Parameters1
Print class parameters (more parameters, more verbose)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.