46#ifndef MUELU_IFPACK2SMOOTHER_DEF_HPP
47#define MUELU_IFPACK2SMOOTHER_DEF_HPP
51#if 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_TpetraBlockCrsMatrix.hpp>
69#include <Xpetra_Matrix.hpp>
70#include <Xpetra_MatrixMatrix.hpp>
71#include <Xpetra_MultiVectorFactory.hpp>
72#include <Xpetra_TpetraMultiVector.hpp>
74#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
78#include "MueLu_Utilities.hpp"
80#include "MueLu_Aggregates.hpp"
83#ifdef HAVE_MUELU_INTREPID2
86#include "Intrepid2_Basis.hpp"
87#include "Kokkos_DynRankView.hpp"
93 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
97 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
98 bool isSupported = Ifpack2::Factory::isSupported<tRowMatrix>(
type_) || (
type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
99 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
100 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
101 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
102 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
103 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
104 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
105 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
106 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
107 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
108 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
109 type_ ==
"LINESMOOTHING_BLOCKRELAXATION" ||
110 type_ ==
"TOPOLOGICAL" ||
111 type_ ==
"AGGREGATE");
117 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
130 std::string prefix = this->
ShortClassName() +
": SetPrecParameters";
132 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
134 paramList.setParameters(list);
141 if(!precList.is_null() && precList->isParameter(
"partitioner: type") && precList->get<std::string>(
"partitioner: type") ==
"linear" &&
142 !precList->isParameter(
"partitioner: local parts")) {
143 LO matrixBlockSize = 1;
144 int lclSize =
A_->getRangeMap()->getLocalNumElements();
145 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_);
146 if (!matA.is_null()) {
147 lclSize = matA->getLocalNumRows();
148 matrixBlockSize = matA->GetFixedBlockSize();
150 precList->set(
"partitioner: local parts", lclSize / matrixBlockSize);
153 prec_->setParameters(*precList);
155 paramList.setParameters(*precList);
158 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
160 this->
Input(currentLevel,
"A");
162 if (
type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
163 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
164 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
165 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
166 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
167 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
168 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
169 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
170 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
171 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
172 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
173 type_ ==
"LINESMOOTHING_BLOCKRELAXATION") {
174 this->
Input(currentLevel,
"CoarseNumZLayers");
175 this->
Input(currentLevel,
"LineDetection_VertLineIds");
177 else if (
type_ ==
"BLOCK RELAXATION" ||
178 type_ ==
"BLOCK_RELAXATION" ||
179 type_ ==
"BLOCKRELAXATION" ||
181 type_ ==
"BANDED_RELAXATION" ||
182 type_ ==
"BANDED RELAXATION" ||
183 type_ ==
"BANDEDRELAXATION" ||
185 type_ ==
"TRIDI_RELAXATION" ||
186 type_ ==
"TRIDI RELAXATION" ||
187 type_ ==
"TRIDIRELAXATION" ||
188 type_ ==
"TRIDIAGONAL_RELAXATION" ||
189 type_ ==
"TRIDIAGONAL RELAXATION" ||
190 type_ ==
"TRIDIAGONALRELAXATION")
194 if(precList.isParameter(
"partitioner: type") &&
195 precList.get<std::string>(
"partitioner: type") ==
"line") {
196 this->
Input(currentLevel,
"Coordinates");
199 else if (
type_ ==
"TOPOLOGICAL")
202 this->
Input(currentLevel,
"pcoarsen: element to node map");
204 else if (
type_ ==
"AGGREGATE")
207 this->
Input(currentLevel,
"Aggregates");
209 else if (
type_ ==
"HIPTMAIR") {
211 this->
Input(currentLevel,
"NodeMatrix");
217 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
221 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
225 if(paramList.isParameter(
"smoother: use blockcrsmatrix storage") && paramList.get<
bool>(
"smoother: use blockcrsmatrix storage")) {
227 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_);
229 blocksize = matA->GetFixedBlockSize();
233 RCP<CrsMatrixWrap> AcrsWrap = rcp_dynamic_cast<CrsMatrixWrap>(
A_);
234 if(AcrsWrap.is_null())
235 throw std::runtime_error(
"Ifpack2Smoother: Cannot convert matrix A to CrsMatrixWrap object.");
236 RCP<CrsMatrix> Acrs = AcrsWrap->getCrsMatrix();
238 throw std::runtime_error(
"Ifpack2Smoother: Cannot extract CrsMatrix from matrix A.");
239 RCP<TpetraCrsMatrix> At = rcp_dynamic_cast<TpetraCrsMatrix>(Acrs);
241 if(!Xpetra::Helpers<Scalar,LO,GO,Node>::isTpetraBlockCrs(matA))
242 throw std::runtime_error(
"Ifpack2Smoother: Cannot extract CrsMatrix or BlockCrsMatrix from matrix A.");
243 this->
GetOStream(
Statistics0) <<
"Ifpack2Smoother: Using (native) BlockCrsMatrix storage with blocksize "<<blocksize<<std::endl;
244 paramList.remove(
"smoother: use blockcrsmatrix storage");
247 RCP<Tpetra::BlockCrsMatrix<Scalar, LO, GO, Node> > blockCrs = Tpetra::convertToBlockCrsMatrix(*At->getTpetra_CrsMatrix(),blocksize);
248 RCP<CrsMatrix> blockCrs_as_crs = rcp(
new TpetraBlockCrsMatrix(blockCrs));
249 RCP<CrsMatrixWrap> blockWrap = rcp(
new CrsMatrixWrap(blockCrs_as_crs));
251 this->
GetOStream(
Statistics0) <<
"Ifpack2Smoother: Using BlockCrsMatrix storage with blocksize "<<blocksize<<std::endl;
253 paramList.remove(
"smoother: use blockcrsmatrix storage");
258 if (
type_ ==
"SCHWARZ")
261 else if (
type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
262 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
263 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
264 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
265 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
266 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
267 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
268 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
269 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
270 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
271 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
272 type_ ==
"LINESMOOTHING_BLOCKRELAXATION")
275 else if (
type_ ==
"BLOCK_RELAXATION" ||
276 type_ ==
"BLOCK RELAXATION" ||
277 type_ ==
"BLOCKRELAXATION" ||
279 type_ ==
"BANDED_RELAXATION" ||
280 type_ ==
"BANDED RELAXATION" ||
281 type_ ==
"BANDEDRELAXATION" ||
283 type_ ==
"TRIDI_RELAXATION" ||
284 type_ ==
"TRIDI RELAXATION" ||
285 type_ ==
"TRIDIRELAXATION" ||
286 type_ ==
"TRIDIAGONAL_RELAXATION" ||
287 type_ ==
"TRIDIAGONAL RELAXATION" ||
288 type_ ==
"TRIDIAGONALRELAXATION")
291 else if (
type_ ==
"CHEBYSHEV")
294 else if (
type_ ==
"TOPOLOGICAL")
296#ifdef HAVE_MUELU_INTREPID2
299 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"'TOPOLOGICAL' smoother choice requires Intrepid2");
302 else if (
type_ ==
"AGGREGATE")
305 else if (
type_ ==
"HIPTMAIR")
318 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
320 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
322 bool reusePreconditioner =
false;
325 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
327 bool isTRowMatrix =
true;
328 RCP<const tRowMatrix> tA;
332 isTRowMatrix =
false;
335 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(
prec_);
336 if (!prec.is_null() && isTRowMatrix) {
338 reusePreconditioner =
true;
340 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
341 "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
345 if (!reusePreconditioner) {
346 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
348 bool isBlockedMatrix =
false;
349 RCP<Matrix> merged2Mat;
351 std::string sublistName =
"subdomain solver parameters";
352 if (paramList.isSublist(sublistName)) {
361 ParameterList& subList = paramList.sublist(sublistName);
363 std::string partName =
"partitioner: type";
369 if (subList.isParameter(partName) && subList.get<std::string>(partName) ==
"vanka user") {
370 isBlockedMatrix =
true;
372 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
374 "Matrix A must be of type BlockedCrsMatrix.");
376 size_t numVels = bA->getMatrix(0,0)->getLocalNumRows();
377 size_t numPres = bA->getMatrix(1,0)->getLocalNumRows();
378 size_t numRows = rcp_dynamic_cast<Matrix>(
A_,
true)->getLocalNumRows();
380 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
382 size_t numBlocks = 0;
383 for (
size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
384 blockSeeds[rowOfB] = numBlocks++;
386 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
388 "Matrix A must be of type BlockedCrsMatrix.");
390 merged2Mat = bA2->Merge();
394 bool haveBoundary =
false;
395 for (LO i = 0; i < boundaryNodes.size(); i++)
396 if (boundaryNodes[i]) {
400 blockSeeds[i] = numBlocks;
406 subList.set(
"partitioner: type",
"user");
407 subList.set(
"partitioner: map", blockSeeds);
408 subList.set(
"partitioner: local parts", as<int>(numBlocks));
411 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
413 isBlockedMatrix =
true;
414 merged2Mat = bA->Merge();
419 RCP<const tRowMatrix> tA;
434 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
437 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
440 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupAggregate(): Setup() has already been called" << std::endl;
441 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupAggregate(): reuse of this type is not available, reverting to full construction" << std::endl;
448 RCP<const LOMultiVector> vertex2AggId = aggregates->GetVertex2AggId();
449 ArrayRCP<LO> aggregate_ids = rcp_const_cast<LOMultiVector>(vertex2AggId)->getDataNonConst(0);
450 ArrayRCP<LO> dof_ids;
454 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_);
456 blocksize = matA->GetFixedBlockSize();
458 dof_ids.resize(aggregate_ids.size() * blocksize);
459 for(LO i=0; i<(LO)aggregate_ids.size(); i++) {
460 for(LO j=0; j<(LO)blocksize; j++)
461 dof_ids[i*blocksize+j] = aggregate_ids[i];
465 dof_ids = aggregate_ids;
469 paramList.set(
"partitioner: map", dof_ids);
470 paramList.set(
"partitioner: type",
"user");
471 paramList.set(
"partitioner: overlap", 0);
472 paramList.set(
"partitioner: local parts", (
int)aggregates->GetNumAggregates());
476 type_ =
"BLOCKRELAXATION";
483#ifdef HAVE_MUELU_INTREPID2
484 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
494 if (this->IsSetup() ==
true) {
495 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
496 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
499 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
501 typedef typename Node::device_type::execution_space ES;
503 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO;
505 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
511 string basisString = paramList.get<
string>(
"pcoarsen: hi basis");
519 string topologyTypeString = paramList.get<
string>(
"smoother: neighborhood type");
521 if (topologyTypeString ==
"node")
523 else if (topologyTypeString ==
"edge")
525 else if (topologyTypeString ==
"face")
527 else if (topologyTypeString ==
"cell")
528 dimension = basis->getBaseCellTopology().getDimension();
530 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
531 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_,
true);
532 vector<vector<LocalOrdinal>> seeds;
537 int myNodeCount = matA->getRowMap()->getLocalNumElements();
538 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
539 int localPartitionNumber = 0;
542 nodeSeeds[seed] = localPartitionNumber++;
545 paramList.remove(
"smoother: neighborhood type");
546 paramList.remove(
"pcoarsen: hi basis");
548 paramList.set(
"partitioner: map", nodeSeeds);
549 paramList.set(
"partitioner: type",
"user");
550 paramList.set(
"partitioner: overlap", 1);
551 paramList.set(
"partitioner: local parts",
int(seeds[dimension].size()));
555 type_ =
"BLOCKRELAXATION";
563 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
566 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
567 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
570 ParameterList& myparamList =
const_cast<ParameterList&
>(this->
GetParameterList());
573 if (CoarseNumZLayers > 0) {
578 for(
size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
579 if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
581 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_,
true);
582 size_t numLocalRows = matA->getLocalNumRows();
585 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
590 int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
591 LO matrixBlockSize = matA->GetFixedBlockSize();
592 if(matrixBlockSize > 1 && actualDofsPerNode > 1)
595 "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
597 else if(matrixBlockSize > 1)
599 actualDofsPerNode = matrixBlockSize;
601 myparamList.set(
"partitioner: PDE equations", actualDofsPerNode);
603 if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
604 myparamList.set(
"partitioner: type",
"user");
605 myparamList.set(
"partitioner: map",TVertLineIdSmoo);
606 myparamList.set(
"partitioner: local parts",maxPart+1);
609 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
612 Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
613 for (
size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
614 for (
size_t dof = 0; dof < numDofsPerNode; dof++)
615 partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
616 myparamList.set(
"partitioner: type",
"user");
617 myparamList.set(
"partitioner: map",partitionerMap);
618 myparamList.set(
"partitioner: local parts",maxPart + 1);
621 if (
type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
622 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
623 type_ ==
"LINESMOOTHING_BANDEDRELAXATION")
624 type_ =
"BANDEDRELAXATION";
625 else if (
type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
626 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
627 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
628 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
629 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
630 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION")
631 type_ =
"TRIDIAGONALRELAXATION";
633 type_ =
"BLOCKRELAXATION";
636 this->
GetOStream(
Runtime0) <<
"Line detection failed: fall back to point-wise relaxation" << std::endl;
637 myparamList.remove(
"partitioner: type",
false);
638 myparamList.remove(
"partitioner: map",
false);
639 myparamList.remove(
"partitioner: local parts",
false);
640 type_ =
"RELAXATION";
651 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
653 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
655 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
661 bool reusePreconditioner =
false;
664 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
666 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(
prec_);
667 if (!prec.is_null()) {
669 reusePreconditioner =
true;
671 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
672 "reverting to full construction" << std::endl;
676 if (!reusePreconditioner) {
677 ParameterList& myparamList =
const_cast<ParameterList&
>(this->
GetParameterList());
679 if(myparamList.isParameter(
"partitioner: type") &&
680 myparamList.get<std::string>(
"partitioner: type") ==
"line") {
681 Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > xCoordinates =
683 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));
685 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_);
686 size_t lclSize =
A_->getRangeMap()->getLocalNumElements();
688 lclSize = matA->getLocalNumRows();
689 size_t numDofsPerNode = lclSize / xCoordinates->getMap()->getLocalNumElements();
690 myparamList.set(
"partitioner: coordinates", coordinates);
691 myparamList.set(
"partitioner: PDE equations", (
int) numDofsPerNode);
702 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
704 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
705 using STS = Teuchos::ScalarTraits<SC>;
706 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
712 bool reusePreconditioner =
false;
716 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): Setup() has already been called, assuming reuse" << std::endl;
718 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(
prec_);
719 if (!prec.is_null()) {
721 reusePreconditioner =
true;
723 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available (failed cast to CanChangeMatrix), "
724 "reverting to full construction" << std::endl;
729 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
730 SC negone = -STS::one();
734 if (!reusePreconditioner) {
749 if (lambdaMax == negone) {
750 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
752 Teuchos::RCP<Ifpack2::Chebyshev<MatrixType> > chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(
prec_);
753 if (chebyPrec != Teuchos::null) {
754 lambdaMax = chebyPrec->getLambdaMaxForApply();
755 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_);
757 matA->SetMaxEigenvalueEstimate(lambdaMax);
758 this->
GetOStream(
Statistics1) <<
"chebyshev: max eigenvalue (calculated by Ifpack2)" <<
" = " << lambdaMax << std::endl;
760 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone,
Exceptions::RuntimeError,
"MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
765 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
767 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
768 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
774 bool reusePreconditioner =
false;
777 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupHiptmair(): Setup() has already been called, assuming reuse" << std::endl;
779 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(
prec_);
780 if (!prec.is_null()) {
782 reusePreconditioner =
true;
784 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupHiptmair(): reuse of this type is not available (failed cast to CanChangeMatrix), "
785 "reverting to full construction" << std::endl;
790 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
791 std::string smoother1 = paramList.get(
"hiptmair: smoother type 1",
"CHEBYSHEV");
792 std::string smoother2 = paramList.get(
"hiptmair: smoother type 2",
"CHEBYSHEV");
795 if(smoother1 ==
"CHEBYSHEV") {
796 ParameterList & list1 = paramList.sublist(
"hiptmair: smoother list 1");
800 if(smoother2 ==
"CHEBYSHEV") {
801 ParameterList & list2 = paramList.sublist(
"hiptmair: smoother list 2");
810 RCP<Operator> NodeMatrix = currentLevel.
Get<RCP<Operator> >(
"NodeMatrix");
811 RCP<Operator> D0 = currentLevel.
Get<RCP<Operator> >(
"D0");
817 Teuchos::ParameterList newlist;
818 newlist.set(
"P",tD0);
819 newlist.set(
"PtAP",tNodeMatrix);
820 if (!reusePreconditioner) {
831 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
834 typedef Teuchos::ScalarTraits<SC> STS;
835 SC negone = -STS::one();
836 RCP<const Matrix> currentA = currentLevel.
Get<RCP<Matrix> >(matrixName);
837 SC lambdaMax = negone;
839 std::string maxEigString =
"chebyshev: max eigenvalue";
840 std::string eigRatioString =
"chebyshev: ratio eigenvalue";
843 if (paramList.isParameter(maxEigString)) {
844 if (paramList.isType<
double>(maxEigString))
845 lambdaMax = paramList.get<
double>(maxEigString);
847 lambdaMax = paramList.get<SC>(maxEigString);
848 this->
GetOStream(
Statistics1) << label << maxEigString <<
" (cached with smoother parameter list) = " << lambdaMax << std::endl;
849 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(
A_);
851 matA->SetMaxEigenvalueEstimate(lambdaMax);
854 lambdaMax = currentA->GetMaxEigenvalueEstimate();
855 if (lambdaMax != negone) {
856 this->
GetOStream(
Statistics1) << label << maxEigString <<
" (cached with matrix) = " << lambdaMax << std::endl;
857 paramList.set(maxEigString, lambdaMax);
862 const SC defaultEigRatio = 20;
864 SC ratio = defaultEigRatio;
865 if (paramList.isParameter(eigRatioString)) {
866 if (paramList.isType<
double>(eigRatioString))
867 ratio = paramList.get<
double>(eigRatioString);
869 ratio = paramList.get<SC>(eigRatioString);
876 RCP<const Matrix> fineA = currentLevel.
GetPreviousLevel()->Get<RCP<Matrix> >(matrixName);
877 size_t nRowsFine = fineA->getGlobalNumRows();
878 size_t nRowsCoarse = currentA->getGlobalNumRows();
880 SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
881 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
886 paramList.set(eigRatioString, ratio);
888 if (paramList.isParameter(
"chebyshev: use rowsumabs diagonal scaling")) {
889 this->
GetOStream(
Runtime1) <<
"chebyshev: using rowsumabs diagonal scaling" << std::endl;
890 bool doScale =
false;
891 doScale = paramList.get<
bool>(
"chebyshev: use rowsumabs diagonal scaling");
892 paramList.remove(
"chebyshev: use rowsumabs diagonal scaling");
893 double chebyReplaceTol = Teuchos::ScalarTraits<Scalar>::eps()*100;
894 std::string paramName =
"chebyshev: rowsumabs diagonal replacement tolerance";
895 if (paramList.isParameter(paramName)) {
896 chebyReplaceTol = paramList.get<
double>(paramName);
897 paramList.remove(paramName);
899 double chebyReplaceVal = Teuchos::ScalarTraits<double>::zero();
900 paramName =
"chebyshev: rowsumabs diagonal replacement value";
901 if (paramList.isParameter(paramName)) {
902 chebyReplaceVal = paramList.get<
double>(paramName);
903 paramList.remove(paramName);
905 bool chebyReplaceSingleEntryRowWithZero =
false;
906 paramName =
"chebyshev: rowsumabs replace single entry row with zero";
907 if (paramList.isParameter(paramName)) {
908 chebyReplaceSingleEntryRowWithZero = paramList.get<
bool>(paramName);
909 paramList.remove(paramName);
911 bool useAverageAbsDiagVal =
false;
912 paramName =
"chebyshev: rowsumabs use automatic diagonal tolerance";
913 if (paramList.isParameter(paramName)) {
914 useAverageAbsDiagVal = paramList.get<
bool>(paramName);
915 paramList.remove(paramName);
918 const bool doReciprocal =
true;
919 RCP<Vector> lumpedDiagonal =
Utilities::GetLumpedMatrixDiagonal(*currentA, doReciprocal, chebyReplaceTol, chebyReplaceVal, chebyReplaceSingleEntryRowWithZero, useAverageAbsDiagVal);
920 const Xpetra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(*lumpedDiagonal);
921 paramList.set(
"chebyshev: operator inv diagonal",tmpVec.getTpetra_Vector());
931 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
933 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
934 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
940 bool reusePreconditioner =
false;
943 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
945 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(
prec_);
946 if (!prec.is_null()) {
948 reusePreconditioner =
true;
950 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
951 "reverting to full construction" << std::endl;
955 if (!reusePreconditioner) {
964 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
974 Teuchos::ParameterList paramList;
975 bool supportInitialGuess =
false;
978 if (
prec_->supportsZeroStartingSolution()) {
979 prec_->setZeroStartingSolution(InitialGuessIsZero);
980 supportInitialGuess =
true;
981 }
else if (
type_ ==
"SCHWARZ") {
982 paramList.set(
"schwarz: zero starting solution", InitialGuessIsZero);
987 prec_->setParameters(paramList);
988 supportInitialGuess =
true;
998 if (InitialGuessIsZero || supportInitialGuess) {
1001 prec_->apply(tpB, tpX);
1003 typedef Teuchos::ScalarTraits<Scalar> TST;
1005 RCP<MultiVector> Residual;
1008 RCP<TimeMonitor> tM = rcp(
new TimeMonitor(*
this, prefix +
"residual calculation",
Timings0));
1012 RCP<MultiVector> Correction = MultiVectorFactory::Build(
A_->getDomainMap(), X.getNumVectors());
1017 prec_->apply(tpB, tpX);
1019 X.update(TST::one(), *Correction, TST::one());
1023 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1030 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1032 std::ostringstream out;
1034 out <<
prec_->description();
1037 out <<
"{type = " <<
type_ <<
"}";
1042 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1047 out0 <<
"Prec. type: " <<
type_ << std::endl;
1050 out0 <<
"Parameter list: " << std::endl;
1051 Teuchos::OSTab tab2(out);
1053 out0 <<
"Overlap: " <<
overlap_ << std::endl;
1057 if (
prec_ != Teuchos::null) {
1058 Teuchos::OSTab tab2(out);
1059 out << *
prec_ << std::endl;
1062 if (verbLevel &
Debug) {
1065 <<
"RCP<prec_>: " <<
prec_ << std::endl;
1069 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1071 typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
1073 RCP<Ifpack2::Relaxation<MatrixType> > pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType> >(
prec_);
1074 if(!pr.is_null())
return pr->getNodeSmootherComplexity();
1076 RCP<Ifpack2::Chebyshev<MatrixType> > pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(
prec_);
1077 if(!pc.is_null())
return pc->getNodeSmootherComplexity();
1079 RCP<Ifpack2::BlockRelaxation<MatrixType> > pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType> >(
prec_);
1080 if(!pb.is_null())
return pb->getNodeSmootherComplexity();
1082 RCP<Ifpack2::ILUT<MatrixType> > pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType> >(
prec_);
1083 if(!pi.is_null())
return pi->getNodeSmootherComplexity();
1085 RCP<Ifpack2::RILUK<MatrixType> > pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType> >(
prec_);
1086 if(!pk.is_null())
return pk->getNodeSmootherComplexity();
1089 return Teuchos::OrdinalTraits<size_t>::invalid();
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
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 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 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.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
RCP< ParameterList > RemoveFactoriesFromList(const ParameterList &list) const
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.
Scalar SetupChebyshevEigenvalues(Level ¤tLevel, const std::string &matrixName, const std::string &label, Teuchos::ParameterList ¶mList) const
void SetupAggregate(Level ¤tLevel)
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 SetupHiptmair(Level ¤tLevel)
LO overlap_
overlap when using the smoother in additive Schwarz mode
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
void SetupTopological(Level ¤tLevel)
RCP< Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node > > prec_
pointer to Ifpack2 preconditioner object
std::string type_
ifpack2-specific key phrase that denote smoother type
friend class Ifpack2Smoother
Constructor.
void DeclareInput(Level ¤tLevel) const
Input.
RCP< Operator > A_
matrix, used in apply if solving residual equation
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.
RCP< Level > & GetPreviousLevel()
Previous level.
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)....
virtual const Teuchos::ParameterList & GetParameterList() const
virtual void SetParameterList(const Teuchos::ParameterList ¶mList)=0
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.
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), bool count_twos_as_dirichlet=false)
static Teuchos::RCP< Vector > GetLumpedMatrixDiagonal(Matrix const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::magnitude(Teuchos::ScalarTraits< Scalar >::zero()), Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false, const bool useAverageAbsDiagVal=false)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
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)
Teuchos::RCP< Intrepid2::Basis< KokkosExecutionSpace, Scalar, Scalar > > BasisFactory(const std::string &name, int °ree)
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.
@ 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).
@ Timings
Print all timing information.
@ Parameters0
Print class parameters.
@ Statistics0
Print statistics that do not involve significant additional computation.
@ Parameters1
Print class parameters (more parameters, more verbose).