10#ifndef MUELU_IFPACK2SMOOTHER_DEF_HPP
11#define MUELU_IFPACK2SMOOTHER_DEF_HPP
14#include "Xpetra_TpetraRowMatrix.hpp"
16#include <Teuchos_ParameterList.hpp>
18#include <Tpetra_RowMatrix.hpp>
20#include <Ifpack2_Chebyshev.hpp>
21#include <Ifpack2_Hiptmair.hpp>
22#include <Ifpack2_RILUK.hpp>
23#include <Ifpack2_Relaxation.hpp>
24#include <Ifpack2_ILUT.hpp>
25#include <Ifpack2_BlockRelaxation.hpp>
26#include <Ifpack2_Factory.hpp>
27#include <Ifpack2_Parameters.hpp>
29#include <Xpetra_BlockedCrsMatrix.hpp>
30#include <Xpetra_CrsMatrix.hpp>
31#include <Xpetra_CrsMatrixWrap.hpp>
32#include <Xpetra_TpetraBlockCrsMatrix.hpp>
33#include <Xpetra_Matrix.hpp>
34#include <Xpetra_MatrixMatrix.hpp>
35#include <Xpetra_MultiVectorFactory.hpp>
36#include <Xpetra_TpetraMultiVector.hpp>
38#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
42#include "MueLu_Utilities.hpp"
44#include "MueLu_Aggregates.hpp"
46#ifdef HAVE_MUELU_INTREPID2
49#include "Intrepid2_Basis.hpp"
50#include "Kokkos_DynRankView.hpp"
55template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
65 bool isSupported = Ifpack2::Factory::isSupported<tRowMatrix>(
type_) || (
type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
66 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
67 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
68 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
69 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
70 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
71 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
72 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
73 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
74 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
75 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
76 type_ ==
"LINESMOOTHING_BLOCKRELAXATION" ||
77 type_ ==
"TOPOLOGICAL" ||
78 type_ ==
"AGGREGATE");
84template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
95template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
97 std::string prefix = this->ShortClassName() +
": SetPrecParameters";
99 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
101 paramList.setParameters(list);
103 RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
106 precList->set(
"timer for apply", this->IsPrint(
Timings));
108 if (!precList.is_null() && precList->isParameter(
"partitioner: type") && precList->get<std::string>(
"partitioner: type") ==
"linear" &&
109 !precList->isParameter(
"partitioner: local parts")) {
110 LO matrixBlockSize = 1;
111 int lclSize = A_->getRangeMap()->getLocalNumElements();
112 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
113 if (!matA.is_null()) {
114 lclSize = matA->getLocalNumRows();
115 matrixBlockSize = matA->GetFixedBlockSize();
117 precList->set(
"partitioner: local parts", lclSize / matrixBlockSize);
120 prec_->setParameters(*precList);
122 paramList.setParameters(*precList);
125template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 this->Input(currentLevel,
"A");
129 if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
130 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
131 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
132 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
133 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
134 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
135 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
136 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
137 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
138 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
139 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
140 type_ ==
"LINESMOOTHING_BLOCKRELAXATION") {
141 this->Input(currentLevel,
"CoarseNumZLayers");
142 this->Input(currentLevel,
"LineDetection_VertLineIds");
143 }
else if (type_ ==
"BLOCK RELAXATION" ||
144 type_ ==
"BLOCK_RELAXATION" ||
145 type_ ==
"BLOCKRELAXATION" ||
147 type_ ==
"BANDED_RELAXATION" ||
148 type_ ==
"BANDED RELAXATION" ||
149 type_ ==
"BANDEDRELAXATION" ||
151 type_ ==
"TRIDI_RELAXATION" ||
152 type_ ==
"TRIDI RELAXATION" ||
153 type_ ==
"TRIDIRELAXATION" ||
154 type_ ==
"TRIDIAGONAL_RELAXATION" ||
155 type_ ==
"TRIDIAGONAL RELAXATION" ||
156 type_ ==
"TRIDIAGONALRELAXATION") {
158 ParameterList precList = this->GetParameterList();
159 if (precList.isParameter(
"partitioner: type") &&
160 precList.get<std::string>(
"partitioner: type") ==
"line") {
161 this->Input(currentLevel,
"Coordinates");
163 }
else if (type_ ==
"TOPOLOGICAL") {
165 this->Input(currentLevel,
"pcoarsen: element to node map");
166 }
else if (type_ ==
"AGGREGATE") {
168 this->Input(currentLevel,
"Aggregates");
169 }
else if (type_ ==
"HIPTMAIR") {
171 this->Input(currentLevel,
"NodeMatrix");
172 this->Input(currentLevel,
"D0");
176template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
179 A_ = Factory::Get<RCP<Operator>>(currentLevel,
"A");
180 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
184 if (paramList.isParameter(
"smoother: use blockcrsmatrix storage") && paramList.get<
bool>(
"smoother: use blockcrsmatrix storage")) {
186 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
188 blocksize = matA->GetFixedBlockSize();
192 RCP<CrsMatrixWrap> AcrsWrap = rcp_dynamic_cast<CrsMatrixWrap>(A_);
193 if (AcrsWrap.is_null())
194 throw std::runtime_error(
"Ifpack2Smoother: Cannot convert matrix A to CrsMatrixWrap object.");
195 RCP<CrsMatrix> Acrs = AcrsWrap->getCrsMatrix();
197 throw std::runtime_error(
"Ifpack2Smoother: Cannot extract CrsMatrix from matrix A.");
198 RCP<TpetraCrsMatrix> At = rcp_dynamic_cast<TpetraCrsMatrix>(Acrs);
200 if (!Xpetra::Helpers<Scalar, LO, GO, Node>::isTpetraBlockCrs(matA))
201 throw std::runtime_error(
"Ifpack2Smoother: Cannot extract CrsMatrix or BlockCrsMatrix from matrix A.");
202 this->GetOStream(
Statistics0) <<
"Ifpack2Smoother: Using (native) BlockCrsMatrix storage with blocksize " << blocksize << std::endl;
203 paramList.remove(
"smoother: use blockcrsmatrix storage");
205 RCP<Tpetra::BlockCrsMatrix<Scalar, LO, GO, Node>> blockCrs = Tpetra::convertToBlockCrsMatrix(*At->getTpetra_CrsMatrix(), blocksize);
206 RCP<CrsMatrix> blockCrs_as_crs = rcp(
new TpetraBlockCrsMatrix(blockCrs));
207 RCP<CrsMatrixWrap> blockWrap = rcp(
new CrsMatrixWrap(blockCrs_as_crs));
209 this->GetOStream(
Statistics0) <<
"Ifpack2Smoother: Using BlockCrsMatrix storage with blocksize " << blocksize << std::endl;
211 paramList.remove(
"smoother: use blockcrsmatrix storage");
216 if (type_ ==
"SCHWARZ")
217 SetupSchwarz(currentLevel);
219 else if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
220 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
221 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
222 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
223 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
224 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION" ||
225 type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
226 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
227 type_ ==
"LINESMOOTHING_BANDEDRELAXATION" ||
228 type_ ==
"LINESMOOTHING_BLOCK_RELAXATION" ||
229 type_ ==
"LINESMOOTHING_BLOCK RELAXATION" ||
230 type_ ==
"LINESMOOTHING_BLOCKRELAXATION")
231 SetupLineSmoothing(currentLevel);
233 else if (type_ ==
"BLOCK_RELAXATION" ||
234 type_ ==
"BLOCK RELAXATION" ||
235 type_ ==
"BLOCKRELAXATION" ||
237 type_ ==
"BANDED_RELAXATION" ||
238 type_ ==
"BANDED RELAXATION" ||
239 type_ ==
"BANDEDRELAXATION" ||
241 type_ ==
"TRIDI_RELAXATION" ||
242 type_ ==
"TRIDI RELAXATION" ||
243 type_ ==
"TRIDIRELAXATION" ||
244 type_ ==
"TRIDIAGONAL_RELAXATION" ||
245 type_ ==
"TRIDIAGONAL RELAXATION" ||
246 type_ ==
"TRIDIAGONALRELAXATION")
247 SetupBlockRelaxation(currentLevel);
249 else if (type_ ==
"CHEBYSHEV")
250 SetupChebyshev(currentLevel);
252 else if (type_ ==
"TOPOLOGICAL") {
253#ifdef HAVE_MUELU_INTREPID2
254 SetupTopological(currentLevel);
256 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"'TOPOLOGICAL' smoother choice requires Intrepid2");
258 }
else if (type_ ==
"AGGREGATE")
259 SetupAggregate(currentLevel);
261 else if (type_ ==
"HIPTMAIR")
262 SetupHiptmair(currentLevel);
265 SetupGeneric(currentLevel);
270 this->GetOStream(
Statistics1) << description() << std::endl;
273template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
275 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
277 bool reusePreconditioner =
false;
278 if (this->IsSetup() ==
true) {
280 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
282 bool isTRowMatrix =
true;
283 RCP<const tRowMatrix> tA;
285 tA = Xpetra::toTpetraRowMatrix(A_);
287 isTRowMatrix =
false;
290 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
291 if (!prec.is_null() && isTRowMatrix) {
293 reusePreconditioner =
true;
295 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
296 "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction"
301 if (!reusePreconditioner) {
302 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
304 bool isBlockedMatrix =
false;
305 RCP<Matrix> merged2Mat;
307 std::string sublistName =
"subdomain solver parameters";
308 if (paramList.isSublist(sublistName)) {
317 ParameterList& subList = paramList.sublist(sublistName);
319 std::string partName =
"partitioner: type";
325 if (subList.isParameter(partName) && subList.get<std::string>(partName) ==
"vanka user") {
326 isBlockedMatrix =
true;
328 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
330 "Matrix A must be of type BlockedCrsMatrix.");
332 size_t numVels = bA->getMatrix(0, 0)->getLocalNumRows();
333 size_t numPres = bA->getMatrix(1, 0)->getLocalNumRows();
334 size_t numRows = rcp_dynamic_cast<Matrix>(A_,
true)->getLocalNumRows();
336 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
338 size_t numBlocks = 0;
339 for (
size_t rowOfB = numVels; rowOfB < numVels + numPres; ++rowOfB)
340 blockSeeds[rowOfB] = numBlocks++;
342 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
344 "Matrix A must be of type BlockedCrsMatrix.");
346 merged2Mat = bA2->Merge();
350 bool haveBoundary =
false;
351 for (LO i = 0; i < boundaryNodes.size(); i++)
352 if (boundaryNodes[i]) {
356 blockSeeds[i] = numBlocks;
362 subList.set(
"partitioner: type",
"user");
363 subList.set(
"partitioner: map", blockSeeds);
364 subList.set(
"partitioner: local parts", as<int>(numBlocks));
367 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
369 isBlockedMatrix =
true;
370 merged2Mat = bA->Merge();
375 RCP<const tRowMatrix> tA;
376 if (isBlockedMatrix ==
true)
377 tA = Xpetra::toTpetraRowMatrix(merged2Mat);
379 tA = Xpetra::toTpetraRowMatrix(A_);
390template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
392 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
394 if (this->IsSetup() ==
true) {
395 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupAggregate(): Setup() has already been called" << std::endl;
396 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupAggregate(): reuse of this type is not available, reverting to full construction" << std::endl;
399 this->GetOStream(
Statistics0) <<
"Ifpack2Smoother: Using Aggregate Smoothing" << std::endl;
401 RCP<Aggregates> aggregates = Factory::Get<RCP<Aggregates>>(currentLevel,
"Aggregates");
403 RCP<const LOMultiVector> vertex2AggId = aggregates->GetVertex2AggId();
404 ArrayRCP<LO> aggregate_ids = rcp_const_cast<LOMultiVector>(vertex2AggId)->getDataNonConst(0);
405 ArrayRCP<LO> dof_ids;
409 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
411 blocksize = matA->GetFixedBlockSize();
413 dof_ids.resize(aggregate_ids.size() * blocksize);
414 for (LO i = 0; i < (LO)aggregate_ids.size(); i++) {
415 for (LO j = 0; j < (LO)blocksize; j++)
416 dof_ids[i * blocksize + j] = aggregate_ids[i];
419 dof_ids = aggregate_ids;
422 paramList.set(
"partitioner: map", dof_ids);
423 paramList.set(
"partitioner: type",
"user");
424 paramList.set(
"partitioner: overlap", 0);
425 paramList.set(
"partitioner: local parts", (
int)aggregates->GetNumAggregates());
427 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Xpetra::toTpetraRowMatrix(A_);
429 type_ =
"BLOCKRELAXATION";
436#ifdef HAVE_MUELU_INTREPID2
437template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
447 if (this->IsSetup() ==
true) {
448 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
449 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
452 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
454 typedef typename Node::device_type::execution_space ES;
456 typedef Kokkos::DynRankView<LocalOrdinal, typename Node::device_type> FCO;
458 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
462 const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO>>(currentLevel,
"pcoarsen: element to node map");
464 string basisString = paramList.get<
string>(
"pcoarsen: hi basis");
470 auto basis = MueLuIntrepid::BasisFactory<double, ES>(basisString, degree);
472 string topologyTypeString = paramList.get<
string>(
"smoother: neighborhood type");
474 if (topologyTypeString ==
"node")
476 else if (topologyTypeString ==
"edge")
478 else if (topologyTypeString ==
"face")
480 else if (topologyTypeString ==
"cell")
481 dimension = basis->getBaseCellTopology().getDimension();
483 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
484 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_,
true);
485 vector<vector<LocalOrdinal>> seeds;
490 int myNodeCount = matA->getRowMap()->getLocalNumElements();
491 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount, lo_invalid);
492 int localPartitionNumber = 0;
494 nodeSeeds[seed] = localPartitionNumber++;
497 paramList.remove(
"smoother: neighborhood type");
498 paramList.remove(
"pcoarsen: hi basis");
500 paramList.set(
"partitioner: map", nodeSeeds);
501 paramList.set(
"partitioner: type",
"user");
502 paramList.set(
"partitioner: overlap", 1);
503 paramList.set(
"partitioner: local parts",
int(seeds[dimension].size()));
505 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Xpetra::toTpetraRowMatrix(A_);
507 type_ =
"BLOCKRELAXATION";
515template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
517 if (this->IsSetup() ==
true) {
518 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
519 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
522 ParameterList& myparamList =
const_cast<ParameterList&
>(this->GetParameterList());
524 LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,
"CoarseNumZLayers");
525 if (CoarseNumZLayers > 0) {
526 Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get<Teuchos::ArrayRCP<LO>>(currentLevel,
"LineDetection_VertLineIds");
530 for (
size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
531 if (maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
533 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_,
true);
534 size_t numLocalRows = matA->getLocalNumRows();
537 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
542 int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
543 LO matrixBlockSize = matA->GetFixedBlockSize();
544 if (matrixBlockSize > 1 && actualDofsPerNode > 1) {
546 "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
547 }
else if (matrixBlockSize > 1) {
548 actualDofsPerNode = matrixBlockSize;
550 myparamList.set(
"partitioner: PDE equations", actualDofsPerNode);
552 if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
553 myparamList.set(
"partitioner: type",
"user");
554 myparamList.set(
"partitioner: map", TVertLineIdSmoo);
555 myparamList.set(
"partitioner: local parts", maxPart + 1);
557 if (myparamList.isParameter(
"partitioner: block size") &&
558 myparamList.get<
int>(
"partitioner: block size") != -1) {
559 int block_size = myparamList.get<
int>(
"partitioner: block size");
561 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the specified block size.");
562 numLocalRows /= block_size;
566 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
569 Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
570 for (
size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
571 for (
size_t dof = 0; dof < numDofsPerNode; dof++)
572 partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
573 myparamList.set(
"partitioner: type",
"user");
574 myparamList.set(
"partitioner: map", partitionerMap);
575 myparamList.set(
"partitioner: local parts", maxPart + 1);
578 if (type_ ==
"LINESMOOTHING_BANDED_RELAXATION" ||
579 type_ ==
"LINESMOOTHING_BANDED RELAXATION" ||
580 type_ ==
"LINESMOOTHING_BANDEDRELAXATION")
581 type_ =
"BANDEDRELAXATION";
582 else if (type_ ==
"LINESMOOTHING_TRIDI_RELAXATION" ||
583 type_ ==
"LINESMOOTHING_TRIDI RELAXATION" ||
584 type_ ==
"LINESMOOTHING_TRIDIRELAXATION" ||
585 type_ ==
"LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
586 type_ ==
"LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
587 type_ ==
"LINESMOOTHING_TRIDIAGONALRELAXATION")
588 type_ =
"TRIDIAGONALRELAXATION";
590 type_ =
"BLOCKRELAXATION";
593 this->GetOStream(
Runtime0) <<
"Line detection failed: fall back to point-wise relaxation" << std::endl;
594 myparamList.remove(
"partitioner: type",
false);
595 myparamList.remove(
"partitioner: map",
false);
596 myparamList.remove(
"partitioner: local parts",
false);
597 type_ =
"RELAXATION";
600 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Xpetra::toTpetraRowMatrix(A_);
608template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
610 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
612 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
616 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
618 bool reusePreconditioner =
false;
619 if (this->IsSetup() ==
true) {
621 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
623 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
624 if (!prec.is_null()) {
626 reusePreconditioner =
true;
628 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
629 "reverting to full construction"
634 if (!reusePreconditioner) {
635 ParameterList& myparamList =
const_cast<ParameterList&
>(this->GetParameterList());
637 if (myparamList.isParameter(
"partitioner: type") &&
638 myparamList.get<std::string>(
"partitioner: type") ==
"line") {
639 Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>> xCoordinates =
640 Factory::Get<Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>>>(currentLevel,
"Coordinates");
641 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));
643 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
644 size_t lclSize = A_->getRangeMap()->getLocalNumElements();
646 lclSize = matA->getLocalNumRows();
647 size_t numDofsPerNode = lclSize / xCoordinates->getMap()->getLocalNumElements();
648 myparamList.set(
"partitioner: coordinates", coordinates);
649 myparamList.set(
"partitioner: PDE equations", (
int)numDofsPerNode);
660template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
662 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
663 using STS = Teuchos::ScalarTraits<SC>;
664 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
668 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
670 bool reusePreconditioner =
false;
672 if (this->IsSetup() ==
true) {
674 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): Setup() has already been called, assuming reuse" << std::endl;
676 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
677 if (!prec.is_null()) {
679 reusePreconditioner =
true;
681 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available (failed cast to CanChangeMatrix), "
682 "reverting to full construction"
688 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
689 SC negone = -STS::one();
690 SC lambdaMax = SetupChebyshevEigenvalues(currentLevel,
"A",
"", paramList);
692 if (!reusePreconditioner) {
707 if (lambdaMax == negone) {
708 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
710 Teuchos::RCP<Ifpack2::Chebyshev<MatrixType>> chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType>>(prec_);
711 if (chebyPrec != Teuchos::null) {
712 lambdaMax = chebyPrec->getLambdaMaxForApply();
713 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
715 matA->SetMaxEigenvalueEstimate(lambdaMax);
716 this->GetOStream(
Statistics1) <<
"chebyshev: max eigenvalue (calculated by Ifpack2)"
717 <<
" = " << lambdaMax << std::endl;
719 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone,
Exceptions::RuntimeError,
"MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
723template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
725 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
726 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
730 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
732 bool reusePreconditioner =
false;
733 if (this->IsSetup() ==
true) {
735 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupHiptmair(): Setup() has already been called, assuming reuse" << std::endl;
737 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
738 if (!prec.is_null()) {
740 reusePreconditioner =
true;
742 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupHiptmair(): reuse of this type is not available (failed cast to CanChangeMatrix), "
743 "reverting to full construction"
749 SC negone = -Teuchos::ScalarTraits<Scalar>::one();
750 ParameterList& paramList =
const_cast<ParameterList&
>(this->GetParameterList());
751 std::string smoother1 = paramList.get(
"hiptmair: smoother type 1",
"CHEBYSHEV");
752 std::string smoother2 = paramList.get(
"hiptmair: smoother type 2",
"CHEBYSHEV");
753 SC lambdaMax11 = negone;
755 if (smoother1 ==
"CHEBYSHEV") {
756 ParameterList& list1 = paramList.sublist(
"hiptmair: smoother list 1");
757 lambdaMax11 = SetupChebyshevEigenvalues(currentLevel,
"A",
"EdgeMatrix ", list1);
759 if (smoother2 ==
"CHEBYSHEV") {
760 ParameterList& list2 = paramList.sublist(
"hiptmair: smoother list 2");
761 SetupChebyshevEigenvalues(currentLevel,
"NodeMatrix",
"NodeMatrix ", list2);
768 RCP<Operator> NodeMatrix = currentLevel.
Get<RCP<Operator>>(
"NodeMatrix");
769 RCP<Operator> D0 = currentLevel.
Get<RCP<Operator>>(
"D0");
771 RCP<tRowMatrix> tNodeMatrix = Xpetra::toTpetraRowMatrix(NodeMatrix);
772 RCP<tRowMatrix> tD0 = Xpetra::toTpetraRowMatrix(D0);
774 Teuchos::ParameterList newlist;
775 newlist.set(
"P", tD0);
776 newlist.set(
"PtAP", tNodeMatrix);
777 if (!reusePreconditioner) {
779 SetPrecParameters(newlist);
786 if (smoother1 ==
"CHEBYSHEV" && lambdaMax11 == negone) {
787 using Teuchos::rcp_dynamic_cast;
788 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
789 auto hiptmairPrec = rcp_dynamic_cast<Ifpack2::Hiptmair<MatrixType>>(prec_);
790 if (hiptmairPrec != Teuchos::null) {
791 auto chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType>>(hiptmairPrec->getPrec1());
792 if (chebyPrec != Teuchos::null) {
793 lambdaMax11 = chebyPrec->getLambdaMaxForApply();
794 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
796 matA->SetMaxEigenvalueEstimate(lambdaMax11);
797 this->GetOStream(
Statistics1) <<
"chebyshev: max eigenvalue (calculated by Ifpack2)"
798 <<
" = " << lambdaMax11 << std::endl;
800 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax11 == negone,
Exceptions::RuntimeError,
"MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
805template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
808 typedef Teuchos::ScalarTraits<SC> STS;
809 SC negone = -STS::one();
810 RCP<Operator> currentA = currentLevel.
Get<RCP<Operator>>(matrixName);
811 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(currentA);
812 SC lambdaMax = negone;
814 std::string maxEigString =
"chebyshev: max eigenvalue";
815 std::string eigRatioString =
"chebyshev: ratio eigenvalue";
818 if (paramList.isParameter(maxEigString)) {
819 if (paramList.isType<
double>(maxEigString))
820 lambdaMax = paramList.get<
double>(maxEigString);
822 lambdaMax = paramList.get<SC>(maxEigString);
823 this->GetOStream(
Statistics1) << label << maxEigString <<
" (cached with smoother parameter list) = " << lambdaMax << std::endl;
825 matA->SetMaxEigenvalueEstimate(lambdaMax);
827 }
else if (!matA.is_null()) {
828 lambdaMax = matA->GetMaxEigenvalueEstimate();
829 if (lambdaMax != negone) {
830 this->GetOStream(
Statistics1) << label << maxEigString <<
" (cached with matrix) = " << lambdaMax << std::endl;
831 paramList.set(maxEigString, lambdaMax);
836 const SC defaultEigRatio = 20;
838 SC ratio = defaultEigRatio;
839 if (paramList.isParameter(eigRatioString)) {
840 if (paramList.isType<
double>(eigRatioString))
841 ratio = paramList.get<
double>(eigRatioString);
843 ratio = paramList.get<SC>(eigRatioString);
850 RCP<const Operator> fineA = currentLevel.
GetPreviousLevel()->Get<RCP<Operator>>(matrixName);
851 size_t nRowsFine = fineA->getDomainMap()->getGlobalNumElements();
852 size_t nRowsCoarse = currentA->getDomainMap()->getGlobalNumElements();
854 SC levelRatio = as<SC>(as<double>(nRowsFine) / nRowsCoarse);
855 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
859 this->GetOStream(
Statistics1) << label << eigRatioString <<
" (computed) = " << ratio << std::endl;
860 paramList.set(eigRatioString, ratio);
862 if (paramList.isParameter(
"chebyshev: use rowsumabs diagonal scaling")) {
863 this->GetOStream(
Runtime1) <<
"chebyshev: using rowsumabs diagonal scaling" << std::endl;
864 bool doScale =
false;
865 doScale = paramList.get<
bool>(
"chebyshev: use rowsumabs diagonal scaling");
866 paramList.remove(
"chebyshev: use rowsumabs diagonal scaling");
867 double chebyReplaceTol = Teuchos::ScalarTraits<Scalar>::eps() * 100;
868 std::string paramName =
"chebyshev: rowsumabs diagonal replacement tolerance";
869 if (paramList.isParameter(paramName)) {
870 chebyReplaceTol = paramList.get<
double>(paramName);
871 paramList.remove(paramName);
873 double chebyReplaceVal = Teuchos::ScalarTraits<double>::zero();
874 paramName =
"chebyshev: rowsumabs diagonal replacement value";
875 if (paramList.isParameter(paramName)) {
876 chebyReplaceVal = paramList.get<
double>(paramName);
877 paramList.remove(paramName);
879 bool chebyReplaceSingleEntryRowWithZero =
false;
880 paramName =
"chebyshev: rowsumabs replace single entry row with zero";
881 if (paramList.isParameter(paramName)) {
882 chebyReplaceSingleEntryRowWithZero = paramList.get<
bool>(paramName);
883 paramList.remove(paramName);
885 bool useAverageAbsDiagVal =
false;
886 paramName =
"chebyshev: rowsumabs use automatic diagonal tolerance";
887 if (paramList.isParameter(paramName)) {
888 useAverageAbsDiagVal = paramList.get<
bool>(paramName);
889 paramList.remove(paramName);
892 TEUCHOS_ASSERT(!matA.is_null());
893 const bool doReciprocal =
true;
894 RCP<Vector> lumpedDiagonal =
Utilities::GetLumpedMatrixDiagonal(*matA, doReciprocal, chebyReplaceTol, chebyReplaceVal, chebyReplaceSingleEntryRowWithZero, useAverageAbsDiagVal);
895 const Xpetra::TpetraVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tmpVec =
dynamic_cast<const Xpetra::TpetraVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>&
>(*lumpedDiagonal);
896 paramList.set(
"chebyshev: operator inv diagonal", tmpVec.getTpetra_Vector());
903template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
905 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
906 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
910 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
912 bool reusePreconditioner =
false;
913 if (this->IsSetup() ==
true) {
915 this->GetOStream(
Runtime1) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
917 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
918 if (!prec.is_null()) {
920 reusePreconditioner =
true;
922 this->GetOStream(
Warnings0) <<
"MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
923 "reverting to full construction"
928 if (!reusePreconditioner) {
937template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
947 bool supportInitialGuess =
false;
949 if (prec_->supportsZeroStartingSolution()) {
950 prec_->setZeroStartingSolution(InitialGuessIsZero);
951 supportInitialGuess =
true;
952 }
else if (type_ ==
"SCHWARZ") {
953 Teuchos::ParameterList paramList;
954 paramList.set(
"schwarz: zero starting solution", InitialGuessIsZero);
959 prec_->setParameters(paramList);
960 supportInitialGuess =
true;
970 if (InitialGuessIsZero || supportInitialGuess) {
971 Tpetra::MultiVector<SC, LO, GO, NO>& tpX = toTpetra(X);
972 const Tpetra::MultiVector<SC, LO, GO, NO>& tpB = toTpetra(B);
973 prec_->apply(tpB, tpX);
975 using TST = Teuchos::ScalarTraits<Scalar>;
977 RCP<MultiVector> Residual;
979 std::string prefix = this->ShortClassName() +
": Apply: ";
980 RCP<TimeMonitor> tM = rcp(
new TimeMonitor(*
this, prefix +
"residual calculation",
Timings0));
984 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
986 Tpetra::MultiVector<SC, LO, GO, NO>& tpX = toTpetra(*Correction);
987 const Tpetra::MultiVector<SC, LO, GO, NO>& tpB = toTpetra(*Residual);
989 prec_->apply(tpB, tpX);
991 X.update(TST::one(), *Correction, TST::one());
995template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
998 smoother->SetParameterList(this->GetParameterList());
1002template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1004 std::ostringstream out;
1006 out << prec_->description();
1009 out <<
"{type = " << type_ <<
"}";
1014template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1019 out0 <<
"Prec. type: " << type_ << std::endl;
1022 out0 <<
"Parameter list: " << std::endl;
1023 Teuchos::OSTab tab2(out);
1024 out << this->GetParameterList();
1025 out0 <<
"Overlap: " << overlap_ << std::endl;
1029 if (prec_ != Teuchos::null) {
1030 Teuchos::OSTab tab2(out);
1031 out << *prec_ << std::endl;
1034 if (verbLevel &
Debug) {
1037 <<
"RCP<prec_>: " << prec_ << std::endl;
1041template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1043 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
1045 RCP<Ifpack2::Relaxation<MatrixType>> pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType>>(prec_);
1046 if (!pr.is_null())
return pr->getNodeSmootherComplexity();
1048 RCP<Ifpack2::Chebyshev<MatrixType>> pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType>>(prec_);
1049 if (!pc.is_null())
return pc->getNodeSmootherComplexity();
1051 RCP<Ifpack2::BlockRelaxation<MatrixType>> pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType>>(prec_);
1052 if (!pb.is_null())
return pb->getNodeSmootherComplexity();
1054 RCP<Ifpack2::ILUT<MatrixType>> pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType>>(prec_);
1055 if (!pi.is_null())
return pi->getNodeSmootherComplexity();
1057 RCP<Ifpack2::RILUK<MatrixType>> pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType>>(prec_);
1058 if (!pk.is_null())
return pk->getNodeSmootherComplexity();
1060 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 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.
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)
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.
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 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, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Magnitude >::zero(), bool count_twos_as_dirichlet=false)
Detect Dirichlet rows.
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)
Extract Matrix Diagonal of lumped matrix.
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
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.
@ 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)