MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Ifpack2Smoother_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_IFPACK2SMOOTHER_DEF_HPP
11#define MUELU_IFPACK2SMOOTHER_DEF_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14#include "Xpetra_TpetraRowMatrix.hpp"
15
16#include <Teuchos_ParameterList.hpp>
17
18#include <Tpetra_RowMatrix.hpp>
19
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>
28
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>
37
38#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
39
41#include "MueLu_Level.hpp"
42#include "MueLu_Utilities.hpp"
43#include "MueLu_Monitor.hpp"
44#include "MueLu_Aggregates.hpp"
45
46#ifdef HAVE_MUELU_INTREPID2
49#include "Intrepid2_Basis.hpp"
50#include "Kokkos_DynRankView.hpp"
51#endif
52
53namespace MueLu {
54
55template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Ifpack2Smoother(const std::string& type, const Teuchos::ParameterList& paramList, const LO& overlap)
57 : type_(type)
58 , overlap_(overlap) {
59 if (!type_.empty()) {
60 // Transform string to "ABCDE" notation
61 std::transform(type_.begin(), type_.end(), type_.begin(), ::toupper);
62 }
63
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");
79 this->declareConstructionOutcome(!isSupported, "Ifpack2 does not provide the smoother '" + type_ + "'.");
80 if (isSupported)
81 SetParameterList(paramList);
82}
83
84template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87
89 // It might be invalid to change parameters after the setup, but it depends entirely on Ifpack implementation.
90 // TODO: I don't know if Ifpack returns an error code or exception or ignore parameters modification in this case...
91 SetPrecParameters();
92 }
93}
94
95template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
97 std::string prefix = this->ShortClassName() + ": SetPrecParameters";
98 RCP<TimeMonitor> tM = rcp(new TimeMonitor(*this, prefix, Timings0));
99 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
100
101 paramList.setParameters(list);
102
103 RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
104
105 // Do we want an Ifpack2 apply timer?
106 precList->set("timer for apply", this->IsPrint(Timings));
107
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();
116 }
117 precList->set("partitioner: local parts", lclSize / matrixBlockSize);
118 }
119
120 prec_->setParameters(*precList);
121
122 paramList.setParameters(*precList); // what about that??
123}
124
125template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
127 this->Input(currentLevel, "A");
128
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"); // necessary for fallback criterion
142 this->Input(currentLevel, "LineDetection_VertLineIds"); // necessary to feed block smoother
143 } else if (type_ == "BLOCK RELAXATION" ||
144 type_ == "BLOCK_RELAXATION" ||
145 type_ == "BLOCKRELAXATION" ||
146 // Banded
147 type_ == "BANDED_RELAXATION" ||
148 type_ == "BANDED RELAXATION" ||
149 type_ == "BANDEDRELAXATION" ||
150 // Tridiagonal
151 type_ == "TRIDI_RELAXATION" ||
152 type_ == "TRIDI RELAXATION" ||
153 type_ == "TRIDIRELAXATION" ||
154 type_ == "TRIDIAGONAL_RELAXATION" ||
155 type_ == "TRIDIAGONAL RELAXATION" ||
156 type_ == "TRIDIAGONALRELAXATION") {
157 // We need to check for the "partitioner type" = "line"
158 ParameterList precList = this->GetParameterList();
159 if (precList.isParameter("partitioner: type") &&
160 precList.get<std::string>("partitioner: type") == "line") {
161 this->Input(currentLevel, "Coordinates");
162 }
163 } else if (type_ == "TOPOLOGICAL") {
164 // for the topological smoother, we require an element to node map:
165 this->Input(currentLevel, "pcoarsen: element to node map");
166 } else if (type_ == "AGGREGATE") {
167 // Aggregate smoothing needs aggregates
168 this->Input(currentLevel, "Aggregates");
169 } else if (type_ == "HIPTMAIR") {
170 // Hiptmair needs D0 and NodeMatrix
171 this->Input(currentLevel, "NodeMatrix");
172 this->Input(currentLevel, "D0");
173 }
174}
175
176template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
178 FactoryMonitor m(*this, "Setup Smoother", currentLevel);
179 A_ = Factory::Get<RCP<Operator>>(currentLevel, "A");
180 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
181
182 // If the user asked us to convert the matrix into BlockCrsMatrix form, we do that now
183
184 if (paramList.isParameter("smoother: use blockcrsmatrix storage") && paramList.get<bool>("smoother: use blockcrsmatrix storage")) {
185 int blocksize = 1;
186 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
187 if (!matA.is_null())
188 blocksize = matA->GetFixedBlockSize();
189 if (blocksize) {
190 // NOTE: Don't think you can move this out of the if block. You can't. The test MueLu_MeshTyingBlocked_SimpleSmoother_2dof_medium_MPI_1 will fail
191
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();
196 if (Acrs.is_null())
197 throw std::runtime_error("Ifpack2Smoother: Cannot extract CrsMatrix from matrix A.");
198 RCP<TpetraCrsMatrix> At = rcp_dynamic_cast<TpetraCrsMatrix>(Acrs);
199 if (At.is_null()) {
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");
204 } else {
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));
208 A_ = blockWrap;
209 this->GetOStream(Statistics0) << "Ifpack2Smoother: Using BlockCrsMatrix storage with blocksize " << blocksize << std::endl;
210
211 paramList.remove("smoother: use blockcrsmatrix storage");
212 }
213 }
214 }
215
216 if (type_ == "SCHWARZ")
217 SetupSchwarz(currentLevel);
218
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);
232
233 else if (type_ == "BLOCK_RELAXATION" ||
234 type_ == "BLOCK RELAXATION" ||
235 type_ == "BLOCKRELAXATION" ||
236 // Banded
237 type_ == "BANDED_RELAXATION" ||
238 type_ == "BANDED RELAXATION" ||
239 type_ == "BANDEDRELAXATION" ||
240 // Tridiagonal
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);
248
249 else if (type_ == "CHEBYSHEV")
250 SetupChebyshev(currentLevel);
251
252 else if (type_ == "TOPOLOGICAL") {
253#ifdef HAVE_MUELU_INTREPID2
254 SetupTopological(currentLevel);
255#else
256 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, "'TOPOLOGICAL' smoother choice requires Intrepid2");
257#endif
258 } else if (type_ == "AGGREGATE")
259 SetupAggregate(currentLevel);
260
261 else if (type_ == "HIPTMAIR")
262 SetupHiptmair(currentLevel);
263
264 else {
265 SetupGeneric(currentLevel);
266 }
267
269
270 this->GetOStream(Statistics1) << description() << std::endl;
271}
272
273template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
275 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
276
277 bool reusePreconditioner = false;
278 if (this->IsSetup() == true) {
279 // Reuse the constructed preconditioner
280 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
281
282 bool isTRowMatrix = true;
283 RCP<const tRowMatrix> tA;
284 try {
285 tA = Xpetra::toTpetraRowMatrix(A_);
286 } catch (Exceptions::BadCast&) {
287 isTRowMatrix = false;
288 }
289
290 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
291 if (!prec.is_null() && isTRowMatrix) {
292 prec->setMatrix(tA);
293 reusePreconditioner = true;
294 } else {
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"
297 << std::endl;
298 }
299 }
300
301 if (!reusePreconditioner) {
302 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
303
304 bool isBlockedMatrix = false;
305 RCP<Matrix> merged2Mat;
306
307 std::string sublistName = "subdomain solver parameters";
308 if (paramList.isSublist(sublistName)) {
309 // If we are doing "user" partitioning, we assume that what the user
310 // really wants to do is make tiny little subdomains with one row
311 // assigned to each subdomain. The rows used for these little
312 // subdomains correspond to those in the 2nd block row. Then,
313 // if we overlap these mini-subdomains, we will do something that
314 // looks like Vanka (grabbing all velocities associated with each
315 // each pressure unknown). In addition, we put all Dirichlet points
316 // as a little mini-domain.
317 ParameterList& subList = paramList.sublist(sublistName);
318
319 std::string partName = "partitioner: type";
320 // Pretty sure no one has been using this. Unfortunately, old if
321 // statement (which checked for equality with "user") prevented
322 // anyone from employing other types of Ifpack2 user partition
323 // options. Leaving this and switching if to "vanka user" just
324 // in case some day someone might want to use this.
325 if (subList.isParameter(partName) && subList.get<std::string>(partName) == "vanka user") {
326 isBlockedMatrix = true;
327
328 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
329 TEUCHOS_TEST_FOR_EXCEPTION(bA.is_null(), Exceptions::BadCast,
330 "Matrix A must be of type BlockedCrsMatrix.");
331
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();
335
336 ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
337
338 size_t numBlocks = 0;
339 for (size_t rowOfB = numVels; rowOfB < numVels + numPres; ++rowOfB)
340 blockSeeds[rowOfB] = numBlocks++;
341
342 RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
343 TEUCHOS_TEST_FOR_EXCEPTION(bA2.is_null(), Exceptions::BadCast,
344 "Matrix A must be of type BlockedCrsMatrix.");
345
346 merged2Mat = bA2->Merge();
347
348 // Add Dirichlet rows to the list of seeds
349 ArrayRCP<const bool> boundaryNodes = Utilities::DetectDirichletRows(*merged2Mat, 0.0);
350 bool haveBoundary = false;
351 for (LO i = 0; i < boundaryNodes.size(); i++)
352 if (boundaryNodes[i]) {
353 // FIXME:
354 // 1. would not this [] overlap with some in the previos blockSeed loop?
355 // 2. do we need to distinguish between pressure and velocity Dirichlet b.c.
356 blockSeeds[i] = numBlocks;
357 haveBoundary = true;
358 }
359 if (haveBoundary)
360 numBlocks++;
361
362 subList.set("partitioner: type", "user");
363 subList.set("partitioner: map", blockSeeds);
364 subList.set("partitioner: local parts", as<int>(numBlocks));
365
366 } else {
367 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
368 if (!bA.is_null()) {
369 isBlockedMatrix = true;
370 merged2Mat = bA->Merge();
371 }
372 }
373 }
374
375 RCP<const tRowMatrix> tA;
376 if (isBlockedMatrix == true)
377 tA = Xpetra::toTpetraRowMatrix(merged2Mat);
378 else
379 tA = Xpetra::toTpetraRowMatrix(A_);
380
381 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
382 SetPrecParameters();
383
384 prec_->initialize();
385 }
386
387 prec_->compute();
388}
389
390template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
392 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
393
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;
397 }
398
399 this->GetOStream(Statistics0) << "Ifpack2Smoother: Using Aggregate Smoothing" << std::endl;
400
401 RCP<Aggregates> aggregates = Factory::Get<RCP<Aggregates>>(currentLevel, "Aggregates");
402
403 RCP<const LOMultiVector> vertex2AggId = aggregates->GetVertex2AggId();
404 ArrayRCP<LO> aggregate_ids = rcp_const_cast<LOMultiVector>(vertex2AggId)->getDataNonConst(0);
405 ArrayRCP<LO> dof_ids;
406
407 // We need to unamalgamate, if the FixedBlockSize > 1
408 LO blocksize = 1;
409 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
410 if (!matA.is_null())
411 blocksize = matA->GetFixedBlockSize();
412 if (blocksize > 1) {
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];
417 }
418 } else {
419 dof_ids = aggregate_ids;
420 }
421
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());
426
427 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Xpetra::toTpetraRowMatrix(A_);
428
429 type_ = "BLOCKRELAXATION";
430 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
431 SetPrecParameters();
432 prec_->initialize();
433 prec_->compute();
434}
435
436#ifdef HAVE_MUELU_INTREPID2
437template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
439 /*
440
441 basic notion:
442
443 Look for user input indicating topo dimension, something like "topological domain type: {node|edge|face}"
444 Call something like what you can find in Poisson example line 1180 to set seeds for a smoother
445
446 */
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;
450 }
451
452 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
453
454 typedef typename Node::device_type::execution_space ES;
455
456 typedef Kokkos::DynRankView<LocalOrdinal, typename Node::device_type> FCO; //
457
458 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
459
460 using namespace std;
461
462 const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO>>(currentLevel, "pcoarsen: element to node map");
463
464 string basisString = paramList.get<string>("pcoarsen: hi basis");
465 int degree;
466 // NOTE: To make sure Stokhos works we only instantiate these guys with double. There's a lot
467 // of stuff in the guts of Intrepid2 that doesn't play well with Stokhos as of yet. Here, we only
468 // care about the assignment of basis ordinals to topological entities, so this code is actually
469 // independent of the Scalar type--hard-coding double here won't hurt us.
470 auto basis = MueLuIntrepid::BasisFactory<double, ES>(basisString, degree);
471
472 string topologyTypeString = paramList.get<string>("smoother: neighborhood type");
473 int dimension;
474 if (topologyTypeString == "node")
475 dimension = 0;
476 else if (topologyTypeString == "edge")
477 dimension = 1;
478 else if (topologyTypeString == "face")
479 dimension = 2;
480 else if (topologyTypeString == "cell")
481 dimension = basis->getBaseCellTopology().getDimension();
482 else
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;
486 MueLuIntrepid::FindGeometricSeedOrdinals(basis, *elemToNode, seeds, *matA->getRowMap(), *matA->getColMap());
487
488 // Ifpack2 wants the seeds in an array of the same length as the number of local elements,
489 // with local partition #s marked for the ones that are seeds, and invalid for the rest
490 int myNodeCount = matA->getRowMap()->getLocalNumElements();
491 ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount, lo_invalid);
492 int localPartitionNumber = 0;
493 for (LocalOrdinal seed : seeds[dimension]) {
494 nodeSeeds[seed] = localPartitionNumber++;
495 }
496
497 paramList.remove("smoother: neighborhood type");
498 paramList.remove("pcoarsen: hi basis");
499
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()));
504
505 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Xpetra::toTpetraRowMatrix(A_);
506
507 type_ = "BLOCKRELAXATION";
508 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
509 SetPrecParameters();
510 prec_->initialize();
511 prec_->compute();
512}
513#endif
514
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;
520 }
521
522 ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
523
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");
527
528 // determine number of local parts
529 LO maxPart = 0;
530 for (size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
531 if (maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
532 }
533 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_, true);
534 size_t numLocalRows = matA->getLocalNumRows();
535
536 TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % TVertLineIdSmoo.size() != 0, Exceptions::RuntimeError,
537 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
538
539 // actualDofsPerNode is the actual number of matrix rows per mesh element.
540 // It is encoded in either the MueLu Level, or in the Xpetra matrix block size.
541 // This value is needed by Ifpack2 to do decoupled block relaxation.
542 int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
543 LO matrixBlockSize = matA->GetFixedBlockSize();
544 if (matrixBlockSize > 1 && actualDofsPerNode > 1) {
545 TEUCHOS_TEST_FOR_EXCEPTION(actualDofsPerNode != matrixBlockSize, Exceptions::RuntimeError,
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;
549 }
550 myparamList.set("partitioner: PDE equations", actualDofsPerNode);
551
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);
556 } else {
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");
560 TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % block_size != 0, Exceptions::RuntimeError,
561 "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the specified block size.");
562 numLocalRows /= block_size;
563 }
564
565 // we assume a constant number of DOFs per node
566 size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
567
568 // Create a new Teuchos::ArrayRCP<LO> of size numLocalRows and fill it with the corresponding information
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);
576 }
577
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";
589 else
590 type_ = "BLOCKRELAXATION";
591 } else {
592 // line detection failed -> fallback to point-wise relaxation
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";
598 }
599
600 RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Xpetra::toTpetraRowMatrix(A_);
601
602 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
603 SetPrecParameters();
604 prec_->initialize();
605 prec_->compute();
606}
607
608template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
610 typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
611
612 RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
613 if (!bA.is_null())
614 A_ = bA->Merge();
615
616 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
617
618 bool reusePreconditioner = false;
619 if (this->IsSetup() == true) {
620 // Reuse the constructed preconditioner
621 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
622
623 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
624 if (!prec.is_null()) {
625 prec->setMatrix(tA);
626 reusePreconditioner = true;
627 } else {
628 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
629 "reverting to full construction"
630 << std::endl;
631 }
632 }
633
634 if (!reusePreconditioner) {
635 ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
636 myparamList.print();
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));
642
643 RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
644 size_t lclSize = A_->getRangeMap()->getLocalNumElements();
645 if (!matA.is_null())
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);
650 }
651
652 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
653 SetPrecParameters();
654 prec_->initialize();
655 }
656
657 prec_->compute();
658}
659
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_);
665 if (!bA.is_null())
666 A_ = bA->Merge();
667
668 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
669
670 bool reusePreconditioner = false;
671
672 if (this->IsSetup() == true) {
673 // Reuse the constructed preconditioner
674 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupChebyshev(): Setup() has already been called, assuming reuse" << std::endl;
675
676 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
677 if (!prec.is_null()) {
678 prec->setMatrix(tA);
679 reusePreconditioner = true;
680 } else {
681 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available (failed cast to CanChangeMatrix), "
682 "reverting to full construction"
683 << std::endl;
684 }
685 }
686
687 // Take care of the eigenvalues
688 ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
689 SC negone = -STS::one();
690 SC lambdaMax = SetupChebyshevEigenvalues(currentLevel, "A", "", paramList);
691
692 if (!reusePreconditioner) {
693 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
694 SetPrecParameters();
695 {
696 SubFactoryMonitor m(*this, "Preconditioner init", currentLevel);
697 prec_->initialize();
698 }
699 } else
700 SetPrecParameters();
701
702 {
703 SubFactoryMonitor m(*this, "Preconditioner compute", currentLevel);
704 prec_->compute();
705 }
706
707 if (lambdaMax == negone) {
708 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
709
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_);
714 if (!matA.is_null())
715 matA->SetMaxEigenvalueEstimate(lambdaMax);
716 this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)"
717 << " = " << lambdaMax << std::endl;
718 }
719 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
720 }
721}
722
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_);
727 if (!bA.is_null())
728 A_ = bA->Merge();
729
730 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
731
732 bool reusePreconditioner = false;
733 if (this->IsSetup() == true) {
734 // Reuse the constructed preconditioner
735 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupHiptmair(): Setup() has already been called, assuming reuse" << std::endl;
736
737 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
738 if (!prec.is_null()) {
739 prec->setMatrix(tA);
740 reusePreconditioner = true;
741 } else {
742 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupHiptmair(): reuse of this type is not available (failed cast to CanChangeMatrix), "
743 "reverting to full construction"
744 << std::endl;
745 }
746 }
747
748 // If we're doing Chebyshev subsmoothing, we'll need to get the eigenvalues
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;
754
755 if (smoother1 == "CHEBYSHEV") {
756 ParameterList& list1 = paramList.sublist("hiptmair: smoother list 1");
757 lambdaMax11 = SetupChebyshevEigenvalues(currentLevel, "A", "EdgeMatrix ", list1);
758 }
759 if (smoother2 == "CHEBYSHEV") {
760 ParameterList& list2 = paramList.sublist("hiptmair: smoother list 2");
761 SetupChebyshevEigenvalues(currentLevel, "NodeMatrix", "NodeMatrix ", list2);
762 }
763
764 // FIXME: Should really add some checks to make sure the eigenvalue calcs worked like in
765 // the regular SetupChebyshev
766
767 // Grab the auxillary matrices and stick them on the list
768 RCP<Operator> NodeMatrix = currentLevel.Get<RCP<Operator>>("NodeMatrix");
769 RCP<Operator> D0 = currentLevel.Get<RCP<Operator>>("D0");
770
771 RCP<tRowMatrix> tNodeMatrix = Xpetra::toTpetraRowMatrix(NodeMatrix);
772 RCP<tRowMatrix> tD0 = Xpetra::toTpetraRowMatrix(D0);
773
774 Teuchos::ParameterList newlist;
775 newlist.set("P", tD0);
776 newlist.set("PtAP", tNodeMatrix);
777 if (!reusePreconditioner) {
778 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
779 SetPrecParameters(newlist);
780 prec_->initialize();
781 }
782
783 prec_->compute();
784
785 // Post-processing the (1,1) eigenvalue, if we have one
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_);
795 if (!matA.is_null())
796 matA->SetMaxEigenvalueEstimate(lambdaMax11);
797 this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)"
798 << " = " << lambdaMax11 << std::endl;
799 }
800 TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax11 == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
801 }
802 }
803}
804
805template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
806Scalar Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::SetupChebyshevEigenvalues(Level& currentLevel, const std::string& matrixName, const std::string& label, ParameterList& paramList) const {
807 // Helper: This gets used for smoothers that want to set up Chebyhev
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;
813
814 std::string maxEigString = "chebyshev: max eigenvalue";
815 std::string eigRatioString = "chebyshev: ratio eigenvalue";
816
817 // Get/calculate the maximum eigenvalue
818 if (paramList.isParameter(maxEigString)) {
819 if (paramList.isType<double>(maxEigString))
820 lambdaMax = paramList.get<double>(maxEigString);
821 else
822 lambdaMax = paramList.get<SC>(maxEigString);
823 this->GetOStream(Statistics1) << label << maxEigString << " (cached with smoother parameter list) = " << lambdaMax << std::endl;
824 if (!matA.is_null())
825 matA->SetMaxEigenvalueEstimate(lambdaMax);
826
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);
832 }
833 }
834
835 // Calculate the eigenvalue ratio
836 const SC defaultEigRatio = 20;
837
838 SC ratio = defaultEigRatio;
839 if (paramList.isParameter(eigRatioString)) {
840 if (paramList.isType<double>(eigRatioString))
841 ratio = paramList.get<double>(eigRatioString);
842 else
843 ratio = paramList.get<SC>(eigRatioString);
844 }
845 if (currentLevel.GetLevelID()) {
846 // Update ratio to be
847 // ratio = max(number of fine DOFs / number of coarse DOFs, defaultValue)
848 //
849 // NOTE: We don't need to request previous level matrix as we know for sure it was constructed
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();
853
854 SC levelRatio = as<SC>(as<double>(nRowsFine) / nRowsCoarse);
855 if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
856 ratio = levelRatio;
857 }
858
859 this->GetOStream(Statistics1) << label << eigRatioString << " (computed) = " << ratio << std::endl;
860 paramList.set(eigRatioString, ratio);
861
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);
872 }
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);
878 }
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);
884 }
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);
890 }
891 if (doScale) {
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());
897 }
898 }
899
900 return lambdaMax;
901}
902
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_);
907 if (!bA.is_null())
908 A_ = bA->Merge();
909
910 RCP<const tRowMatrix> tA = Xpetra::toTpetraRowMatrix(A_);
911
912 bool reusePreconditioner = false;
913 if (this->IsSetup() == true) {
914 // Reuse the constructed preconditioner
915 this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
916
917 RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix>> prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix>>(prec_);
918 if (!prec.is_null()) {
919 prec->setMatrix(tA);
920 reusePreconditioner = true;
921 } else {
922 this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
923 "reverting to full construction"
924 << std::endl;
925 }
926 }
927
928 if (!reusePreconditioner) {
929 prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
930 SetPrecParameters();
931 prec_->initialize();
932 }
933
934 prec_->compute();
935}
936
937template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
938void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
939 TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Apply(): Setup() has not been called");
940
941 // Forward the InitialGuessIsZero option to Ifpack2
942 // TODO: It might be nice to switch back the internal
943 // "zero starting solution" option of the ifpack2 object prec_ to his
944 // initial value at the end but there is no way right now to get
945 // the current value of the "zero starting solution" in ifpack2.
946 // It's not really an issue, as prec_ can only be used by this method.
947 bool supportInitialGuess = false;
948
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);
955 // Because additive Schwarz has "delta" semantics, it's sufficient to
956 // toggle only the zero initial guess flag, and not pass in already
957 // set parameters. If we call SetPrecParameters, the subdomain solver
958 // will be destroyed.
959 prec_->setParameters(paramList);
960 supportInitialGuess = true;
961 }
962
963 // TODO JJH 30Apr2014 Calling SetPrecParameters(paramList) when the smoother
964 // is Ifpack2::AdditiveSchwarz::setParameterList() will destroy the subdomain
965 //(aka inner) solver. This behavior is documented but a departure from what
966 // it previously did, and what other Ifpack2 solvers currently do. So I have
967 // moved SetPrecParameters(paramList) into the if-else block above.
968
969 // Apply
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);
974 } else {
975 using TST = Teuchos::ScalarTraits<Scalar>;
976
977 RCP<MultiVector> Residual;
978 {
979 std::string prefix = this->ShortClassName() + ": Apply: ";
980 RCP<TimeMonitor> tM = rcp(new TimeMonitor(*this, prefix + "residual calculation", Timings0));
981 Residual = Utilities::Residual(*A_, X, B);
982 }
983
984 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
985
986 Tpetra::MultiVector<SC, LO, GO, NO>& tpX = toTpetra(*Correction);
987 const Tpetra::MultiVector<SC, LO, GO, NO>& tpB = toTpetra(*Residual);
988
989 prec_->apply(tpB, tpX);
990
991 X.update(TST::one(), *Correction, TST::one());
992 }
993}
994
995template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
996RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
997 RCP<Ifpack2Smoother> smoother = rcp(new Ifpack2Smoother(*this));
998 smoother->SetParameterList(this->GetParameterList());
999 return smoother;
1000}
1001
1002template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1004 std::ostringstream out;
1006 out << prec_->description();
1007 } else {
1009 out << "{type = " << type_ << "}";
1010 }
1011 return out.str();
1012}
1013
1014template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1015void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream& out, const VerbLevel verbLevel) const {
1017
1018 if (verbLevel & Parameters0)
1019 out0 << "Prec. type: " << type_ << std::endl;
1020
1021 if (verbLevel & Parameters1) {
1022 out0 << "Parameter list: " << std::endl;
1023 Teuchos::OSTab tab2(out);
1024 out << this->GetParameterList();
1025 out0 << "Overlap: " << overlap_ << std::endl;
1026 }
1027
1028 if (verbLevel & External)
1029 if (prec_ != Teuchos::null) {
1030 Teuchos::OSTab tab2(out);
1031 out << *prec_ << std::endl;
1032 }
1033
1034 if (verbLevel & Debug) {
1035 out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
1036 << "-" << std::endl
1037 << "RCP<prec_>: " << prec_ << std::endl;
1038 }
1039}
1040
1041template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1043 typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
1044 // NOTE: Only works for a subset of Ifpack2's smoothers
1045 RCP<Ifpack2::Relaxation<MatrixType>> pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType>>(prec_);
1046 if (!pr.is_null()) return pr->getNodeSmootherComplexity();
1047
1048 RCP<Ifpack2::Chebyshev<MatrixType>> pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType>>(prec_);
1049 if (!pc.is_null()) return pc->getNodeSmootherComplexity();
1050
1051 RCP<Ifpack2::BlockRelaxation<MatrixType>> pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType>>(prec_);
1052 if (!pb.is_null()) return pb->getNodeSmootherComplexity();
1053
1054 RCP<Ifpack2::ILUT<MatrixType>> pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType>>(prec_);
1055 if (!pi.is_null()) return pi->getNodeSmootherComplexity();
1056
1057 RCP<Ifpack2::RILUK<MatrixType>> pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType>>(prec_);
1058 if (!pk.is_null()) return pk->getNodeSmootherComplexity();
1059
1060 return Teuchos::OrdinalTraits<size_t>::invalid();
1061}
1062
1063} // namespace MueLu
1064
1065#endif // MUELU_IFPACK2SMOOTHER_DEF_HPP
#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 &currentLevel)
void Setup(Level &currentLevel)
Set up the smoother.
Scalar SetupChebyshevEigenvalues(Level &currentLevel, const std::string &matrixName, const std::string &label, Teuchos::ParameterList &paramList) const
void SetupAggregate(Level &currentLevel)
void SetupBlockRelaxation(Level &currentLevel)
void SetupChebyshev(Level &currentLevel)
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetParameterList(const Teuchos::ParameterList &paramList)
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 &currentLevel)
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
void SetupTopological(Level &currentLevel)
std::string type_
ifpack2-specific key phrase that denote smoother type
friend class Ifpack2Smoother
Constructor.
void DeclareInput(Level &currentLevel) const
Input.
void SetupLineSmoothing(Level &currentLevel)
void SetupGeneric(Level &currentLevel)
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 &paramList)=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)