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_Relaxation.hpp>
59 #include <Ifpack2_ILUT.hpp>
60 #include <Ifpack2_BlockRelaxation.hpp>
61 #include <Ifpack2_Factory.hpp>
62 #include <Ifpack2_Parameters.hpp>
64 #include <Xpetra_BlockedCrsMatrix.hpp>
65 #include <Xpetra_CrsMatrix.hpp>
66 #include <Xpetra_CrsMatrixWrap.hpp>
67 #include <Xpetra_Matrix.hpp>
68 #include <Xpetra_MultiVectorFactory.hpp>
69 #include <Xpetra_TpetraMultiVector.hpp>
74 #include "MueLu_Utilities.hpp"
77 #ifdef HAVE_MUELU_INTREPID2
80 #include "Intrepid2_Basis.hpp"
88 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
90 : type_(type), overlap_(overlap)
95 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
106 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
108 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
109 paramList.setParameters(list);
111 RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
113 prec_->setParameters(*precList);
115 paramList.setParameters(*precList);
118 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
120 this->Input(currentLevel,
"A");
122 if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
123 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
124 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
125 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
126 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
127 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
128 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
129 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
130 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
131 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
132 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
133 type_ ==
"LINESMOOTHING_BLOCKRELAXATION") {
134 this->Input(currentLevel,
"CoarseNumZLayers");
135 this->Input(currentLevel,
"LineDetection_VertLineIds");
137 else if (type_ ==
"BLOCK RELAXATION" ||
138 type_ ==
"BLOCK_RELAXATION" ||
139 type_ ==
"BLOCKRELAXATION" ||
141 type_ ==
"BANDED_RELAXATION" ||
142 type_ ==
"BANDED RELAXATION" ||
143 type_ ==
"BANDEDRELAXATION" ||
145 type_ ==
"TRIDI_RELAXATION" ||
146 type_ ==
"TRIDI RELAXATION" ||
147 type_ ==
"TRIDIRELAXATION" ||
148 type_ ==
"TRIDIAGONAL_RELAXATION" ||
149 type_ ==
"TRIDIAGONAL RELAXATION" ||
150 type_ ==
"TRIDIAGONALRELAXATION")
153 ParameterList precList = this->GetParameterList();
154 if(precList.isParameter(
"partitioner: type") &&
155 precList.get<std::string>(
"partitioner: type") ==
"line") {
156 this->Input(currentLevel,
"Coordinates");
159 else if (type_ ==
"TOPOLOGICAL")
162 this->Input(currentLevel,
"pcoarsen: element to node map");
166 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
170 A_ = Factory::Get< RCP<Matrix> >(currentLevel,
"A");
172 if (type_ ==
"SCHWARZ")
173 SetupSchwarz(currentLevel);
175 else if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
176 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
177 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
178 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
179 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
180 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
181 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
182 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
183 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
184 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
185 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
186 type_ ==
"LINESMOOTHING_BLOCKRELAXATION")
187 SetupLineSmoothing(currentLevel);
189 else if (type_ ==
"BLOCK_RELAXATION" ||
190 type_ ==
"BLOCK RELAXATION" ||
191 type_ ==
"BLOCKRELAXATION" ||
193 type_ ==
"BANDED_RELAXATION" ||
194 type_ ==
"BANDED RELAXATION" ||
195 type_ ==
"BANDEDRELAXATION" ||
197 type_ ==
"TRIDI_RELAXATION" ||
198 type_ ==
"TRIDI RELAXATION" ||
199 type_ ==
"TRIDIRELAXATION" ||
200 type_ ==
"TRIDIAGONAL_RELAXATION" ||
201 type_ ==
"TRIDIAGONAL RELAXATION" ||
202 type_ ==
"TRIDIAGONALRELAXATION")
203 SetupBlockRelaxation(currentLevel);
205 else if (type_ ==
"CHEBYSHEV")
206 SetupChebyshev(currentLevel);
208 else if (type_ ==
"TOPOLOGICAL")
210 #ifdef HAVE_MUELU_INTREPID2
211 SetupTopological(currentLevel);
213 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"'TOPOLOGICAL' smoother choice requires Intrepid2");
218 SetupGeneric(currentLevel);
223 this->GetOStream(
Statistics1) << description() << std::endl;
226 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
228 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
230 bool reusePreconditioner =
false;
231 if (this->IsSetup() ==
true) {
233 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
235 bool isTRowMatrix =
true;
236 RCP<const tRowMatrix> tA;
240 isTRowMatrix =
false;
243 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
244 if (!prec.is_null() && isTRowMatrix) {
245 #ifdef IFPACK2_HAS_PROPER_REUSE
246 prec->resetMatrix(tA);
247 reusePreconditioner =
true;
249 this->GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
253 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
254 "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
258 if (!reusePreconditioner) {
259 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
261 bool isBlockedMatrix =
false;
262 RCP<Matrix> merged2Mat;
264 std::string sublistName =
"subdomain solver parameters";
265 if (paramList.isSublist(sublistName)) {
274 ParameterList& subList = paramList.sublist(sublistName);
276 std::string partName =
"partitioner: type";
277 if (subList.isParameter(partName) && subList.get<std::string>(partName) ==
"user") {
278 isBlockedMatrix =
true;
280 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
282 "Matrix A must be of type BlockedCrsMatrix.");
284 size_t numVels = bA->getMatrix(0,0)->getNodeNumRows();
285 size_t numPres = bA->getMatrix(1,0)->getNodeNumRows();
286 size_t numRows = A_->getNodeNumRows();
288 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
290 size_t numBlocks = 0;
291 for (
size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
292 blockSeeds[rowOfB] = numBlocks++;
294 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
296 "Matrix A must be of type BlockedCrsMatrix.");
298 merged2Mat = bA2->Merge();
302 bool haveBoundary =
false;
303 for (LO i = 0; i < boundaryNodes.size(); i++)
304 if (boundaryNodes[i]) {
308 blockSeeds[i] = numBlocks;
314 subList.set(
"partitioner: map", blockSeeds);
315 subList.set(
"partitioner: local parts", as<int>(numBlocks));
318 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
320 isBlockedMatrix =
true;
321 merged2Mat = bA->Merge();
326 RCP<const tRowMatrix> tA;
339 #ifdef HAVE_MUELU_INTREPID2
340 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
350 if (this->IsSetup() ==
true) {
351 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
352 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
355 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
357 typedef typename Node::device_type::execution_space ES;
359 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO;
361 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
365 const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO> >(currentLevel,
"pcoarsen: element to node map");
367 string basisString = paramList.get<
string>(
"pcoarsen: hi basis");
373 auto basis = MueLuIntrepid::BasisFactory<double,ES>(basisString, degree);
375 string topologyTypeString = paramList.get<
string>(
"smoother: neighborhood type");
377 if (topologyTypeString ==
"node")
379 else if (topologyTypeString ==
"edge")
381 else if (topologyTypeString ==
"face")
383 else if (topologyTypeString ==
"cell")
384 dimension = basis->getBaseCellTopology().getDimension();
386 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
387 vector<vector<LocalOrdinal>> seeds;
392 int myNodeCount = A_->getRowMap()->getNodeNumElements();
393 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
394 int localPartitionNumber = 0;
397 nodeSeeds[seed] = localPartitionNumber++;
400 paramList.remove(
"smoother: neighborhood type");
401 paramList.remove(
"pcoarsen: hi basis");
403 paramList.set(
"partitioner: map", nodeSeeds);
404 paramList.set(
"partitioner: type",
"user");
405 paramList.set(
"partitioner: overlap", 1);
406 paramList.set(
"partitioner: local parts",
int(seeds[dimension].size()));
410 type_ =
"BLOCKRELAXATION";
418 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
420 if (this->IsSetup() ==
true) {
421 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
422 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
425 ParameterList& myparamList =
const_cast<ParameterList&
>(this->GetParameterList());
427 LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,
"CoarseNumZLayers");
428 if (CoarseNumZLayers > 0) {
429 Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get< Teuchos::ArrayRCP<LO> >(currentLevel,
"LineDetection_VertLineIds");
433 for(
size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
434 if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
436 size_t numLocalRows = A_->getNodeNumRows();
439 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
444 int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
445 LO matrixBlockSize = A_->GetFixedBlockSize();
446 if(matrixBlockSize > 1 && actualDofsPerNode > 1)
449 "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
451 else if(matrixBlockSize > 1)
453 actualDofsPerNode = A_->GetFixedBlockSize();
455 myparamList.set(
"partitioner: PDE equations", actualDofsPerNode);
457 if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
458 myparamList.set(
"partitioner: type",
"user");
459 myparamList.set(
"partitioner: map",TVertLineIdSmoo);
460 myparamList.set(
"partitioner: local parts",maxPart+1);
463 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
466 Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
467 for (
size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
468 for (
size_t dof = 0; dof < numDofsPerNode; dof++)
469 partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
470 myparamList.set(
"partitioner: type",
"user");
471 myparamList.set(
"partitioner: map",partitionerMap);
472 myparamList.set(
"partitioner: local parts",maxPart + 1);
475 if (type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
476 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
477 type_ ==
"LINESMOOTHING_BANDEDRELAXATION")
478 type_ =
"BANDEDRELAXATION";
479 else if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
480 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
481 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
482 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
483 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
484 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION")
485 type_ =
"TRIDIAGONALRELAXATION";
487 type_ =
"BLOCKRELAXATION";
490 this->GetOStream(
Runtime0) <<
"Line detection failed: fall back to point-wise relaxation" << std::endl;
491 myparamList.remove(
"partitioner: type",
false);
492 myparamList.remove(
"partitioner: map",
false);
493 myparamList.remove(
"partitioner: local parts",
false);
494 type_ =
"RELAXATION";
505 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
507 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
509 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
515 bool reusePreconditioner =
false;
516 if (this->IsSetup() ==
true) {
518 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
520 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
521 if (!prec.is_null()) {
522 #ifdef IFPACK2_HAS_PROPER_REUSE
523 prec->resetMatrix(tA);
524 reusePreconditioner =
true;
526 this->GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
530 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
531 "reverting to full construction" << std::endl;
535 if (!reusePreconditioner) {
536 ParameterList& myparamList =
const_cast<ParameterList&
>(this->GetParameterList());
538 if(myparamList.isParameter(
"partitioner: type") &&
539 myparamList.get<std::string>(
"partitioner: type") ==
"line") {
540 Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > xCoordinates =
541 Factory::Get<Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > >(currentLevel,
"Coordinates");
542 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));
544 size_t numDofsPerNode = A_->getNodeNumRows() / xCoordinates->getMap()->getNodeNumElements();
545 myparamList.set(
"partitioner: coordinates", coordinates);
546 myparamList.set(
"partitioner: PDE equations", (
int) numDofsPerNode);
557 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
559 if (this->IsSetup() ==
true) {
560 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): SetupChebyshev() has already been called" << std::endl;
561 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available, reverting to full construction" << std::endl;
564 typedef Teuchos::ScalarTraits<SC> STS;
565 SC negone = -STS::one();
567 SC lambdaMax = negone;
569 std::string maxEigString =
"chebyshev: max eigenvalue";
570 std::string eigRatioString =
"chebyshev: ratio eigenvalue";
572 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
575 if (paramList.isParameter(maxEigString)) {
576 if (paramList.isType<
double>(maxEigString))
577 lambdaMax = paramList.get<
double>(maxEigString);
579 lambdaMax = paramList.get<SC>(maxEigString);
580 this->GetOStream(
Statistics1) << maxEigString <<
" (cached with smoother parameter list) = " << lambdaMax << std::endl;
583 lambdaMax = A_->GetMaxEigenvalueEstimate();
584 if (lambdaMax != negone) {
585 this->GetOStream(
Statistics1) << maxEigString <<
" (cached with matrix) = " << lambdaMax << std::endl;
586 paramList.set(maxEigString, lambdaMax);
591 const SC defaultEigRatio = 20;
593 SC ratio = defaultEigRatio;
594 if (paramList.isParameter(eigRatioString)) {
595 if (paramList.isType<
double>(eigRatioString))
596 ratio = paramList.get<
double>(eigRatioString);
598 ratio = paramList.get<SC>(eigRatioString);
605 RCP<const Matrix> fineA = currentLevel.
GetPreviousLevel()->Get<RCP<Matrix> >(
"A");
606 size_t nRowsFine = fineA->getGlobalNumRows();
607 size_t nRowsCoarse = A_->getGlobalNumRows();
609 SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
610 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
614 this->GetOStream(
Statistics1) << eigRatioString <<
" (computed) = " << ratio << std::endl;
615 paramList.set(eigRatioString, ratio);
631 if (lambdaMax == negone) {
632 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
634 Teuchos::RCP<Ifpack2::Chebyshev<MatrixType> > chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
635 if (chebyPrec != Teuchos::null) {
636 lambdaMax = chebyPrec->getLambdaMaxForApply();
637 A_->SetMaxEigenvalueEstimate(lambdaMax);
638 this->GetOStream(
Statistics1) <<
"chebyshev: max eigenvalue (calculated by Ifpack2)" <<
" = " << lambdaMax << std::endl;
640 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone,
Exceptions::RuntimeError,
"MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
644 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
646 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
648 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
654 bool reusePreconditioner =
false;
655 if (this->IsSetup() ==
true) {
657 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
659 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
660 if (!prec.is_null()) {
661 #ifdef IFPACK2_HAS_PROPER_REUSE
662 prec->resetMatrix(tA);
663 reusePreconditioner =
true;
665 this->GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
669 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
670 "reverting to full construction" << std::endl;
674 if (!reusePreconditioner) {
683 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
696 Teuchos::ParameterList paramList;
697 bool supportInitialGuess =
false;
698 if (type_ ==
"CHEBYSHEV") {
699 paramList.set(
"chebyshev: zero starting solution", InitialGuessIsZero);
700 SetPrecParameters(paramList);
701 supportInitialGuess =
true;
703 }
else if (type_ ==
"RELAXATION") {
704 paramList.set(
"relaxation: zero starting solution", InitialGuessIsZero);
705 SetPrecParameters(paramList);
706 supportInitialGuess =
true;
708 }
else if (type_ ==
"KRYLOV") {
709 paramList.set(
"krylov: zero starting solution", InitialGuessIsZero);
710 SetPrecParameters(paramList);
711 supportInitialGuess =
true;
713 }
else if (type_ ==
"SCHWARZ") {
714 paramList.set(
"schwarz: zero starting solution", InitialGuessIsZero);
719 prec_->setParameters(paramList);
720 supportInitialGuess =
true;
730 if (InitialGuessIsZero || supportInitialGuess) {
733 prec_->apply(tpB, tpX);
735 typedef Teuchos::ScalarTraits<Scalar> TST;
737 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
742 prec_->apply(tpB, tpX);
744 X.update(TST::one(), *Correction, TST::one());
748 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
751 smoother->SetParameterList(this->GetParameterList());
755 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
757 std::ostringstream out;
759 out << prec_->description();
762 out <<
"{type = " << type_ <<
"}";
767 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
772 out0 <<
"Prec. type: " << type_ << std::endl;
775 out0 <<
"Parameter list: " << std::endl;
776 Teuchos::OSTab tab2(out);
777 out << this->GetParameterList();
778 out0 <<
"Overlap: " << overlap_ << std::endl;
782 if (prec_ != Teuchos::null) {
783 Teuchos::OSTab tab2(out);
784 out << *prec_ << std::endl;
787 if (verbLevel &
Debug) {
790 <<
"RCP<prec_>: " << prec_ << std::endl;
794 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
796 typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
798 RCP<Ifpack2::Relaxation<MatrixType> > pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType> >(prec_);
799 if(!pr.is_null())
return pr->getNodeSmootherComplexity();
801 RCP<Ifpack2::Chebyshev<MatrixType> > pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
802 if(!pc.is_null())
return pc->getNodeSmootherComplexity();
804 RCP<Ifpack2::BlockRelaxation<MatrixType> > pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType> >(prec_);
805 if(!pb.is_null())
return pb->getNodeSmootherComplexity();
807 RCP<Ifpack2::ILUT<MatrixType> > pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType> >(prec_);
808 if(!pi.is_null())
return pi->getNodeSmootherComplexity();
810 RCP<Ifpack2::RILUK<MatrixType> > pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType> >(prec_);
811 if(!pk.is_null())
return pk->getNodeSmootherComplexity();
814 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)
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.
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.