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> 73 #include "MueLu_Utilities.hpp" 76 #ifdef HAVE_MUELU_INTREPID2 79 #include "Intrepid2_Basis.hpp" 81 #ifdef HAVE_MUELU_INTREPID2_REFACTOR 84 #include "Intrepid2_FieldContainer.hpp" 92 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
94 : type_(type), overlap_(overlap)
99 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
112 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
113 paramList.setParameters(list);
117 prec_->setParameters(*precList);
119 paramList.setParameters(*precList);
122 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
124 this->
Input(currentLevel,
"A");
126 if (
type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
127 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
128 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
129 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
130 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
131 type_ ==
"LINESMOOTHING_BLOCKRELAXATION") {
132 this->
Input(currentLevel,
"CoarseNumZLayers");
133 this->
Input(currentLevel,
"LineDetection_VertLineIds");
135 else if (
type_ ==
"TOPOLOGICAL")
138 this->
Input(currentLevel,
"pcoarsen: element to node map");
142 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
146 A_ = Factory::Get< RCP<Matrix> >(currentLevel,
"A");
148 if (
type_ ==
"SCHWARZ")
151 else if (
type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
152 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
153 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
154 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
155 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
156 type_ ==
"LINESMOOTHING_BLOCKRELAXATION")
159 else if (
type_ ==
"CHEBYSHEV")
162 else if (
type_ ==
"TOPOLOGICAL")
164 #ifdef HAVE_MUELU_INTREPID2 167 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"'TOPOLOGICAL' smoother choice requires Intrepid2");
180 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
182 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
184 bool reusePreconditioner =
false;
187 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
189 bool isTRowMatrix =
true;
190 RCP<const tRowMatrix> tA;
194 isTRowMatrix =
false;
198 if (!prec.is_null() && isTRowMatrix) {
199 #ifdef IFPACK2_HAS_PROPER_REUSE 200 prec->resetMatrix(tA);
201 reusePreconditioner =
true;
203 this->
GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
207 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available " 208 "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
212 if (!reusePreconditioner) {
213 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
215 bool isBlockedMatrix =
false;
216 RCP<Matrix> merged2Mat;
218 std::string sublistName =
"subdomain solver parameters";
219 if (paramList.isSublist(sublistName)) {
228 ParameterList& subList = paramList.sublist(sublistName);
230 std::string partName =
"partitioner: type";
231 if (subList.isParameter(partName) && subList.get<std::string>(partName) ==
"user") {
232 isBlockedMatrix =
true;
234 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
236 "Matrix A must be of type BlockedCrsMatrix.");
238 size_t numVels = bA->getMatrix(0,0)->getNodeNumRows();
239 size_t numPres = bA->getMatrix(1,0)->getNodeNumRows();
240 size_t numRows =
A_->getNodeNumRows();
242 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
244 size_t numBlocks = 0;
245 for (
size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
246 blockSeeds[rowOfB] = numBlocks++;
248 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
250 "Matrix A must be of type BlockedCrsMatrix.");
252 merged2Mat = bA2->Merge();
256 bool haveBoundary =
false;
257 for (LO i = 0; i < boundaryNodes.size(); i++)
258 if (boundaryNodes[i]) {
262 blockSeeds[i] = numBlocks;
268 subList.set(
"partitioner: map", blockSeeds);
269 subList.set(
"partitioner: local parts", as<int>(numBlocks));
272 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
274 isBlockedMatrix =
true;
275 merged2Mat = bA->Merge();
280 RCP<const tRowMatrix> tA;
293 #ifdef HAVE_MUELU_INTREPID2 294 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
305 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
306 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
309 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
311 typedef typename Node::device_type::execution_space ES;
313 #ifdef HAVE_MUELU_INTREPID2_REFACTOR 314 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO;
316 typedef Intrepid2::FieldContainer<LocalOrdinal> FCO;
319 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
323 const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO> >(currentLevel,
"pcoarsen: element to node map");
325 string basisString = paramList.get<
string>(
"pcoarsen: hi basis");
331 #ifdef HAVE_MUELU_INTREPID2_REFACTOR 332 auto basis = MueLuIntrepid::BasisFactory<double,ES>(basisString, degree);
334 auto basis = MueLuIntrepid::BasisFactory<double>(basisString, degree);
337 string topologyTypeString = paramList.get<
string>(
"smoother: neighborhood type");
339 if (topologyTypeString ==
"node")
341 else if (topologyTypeString ==
"edge")
343 else if (topologyTypeString ==
"face")
345 else if (topologyTypeString ==
"cell")
346 dimension = basis->getBaseCellTopology().getDimension();
348 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
349 vector<vector<LocalOrdinal>> seeds;
354 int myNodeCount =
A_->getRowMap()->getNodeNumElements();
355 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
356 int localPartitionNumber = 0;
357 for (LocalOrdinal seed : seeds[dimension])
359 nodeSeeds[seed] = localPartitionNumber++;
362 paramList.remove(
"smoother: neighborhood type");
363 paramList.remove(
"pcoarsen: hi basis");
365 paramList.set(
"partitioner: map", nodeSeeds);
366 paramList.set(
"partitioner: type",
"user");
367 paramList.set(
"partitioner: overlap", 1);
368 paramList.set(
"partitioner: local parts",
int(seeds[dimension].size()));
372 type_ =
"BLOCKRELAXATION";
380 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
383 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
384 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
387 ParameterList& myparamList =
const_cast<ParameterList&
>(this->
GetParameterList());
389 LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,
"CoarseNumZLayers");
390 if (CoarseNumZLayers > 0) {
391 Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get< Teuchos::ArrayRCP<LO> >(currentLevel,
"LineDetection_VertLineIds");
395 for(
size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
396 if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
399 size_t numLocalRows =
A_->getNodeNumRows();
401 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
403 if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
404 myparamList.set(
"partitioner: type",
"user");
405 myparamList.set(
"partitioner: map",TVertLineIdSmoo);
406 myparamList.set(
"partitioner: local parts",maxPart+1);
409 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
412 Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
413 for (
size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
414 for (
size_t dof = 0; dof < numDofsPerNode; dof++)
415 partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
416 myparamList.set(
"partitioner: type",
"user");
417 myparamList.set(
"partitioner: map",partitionerMap);
418 myparamList.set(
"partitioner: local parts",maxPart + 1);
421 if (
type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
422 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
423 type_ ==
"LINESMOOTHING_BANDEDRELAXATION")
424 type_ =
"BANDEDRELAXATION";
426 type_ =
"BLOCKRELAXATION";
429 this->
GetOStream(
Runtime0) <<
"Line detection failed: fall back to point-wise relaxation" << std::endl;
430 myparamList.remove(
"partitioner: type",
false);
431 myparamList.remove(
"partitioner: map",
false);
432 myparamList.remove(
"partitioner: local parts",
false);
433 type_ =
"RELAXATION";
444 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
447 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): SetupChebyshev() has already been called" << std::endl;
448 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available, reverting to full construction" << std::endl;
451 typedef Teuchos::ScalarTraits<SC> STS;
452 SC negone = -STS::one();
454 SC lambdaMax = negone;
456 std::string maxEigString =
"chebyshev: max eigenvalue";
457 std::string eigRatioString =
"chebyshev: ratio eigenvalue";
459 ParameterList& paramList =
const_cast<ParameterList&
>(this->
GetParameterList());
462 if (paramList.isParameter(maxEigString)) {
463 if (paramList.isType<
double>(maxEigString))
464 lambdaMax = paramList.get<
double>(maxEigString);
466 lambdaMax = paramList.get<SC>(maxEigString);
467 this->
GetOStream(
Statistics1) << maxEigString <<
" (cached with smoother parameter list) = " << lambdaMax << std::endl;
470 lambdaMax =
A_->GetMaxEigenvalueEstimate();
471 if (lambdaMax != negone) {
472 this->
GetOStream(
Statistics1) << maxEigString <<
" (cached with matrix) = " << lambdaMax << std::endl;
473 paramList.set(maxEigString, lambdaMax);
478 const SC defaultEigRatio = 20;
480 SC ratio = defaultEigRatio;
481 if (paramList.isParameter(eigRatioString)) {
482 if (paramList.isType<
double>(eigRatioString))
483 ratio = paramList.get<
double>(eigRatioString);
485 ratio = paramList.get<SC>(eigRatioString);
492 RCP<const Matrix> fineA = currentLevel.
GetPreviousLevel()->Get<RCP<Matrix> >(
"A");
493 size_t nRowsFine = fineA->getGlobalNumRows();
494 size_t nRowsCoarse =
A_->getGlobalNumRows();
496 SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
497 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
502 paramList.set(eigRatioString, ratio);
512 if (lambdaMax == negone) {
513 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
516 if (chebyPrec != Teuchos::null) {
518 A_->SetMaxEigenvalueEstimate(lambdaMax);
519 this->
GetOStream(
Statistics1) <<
"chebyshev: max eigenvalue (calculated by Ifpack2)" <<
" = " << lambdaMax << std::endl;
521 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone,
Exceptions::RuntimeError,
"MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
525 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
527 typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
529 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(
A_);
535 bool reusePreconditioner =
false;
538 this->
GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
541 if (!prec.is_null()) {
542 #ifdef IFPACK2_HAS_PROPER_REUSE 543 prec->resetMatrix(tA);
544 reusePreconditioner =
true;
546 this->
GetOStream(
Errors) <<
"Ifpack2 does not have proper reuse yet." << std::endl;
550 this->
GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available (failed cast to CanChangeMatrix), " 551 "reverting to full construction" << std::endl;
555 if (!reusePreconditioner) {
564 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
577 Teuchos::ParameterList paramList;
578 bool supportInitialGuess =
false;
579 if (
type_ ==
"CHEBYSHEV") {
580 paramList.set(
"chebyshev: zero starting solution", InitialGuessIsZero);
582 supportInitialGuess =
true;
584 }
else if (
type_ ==
"RELAXATION") {
585 paramList.set(
"relaxation: zero starting solution", InitialGuessIsZero);
587 supportInitialGuess =
true;
589 }
else if (
type_ ==
"KRYLOV") {
590 paramList.set(
"krylov: zero starting solution", InitialGuessIsZero);
592 supportInitialGuess =
true;
594 }
else if (
type_ ==
"SCHWARZ") {
595 paramList.set(
"schwarz: zero starting solution", InitialGuessIsZero);
600 prec_->setParameters(paramList);
601 supportInitialGuess =
true;
611 if (InitialGuessIsZero || supportInitialGuess) {
614 prec_->apply(tpB, tpX);
616 typedef Teuchos::ScalarTraits<Scalar> TST;
618 RCP<MultiVector> Correction = MultiVectorFactory::Build(
A_->getDomainMap(), X.getNumVectors());
623 prec_->apply(tpB, tpX);
625 X.update(TST::one(), *Correction, TST::one());
629 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
636 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
638 std::ostringstream out;
640 out <<
prec_->description();
643 out <<
"{type = " <<
type_ <<
"}";
648 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
653 out0 <<
"Prec. type: " <<
type_ << std::endl;
656 out0 <<
"Parameter list: " << std::endl;
657 Teuchos::OSTab tab2(out);
659 out0 <<
"Overlap: " <<
overlap_ << std::endl;
663 if (
prec_ != Teuchos::null) {
664 Teuchos::OSTab tab2(out);
665 out << *
prec_ << std::endl;
668 if (verbLevel &
Debug) {
671 <<
"RCP<prec_>: " <<
prec_ << std::endl;
675 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
677 typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
692 return Teuchos::OrdinalTraits<size_t>::invalid();
698 #endif // HAVE_MUELU_TPETRA && HAVE_MUELU_IFPACK2 699 #endif // MUELU_IFPACK2SMOOTHER_DEF_HPP Important warning messages (one line)
void SetupGeneric(Level ¤tLevel)
RCP< SmootherPrototype > Copy() const
Exception indicating invalid cast attempted.
RCP< Level > & GetPreviousLevel()
Previous level.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
size_t getNodeSmootherComplexity() const
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)
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.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
std::string description() const
Return a simple one-line description of this object.
Print external lib objects.
bool IsSetup() const
Get the state of a smoother prototype.
friend class Ifpack2Smoother
Constructor.
Timer to be used in factories. Similar to Monitor but with additional timers.
RCP< Matrix > A_
matrix, used in apply if solving residual equation
RCP< ParameterList > RemoveFactoriesFromList(const ParameterList &list) const
Print additional debugging information.
One-liner description of what is happening.
void SetParameterList(const Teuchos::ParameterList ¶mList)
Namespace for MueLu classes and methods.
LO overlap_
overlap when using the smoother in additive Schwarz mode
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
int GetLevelID() const
Return level number.
std::string type_
ifpack2-specific key phrase that denote smoother type
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetupSchwarz(Level ¤tLevel)
void SetupLineSmoothing(Level ¤tLevel)
MatrixType::scalar_type getLambdaMaxForApply() const
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)
virtual void SetParameterList(const ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
Class that holds all level-specific information.
void Setup(Level ¤tLevel)
Set up the smoother.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
void SetupTopological(Level ¤tLevel)
size_t getNodeSmootherComplexity() const
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)
void DeclareInput(Level ¤tLevel) const
Input.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
virtual const Teuchos::ParameterList & GetParameterList() const
void Input(Level &level, const std::string &varName) const
size_t getNodeSmootherComplexity() const
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::magnitude(0.))
Print class parameters (more parameters, more verbose)
size_t getNodeSmootherComplexity() const
Exception throws to report errors in the internal logical of the program.
void SetupChebyshev(Level ¤tLevel)
Description of what is happening (more verbose)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
virtual std::string description() const
Return a simple one-line description of this object.
RCP< Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node > > prec_
pointer to Ifpack2 preconditioner object