MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
15 #if defined(HAVE_MUELU_IFPACK2)
16 
18 
19 #include <Tpetra_RowMatrix.hpp>
20 
21 #include <Ifpack2_Chebyshev.hpp>
22 #include <Ifpack2_Hiptmair.hpp>
23 #include <Ifpack2_RILUK.hpp>
24 #include <Ifpack2_Relaxation.hpp>
25 #include <Ifpack2_ILUT.hpp>
26 #include <Ifpack2_BlockRelaxation.hpp>
27 #include <Ifpack2_Factory.hpp>
28 #include <Ifpack2_Parameters.hpp>
29 
31 #include <Xpetra_CrsMatrix.hpp>
32 #include <Xpetra_CrsMatrixWrap.hpp>
33 #include <Xpetra_TpetraBlockCrsMatrix.hpp>
34 #include <Xpetra_Matrix.hpp>
35 #include <Xpetra_MatrixMatrix.hpp>
36 #include <Xpetra_MultiVectorFactory.hpp>
37 #include <Xpetra_TpetraMultiVector.hpp>
38 
39 #include <Tpetra_BlockCrsMatrix_Helpers.hpp>
40 
42 #include "MueLu_Level.hpp"
43 #include "MueLu_Utilities.hpp"
44 #include "MueLu_Monitor.hpp"
45 #include "MueLu_Aggregates.hpp"
46 
47 #ifdef HAVE_MUELU_INTREPID2
50 #include "Intrepid2_Basis.hpp"
51 #include "Kokkos_DynRankView.hpp"
52 #endif
53 
54 namespace MueLu {
55 
56 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
58  : type_(type)
59  , overlap_(overlap) {
60  typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
61  bool isSupported = Ifpack2::Factory::isSupported<tRowMatrix>(type_) || (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
62  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
63  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
64  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
65  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
66  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
67  type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
68  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
69  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
70  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
71  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
72  type_ == "LINESMOOTHING_BLOCKRELAXATION" ||
73  type_ == "TOPOLOGICAL" ||
74  type_ == "AGGREGATE");
75  this->declareConstructionOutcome(!isSupported, "Ifpack2 does not provide the smoother '" + type_ + "'.");
76  if (isSupported)
77  SetParameterList(paramList);
78 }
79 
80 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82  Factory::SetParameterList(paramList);
83 
85  // It might be invalid to change parameters after the setup, but it depends entirely on Ifpack implementation.
86  // TODO: I don't know if Ifpack returns an error code or exception or ignore parameters modification in this case...
87  SetPrecParameters();
88  }
89 }
90 
91 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
93  std::string prefix = this->ShortClassName() + ": SetPrecParameters";
94  RCP<TimeMonitor> tM = rcp(new TimeMonitor(*this, prefix, Timings0));
95  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
96 
97  paramList.setParameters(list);
98 
99  RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
100 
101  // Do we want an Ifpack2 apply timer?
102  precList->set("timer for apply", this->IsPrint(Timings));
103 
104  if (!precList.is_null() && precList->isParameter("partitioner: type") && precList->get<std::string>("partitioner: type") == "linear" &&
105  !precList->isParameter("partitioner: local parts")) {
106  LO matrixBlockSize = 1;
107  int lclSize = A_->getRangeMap()->getLocalNumElements();
108  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
109  if (!matA.is_null()) {
110  lclSize = matA->getLocalNumRows();
111  matrixBlockSize = matA->GetFixedBlockSize();
112  }
113  precList->set("partitioner: local parts", lclSize / matrixBlockSize);
114  }
115 
116  prec_->setParameters(*precList);
117 
118  paramList.setParameters(*precList); // what about that??
119 }
120 
121 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
123  this->Input(currentLevel, "A");
124 
125  if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
126  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
127  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
128  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
129  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
130  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
131  type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
132  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
133  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
134  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
135  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
136  type_ == "LINESMOOTHING_BLOCKRELAXATION") {
137  this->Input(currentLevel, "CoarseNumZLayers"); // necessary for fallback criterion
138  this->Input(currentLevel, "LineDetection_VertLineIds"); // necessary to feed block smoother
139  } else if (type_ == "BLOCK RELAXATION" ||
140  type_ == "BLOCK_RELAXATION" ||
141  type_ == "BLOCKRELAXATION" ||
142  // Banded
143  type_ == "BANDED_RELAXATION" ||
144  type_ == "BANDED RELAXATION" ||
145  type_ == "BANDEDRELAXATION" ||
146  // Tridiagonal
147  type_ == "TRIDI_RELAXATION" ||
148  type_ == "TRIDI RELAXATION" ||
149  type_ == "TRIDIRELAXATION" ||
150  type_ == "TRIDIAGONAL_RELAXATION" ||
151  type_ == "TRIDIAGONAL RELAXATION" ||
152  type_ == "TRIDIAGONALRELAXATION") {
153  // We need to check for the "partitioner type" = "line"
154  ParameterList precList = this->GetParameterList();
155  if (precList.isParameter("partitioner: type") &&
156  precList.get<std::string>("partitioner: type") == "line") {
157  this->Input(currentLevel, "Coordinates");
158  }
159  } else if (type_ == "TOPOLOGICAL") {
160  // for the topological smoother, we require an element to node map:
161  this->Input(currentLevel, "pcoarsen: element to node map");
162  } else if (type_ == "AGGREGATE") {
163  // Aggregate smoothing needs aggregates
164  this->Input(currentLevel, "Aggregates");
165  } else if (type_ == "HIPTMAIR") {
166  // Hiptmair needs D0 and NodeMatrix
167  this->Input(currentLevel, "NodeMatrix");
168  this->Input(currentLevel, "D0");
169  }
170 }
171 
172 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174  FactoryMonitor m(*this, "Setup Smoother", currentLevel);
175  A_ = Factory::Get<RCP<Operator>>(currentLevel, "A");
176  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
177 
178  // If the user asked us to convert the matrix into BlockCrsMatrix form, we do that now
179 
180  if (paramList.isParameter("smoother: use blockcrsmatrix storage") && paramList.get<bool>("smoother: use blockcrsmatrix storage")) {
181  int blocksize = 1;
182  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
183  if (!matA.is_null())
184  blocksize = matA->GetFixedBlockSize();
185  if (blocksize) {
186  // 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
187 
188  RCP<CrsMatrixWrap> AcrsWrap = rcp_dynamic_cast<CrsMatrixWrap>(A_);
189  if (AcrsWrap.is_null())
190  throw std::runtime_error("Ifpack2Smoother: Cannot convert matrix A to CrsMatrixWrap object.");
191  RCP<CrsMatrix> Acrs = AcrsWrap->getCrsMatrix();
192  if (Acrs.is_null())
193  throw std::runtime_error("Ifpack2Smoother: Cannot extract CrsMatrix from matrix A.");
194  RCP<TpetraCrsMatrix> At = rcp_dynamic_cast<TpetraCrsMatrix>(Acrs);
195  if (At.is_null()) {
197  throw std::runtime_error("Ifpack2Smoother: Cannot extract CrsMatrix or BlockCrsMatrix from matrix A.");
198  this->GetOStream(Statistics0) << "Ifpack2Smoother: Using (native) BlockCrsMatrix storage with blocksize " << blocksize << std::endl;
199  paramList.remove("smoother: use blockcrsmatrix storage");
200  } else {
201  RCP<Tpetra::BlockCrsMatrix<Scalar, LO, GO, Node>> blockCrs = Tpetra::convertToBlockCrsMatrix(*At->getTpetra_CrsMatrix(), blocksize);
202  RCP<CrsMatrix> blockCrs_as_crs = rcp(new TpetraBlockCrsMatrix(blockCrs));
203  RCP<CrsMatrixWrap> blockWrap = rcp(new CrsMatrixWrap(blockCrs_as_crs));
204  A_ = blockWrap;
205  this->GetOStream(Statistics0) << "Ifpack2Smoother: Using BlockCrsMatrix storage with blocksize " << blocksize << std::endl;
206 
207  paramList.remove("smoother: use blockcrsmatrix storage");
208  }
209  }
210  }
211 
212  if (type_ == "SCHWARZ")
213  SetupSchwarz(currentLevel);
214 
215  else if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
216  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
217  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
218  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
219  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
220  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
221  type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
222  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
223  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
224  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
225  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
226  type_ == "LINESMOOTHING_BLOCKRELAXATION")
227  SetupLineSmoothing(currentLevel);
228 
229  else if (type_ == "BLOCK_RELAXATION" ||
230  type_ == "BLOCK RELAXATION" ||
231  type_ == "BLOCKRELAXATION" ||
232  // Banded
233  type_ == "BANDED_RELAXATION" ||
234  type_ == "BANDED RELAXATION" ||
235  type_ == "BANDEDRELAXATION" ||
236  // Tridiagonal
237  type_ == "TRIDI_RELAXATION" ||
238  type_ == "TRIDI RELAXATION" ||
239  type_ == "TRIDIRELAXATION" ||
240  type_ == "TRIDIAGONAL_RELAXATION" ||
241  type_ == "TRIDIAGONAL RELAXATION" ||
242  type_ == "TRIDIAGONALRELAXATION")
243  SetupBlockRelaxation(currentLevel);
244 
245  else if (type_ == "CHEBYSHEV")
246  SetupChebyshev(currentLevel);
247 
248  else if (type_ == "TOPOLOGICAL") {
249 #ifdef HAVE_MUELU_INTREPID2
250  SetupTopological(currentLevel);
251 #else
252  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, "'TOPOLOGICAL' smoother choice requires Intrepid2");
253 #endif
254  } else if (type_ == "AGGREGATE")
255  SetupAggregate(currentLevel);
256 
257  else if (type_ == "HIPTMAIR")
258  SetupHiptmair(currentLevel);
259 
260  else {
261  SetupGeneric(currentLevel);
262  }
263 
265 
266  this->GetOStream(Statistics1) << description() << std::endl;
267 }
268 
269 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271  typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
272 
273  bool reusePreconditioner = false;
274  if (this->IsSetup() == true) {
275  // Reuse the constructed preconditioner
276  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
277 
278  bool isTRowMatrix = true;
280  try {
282  } catch (Exceptions::BadCast&) {
283  isTRowMatrix = false;
284  }
285 
287  if (!prec.is_null() && isTRowMatrix) {
288  prec->setMatrix(tA);
289  reusePreconditioner = true;
290  } else {
291  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
292  "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction"
293  << std::endl;
294  }
295  }
296 
297  if (!reusePreconditioner) {
298  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
299 
300  bool isBlockedMatrix = false;
301  RCP<Matrix> merged2Mat;
302 
303  std::string sublistName = "subdomain solver parameters";
304  if (paramList.isSublist(sublistName)) {
305  // If we are doing "user" partitioning, we assume that what the user
306  // really wants to do is make tiny little subdomains with one row
307  // assigned to each subdomain. The rows used for these little
308  // subdomains correspond to those in the 2nd block row. Then,
309  // if we overlap these mini-subdomains, we will do something that
310  // looks like Vanka (grabbing all velocities associated with each
311  // each pressure unknown). In addition, we put all Dirichlet points
312  // as a little mini-domain.
313  ParameterList& subList = paramList.sublist(sublistName);
314 
315  std::string partName = "partitioner: type";
316  // Pretty sure no one has been using this. Unfortunately, old if
317  // statement (which checked for equality with "user") prevented
318  // anyone from employing other types of Ifpack2 user partition
319  // options. Leaving this and switching if to "vanka user" just
320  // in case some day someone might want to use this.
321  if (subList.isParameter(partName) && subList.get<std::string>(partName) == "vanka user") {
322  isBlockedMatrix = true;
323 
324  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
326  "Matrix A must be of type BlockedCrsMatrix.");
327 
328  size_t numVels = bA->getMatrix(0, 0)->getLocalNumRows();
329  size_t numPres = bA->getMatrix(1, 0)->getLocalNumRows();
330  size_t numRows = rcp_dynamic_cast<Matrix>(A_, true)->getLocalNumRows();
331 
333 
334  size_t numBlocks = 0;
335  for (size_t rowOfB = numVels; rowOfB < numVels + numPres; ++rowOfB)
336  blockSeeds[rowOfB] = numBlocks++;
337 
338  RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
340  "Matrix A must be of type BlockedCrsMatrix.");
341 
342  merged2Mat = bA2->Merge();
343 
344  // Add Dirichlet rows to the list of seeds
345  ArrayRCP<const bool> boundaryNodes = Utilities::DetectDirichletRows(*merged2Mat, 0.0);
346  bool haveBoundary = false;
347  for (LO i = 0; i < boundaryNodes.size(); i++)
348  if (boundaryNodes[i]) {
349  // FIXME:
350  // 1. would not this [] overlap with some in the previos blockSeed loop?
351  // 2. do we need to distinguish between pressure and velocity Dirichlet b.c.
352  blockSeeds[i] = numBlocks;
353  haveBoundary = true;
354  }
355  if (haveBoundary)
356  numBlocks++;
357 
358  subList.set("partitioner: type", "user");
359  subList.set("partitioner: map", blockSeeds);
360  subList.set("partitioner: local parts", as<int>(numBlocks));
361 
362  } else {
363  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
364  if (!bA.is_null()) {
365  isBlockedMatrix = true;
366  merged2Mat = bA->Merge();
367  }
368  }
369  }
370 
372  if (isBlockedMatrix == true)
373  tA = Utilities::Op2NonConstTpetraRow(merged2Mat);
374  else
376 
377  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
378  SetPrecParameters();
379 
380  prec_->initialize();
381  }
382 
383  prec_->compute();
384 }
385 
386 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
388  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
389 
390  if (this->IsSetup() == true) {
391  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupAggregate(): Setup() has already been called" << std::endl;
392  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupAggregate(): reuse of this type is not available, reverting to full construction" << std::endl;
393  }
394 
395  this->GetOStream(Statistics0) << "Ifpack2Smoother: Using Aggregate Smoothing" << std::endl;
396 
397  RCP<Aggregates> aggregates = Factory::Get<RCP<Aggregates>>(currentLevel, "Aggregates");
398 
399  RCP<const LOMultiVector> vertex2AggId = aggregates->GetVertex2AggId();
400  ArrayRCP<LO> aggregate_ids = rcp_const_cast<LOMultiVector>(vertex2AggId)->getDataNonConst(0);
401  ArrayRCP<LO> dof_ids;
402 
403  // We need to unamalgamate, if the FixedBlockSize > 1
404  LO blocksize = 1;
405  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
406  if (!matA.is_null())
407  blocksize = matA->GetFixedBlockSize();
408  if (blocksize > 1) {
409  dof_ids.resize(aggregate_ids.size() * blocksize);
410  for (LO i = 0; i < (LO)aggregate_ids.size(); i++) {
411  for (LO j = 0; j < (LO)blocksize; j++)
412  dof_ids[i * blocksize + j] = aggregate_ids[i];
413  }
414  } else {
415  dof_ids = aggregate_ids;
416  }
417 
418  paramList.set("partitioner: map", dof_ids);
419  paramList.set("partitioner: type", "user");
420  paramList.set("partitioner: overlap", 0);
421  paramList.set("partitioner: local parts", (int)aggregates->GetNumAggregates());
422 
424 
425  type_ = "BLOCKRELAXATION";
426  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
427  SetPrecParameters();
428  prec_->initialize();
429  prec_->compute();
430 }
431 
432 #ifdef HAVE_MUELU_INTREPID2
433 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
435  /*
436 
437  basic notion:
438 
439  Look for user input indicating topo dimension, something like "topological domain type: {node|edge|face}"
440  Call something like what you can find in Poisson example line 1180 to set seeds for a smoother
441 
442  */
443  if (this->IsSetup() == true) {
444  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
445  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
446  }
447 
448  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
449 
450  typedef typename Node::device_type::execution_space ES;
451 
452  typedef Kokkos::DynRankView<LocalOrdinal, typename Node::device_type> FCO; //
453 
455 
456  using namespace std;
457 
458  const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO>>(currentLevel, "pcoarsen: element to node map");
459 
460  string basisString = paramList.get<string>("pcoarsen: hi basis");
461  int degree;
462  // NOTE: To make sure Stokhos works we only instantiate these guys with double. There's a lot
463  // of stuff in the guts of Intrepid2 that doesn't play well with Stokhos as of yet. Here, we only
464  // care about the assignment of basis ordinals to topological entities, so this code is actually
465  // independent of the Scalar type--hard-coding double here won't hurt us.
466  auto basis = MueLuIntrepid::BasisFactory<double, ES>(basisString, degree);
467 
468  string topologyTypeString = paramList.get<string>("smoother: neighborhood type");
469  int dimension;
470  if (topologyTypeString == "node")
471  dimension = 0;
472  else if (topologyTypeString == "edge")
473  dimension = 1;
474  else if (topologyTypeString == "face")
475  dimension = 2;
476  else if (topologyTypeString == "cell")
477  dimension = basis->getBaseCellTopology().getDimension();
478  else
479  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
480  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_, true);
481  vector<vector<LocalOrdinal>> seeds;
482  MueLuIntrepid::FindGeometricSeedOrdinals(basis, *elemToNode, seeds, *matA->getRowMap(), *matA->getColMap());
483 
484  // Ifpack2 wants the seeds in an array of the same length as the number of local elements,
485  // with local partition #s marked for the ones that are seeds, and invalid for the rest
486  int myNodeCount = matA->getRowMap()->getLocalNumElements();
487  ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount, lo_invalid);
488  int localPartitionNumber = 0;
489  for (LocalOrdinal seed : seeds[dimension]) {
490  nodeSeeds[seed] = localPartitionNumber++;
491  }
492 
493  paramList.remove("smoother: neighborhood type");
494  paramList.remove("pcoarsen: hi basis");
495 
496  paramList.set("partitioner: map", nodeSeeds);
497  paramList.set("partitioner: type", "user");
498  paramList.set("partitioner: overlap", 1);
499  paramList.set("partitioner: local parts", int(seeds[dimension].size()));
500 
501  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO>> tA = Utilities::Op2NonConstTpetraRow(A_);
502 
503  type_ = "BLOCKRELAXATION";
504  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
505  SetPrecParameters();
506  prec_->initialize();
507  prec_->compute();
508 }
509 #endif
510 
511 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
513  if (this->IsSetup() == true) {
514  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
515  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
516  }
517 
518  ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
519 
520  LO CoarseNumZLayers = Factory::Get<LO>(currentLevel, "CoarseNumZLayers");
521  if (CoarseNumZLayers > 0) {
522  Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get<Teuchos::ArrayRCP<LO>>(currentLevel, "LineDetection_VertLineIds");
523 
524  // determine number of local parts
525  LO maxPart = 0;
526  for (size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
527  if (maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
528  }
529  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_, true);
530  size_t numLocalRows = matA->getLocalNumRows();
531 
532  TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % TVertLineIdSmoo.size() != 0, Exceptions::RuntimeError,
533  "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
534 
535  // actualDofsPerNode is the actual number of matrix rows per mesh element.
536  // It is encoded in either the MueLu Level, or in the Xpetra matrix block size.
537  // This value is needed by Ifpack2 to do decoupled block relaxation.
538  int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
539  LO matrixBlockSize = matA->GetFixedBlockSize();
540  if (matrixBlockSize > 1 && actualDofsPerNode > 1) {
541  TEUCHOS_TEST_FOR_EXCEPTION(actualDofsPerNode != matrixBlockSize, Exceptions::RuntimeError,
542  "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
543  } else if (matrixBlockSize > 1) {
544  actualDofsPerNode = matrixBlockSize;
545  }
546  myparamList.set("partitioner: PDE equations", actualDofsPerNode);
547 
548  if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
549  myparamList.set("partitioner: type", "user");
550  myparamList.set("partitioner: map", TVertLineIdSmoo);
551  myparamList.set("partitioner: local parts", maxPart + 1);
552  } else {
553  if (myparamList.isParameter("partitioner: block size") &&
554  myparamList.get<int>("partitioner: block size") != -1) {
555  int block_size = myparamList.get<int>("partitioner: block size");
556  TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % block_size != 0, Exceptions::RuntimeError,
557  "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the specified block size.");
558  numLocalRows /= block_size;
559  }
560 
561  // we assume a constant number of DOFs per node
562  size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
563 
564  // Create a new Teuchos::ArrayRCP<LO> of size numLocalRows and fill it with the corresponding information
566  for (size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
567  for (size_t dof = 0; dof < numDofsPerNode; dof++)
568  partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
569  myparamList.set("partitioner: type", "user");
570  myparamList.set("partitioner: map", partitionerMap);
571  myparamList.set("partitioner: local parts", maxPart + 1);
572  }
573 
574  if (type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
575  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
576  type_ == "LINESMOOTHING_BANDEDRELAXATION")
577  type_ = "BANDEDRELAXATION";
578  else if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
579  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
580  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
581  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
582  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
583  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION")
584  type_ = "TRIDIAGONALRELAXATION";
585  else
586  type_ = "BLOCKRELAXATION";
587  } else {
588  // line detection failed -> fallback to point-wise relaxation
589  this->GetOStream(Runtime0) << "Line detection failed: fall back to point-wise relaxation" << std::endl;
590  myparamList.remove("partitioner: type", false);
591  myparamList.remove("partitioner: map", false);
592  myparamList.remove("partitioner: local parts", false);
593  type_ = "RELAXATION";
594  }
595 
597 
598  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
599  SetPrecParameters();
600  prec_->initialize();
601  prec_->compute();
602 }
603 
604 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
606  typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
607 
608  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
609  if (!bA.is_null())
610  A_ = bA->Merge();
611 
613 
614  bool reusePreconditioner = false;
615  if (this->IsSetup() == true) {
616  // Reuse the constructed preconditioner
617  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
618 
620  if (!prec.is_null()) {
621  prec->setMatrix(tA);
622  reusePreconditioner = true;
623  } else {
624  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
625  "reverting to full construction"
626  << std::endl;
627  }
628  }
629 
630  if (!reusePreconditioner) {
631  ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
632  myparamList.print();
633  if (myparamList.isParameter("partitioner: type") &&
634  myparamList.get<std::string>("partitioner: type") == "line") {
636  Factory::Get<Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LO, GO, NO>>>(currentLevel, "Coordinates");
638 
639  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
640  size_t lclSize = A_->getRangeMap()->getLocalNumElements();
641  if (!matA.is_null())
642  lclSize = matA->getLocalNumRows();
643  size_t numDofsPerNode = lclSize / xCoordinates->getMap()->getLocalNumElements();
644  myparamList.set("partitioner: coordinates", coordinates);
645  myparamList.set("partitioner: PDE equations", (int)numDofsPerNode);
646  }
647 
648  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
649  SetPrecParameters();
650  prec_->initialize();
651  }
652 
653  prec_->compute();
654 }
655 
656 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
658  typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
659  using STS = Teuchos::ScalarTraits<SC>;
660  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
661  if (!bA.is_null())
662  A_ = bA->Merge();
663 
665 
666  bool reusePreconditioner = false;
667 
668  if (this->IsSetup() == true) {
669  // Reuse the constructed preconditioner
670  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupChebyshev(): Setup() has already been called, assuming reuse" << std::endl;
671 
673  if (!prec.is_null()) {
674  prec->setMatrix(tA);
675  reusePreconditioner = true;
676  } else {
677  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available (failed cast to CanChangeMatrix), "
678  "reverting to full construction"
679  << std::endl;
680  }
681  }
682 
683  // Take care of the eigenvalues
684  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
685  SC negone = -STS::one();
686  SC lambdaMax = SetupChebyshevEigenvalues(currentLevel, "A", "", paramList);
687 
688  if (!reusePreconditioner) {
689  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
690  SetPrecParameters();
691  {
692  SubFactoryMonitor m(*this, "Preconditioner init", currentLevel);
693  prec_->initialize();
694  }
695  } else
696  SetPrecParameters();
697 
698  {
699  SubFactoryMonitor m(*this, "Preconditioner compute", currentLevel);
700  prec_->compute();
701  }
702 
703  if (lambdaMax == negone) {
704  typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
705 
707  if (chebyPrec != Teuchos::null) {
708  lambdaMax = chebyPrec->getLambdaMaxForApply();
709  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
710  if (!matA.is_null())
711  matA->SetMaxEigenvalueEstimate(lambdaMax);
712  this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)"
713  << " = " << lambdaMax << std::endl;
714  }
715  TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
716  }
717 }
718 
719 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
721  typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
722  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
723  if (!bA.is_null())
724  A_ = bA->Merge();
725 
727 
728  bool reusePreconditioner = false;
729  if (this->IsSetup() == true) {
730  // Reuse the constructed preconditioner
731  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupHiptmair(): Setup() has already been called, assuming reuse" << std::endl;
732 
734  if (!prec.is_null()) {
735  prec->setMatrix(tA);
736  reusePreconditioner = true;
737  } else {
738  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupHiptmair(): reuse of this type is not available (failed cast to CanChangeMatrix), "
739  "reverting to full construction"
740  << std::endl;
741  }
742  }
743 
744  // If we're doing Chebyshev subsmoothing, we'll need to get the eigenvalues
746  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
747  std::string smoother1 = paramList.get("hiptmair: smoother type 1", "CHEBYSHEV");
748  std::string smoother2 = paramList.get("hiptmair: smoother type 2", "CHEBYSHEV");
749  SC lambdaMax11 = negone;
750 
751  if (smoother1 == "CHEBYSHEV") {
752  ParameterList& list1 = paramList.sublist("hiptmair: smoother list 1");
753  lambdaMax11 = SetupChebyshevEigenvalues(currentLevel, "A", "EdgeMatrix ", list1);
754  }
755  if (smoother2 == "CHEBYSHEV") {
756  ParameterList& list2 = paramList.sublist("hiptmair: smoother list 2");
757  SetupChebyshevEigenvalues(currentLevel, "NodeMatrix", "NodeMatrix ", list2);
758  }
759 
760  // FIXME: Should really add some checks to make sure the eigenvalue calcs worked like in
761  // the regular SetupChebyshev
762 
763  // Grab the auxillary matrices and stick them on the list
764  RCP<Operator> NodeMatrix = currentLevel.Get<RCP<Operator>>("NodeMatrix");
765  RCP<Operator> D0 = currentLevel.Get<RCP<Operator>>("D0");
766 
767  RCP<tRowMatrix> tNodeMatrix = Utilities::Op2NonConstTpetraRow(NodeMatrix);
769 
770  Teuchos::ParameterList newlist;
771  newlist.set("P", tD0);
772  newlist.set("PtAP", tNodeMatrix);
773  if (!reusePreconditioner) {
774  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
775  SetPrecParameters(newlist);
776  prec_->initialize();
777  }
778 
779  prec_->compute();
780 
781  // Post-processing the (1,1) eigenvalue, if we have one
782  if (smoother1 == "CHEBYSHEV" && lambdaMax11 == negone) {
783  using Teuchos::rcp_dynamic_cast;
784  typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
785  auto hiptmairPrec = rcp_dynamic_cast<Ifpack2::Hiptmair<MatrixType>>(prec_);
786  if (hiptmairPrec != Teuchos::null) {
787  auto chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType>>(hiptmairPrec->getPrec1());
788  if (chebyPrec != Teuchos::null) {
789  lambdaMax11 = chebyPrec->getLambdaMaxForApply();
790  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(A_);
791  if (!matA.is_null())
792  matA->SetMaxEigenvalueEstimate(lambdaMax11);
793  this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)"
794  << " = " << lambdaMax11 << std::endl;
795  }
796  TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax11 == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
797  }
798  }
799 }
800 
801 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
802 Scalar Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::SetupChebyshevEigenvalues(Level& currentLevel, const std::string& matrixName, const std::string& label, ParameterList& paramList) const {
803  // Helper: This gets used for smoothers that want to set up Chebyhev
804  typedef Teuchos::ScalarTraits<SC> STS;
805  SC negone = -STS::one();
806  RCP<Operator> currentA = currentLevel.Get<RCP<Operator>>(matrixName);
807  RCP<Matrix> matA = rcp_dynamic_cast<Matrix>(currentA);
808  SC lambdaMax = negone;
809 
810  std::string maxEigString = "chebyshev: max eigenvalue";
811  std::string eigRatioString = "chebyshev: ratio eigenvalue";
812 
813  // Get/calculate the maximum eigenvalue
814  if (paramList.isParameter(maxEigString)) {
815  if (paramList.isType<double>(maxEigString))
816  lambdaMax = paramList.get<double>(maxEigString);
817  else
818  lambdaMax = paramList.get<SC>(maxEigString);
819  this->GetOStream(Statistics1) << label << maxEigString << " (cached with smoother parameter list) = " << lambdaMax << std::endl;
820  if (!matA.is_null())
821  matA->SetMaxEigenvalueEstimate(lambdaMax);
822 
823  } else if (!matA.is_null()) {
824  lambdaMax = matA->GetMaxEigenvalueEstimate();
825  if (lambdaMax != negone) {
826  this->GetOStream(Statistics1) << label << maxEigString << " (cached with matrix) = " << lambdaMax << std::endl;
827  paramList.set(maxEigString, lambdaMax);
828  }
829  }
830 
831  // Calculate the eigenvalue ratio
832  const SC defaultEigRatio = 20;
833 
834  SC ratio = defaultEigRatio;
835  if (paramList.isParameter(eigRatioString)) {
836  if (paramList.isType<double>(eigRatioString))
837  ratio = paramList.get<double>(eigRatioString);
838  else
839  ratio = paramList.get<SC>(eigRatioString);
840  }
841  if (currentLevel.GetLevelID()) {
842  // Update ratio to be
843  // ratio = max(number of fine DOFs / number of coarse DOFs, defaultValue)
844  //
845  // NOTE: We don't need to request previous level matrix as we know for sure it was constructed
846  RCP<const Operator> fineA = currentLevel.GetPreviousLevel()->Get<RCP<Operator>>(matrixName);
847  size_t nRowsFine = fineA->getDomainMap()->getGlobalNumElements();
848  size_t nRowsCoarse = currentA->getDomainMap()->getGlobalNumElements();
849 
850  SC levelRatio = as<SC>(as<double>(nRowsFine) / nRowsCoarse);
851  if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
852  ratio = levelRatio;
853  }
854 
855  this->GetOStream(Statistics1) << label << eigRatioString << " (computed) = " << ratio << std::endl;
856  paramList.set(eigRatioString, ratio);
857 
858  if (paramList.isParameter("chebyshev: use rowsumabs diagonal scaling")) {
859  this->GetOStream(Runtime1) << "chebyshev: using rowsumabs diagonal scaling" << std::endl;
860  bool doScale = false;
861  doScale = paramList.get<bool>("chebyshev: use rowsumabs diagonal scaling");
862  paramList.remove("chebyshev: use rowsumabs diagonal scaling");
863  double chebyReplaceTol = Teuchos::ScalarTraits<Scalar>::eps() * 100;
864  std::string paramName = "chebyshev: rowsumabs diagonal replacement tolerance";
865  if (paramList.isParameter(paramName)) {
866  chebyReplaceTol = paramList.get<double>(paramName);
867  paramList.remove(paramName);
868  }
869  double chebyReplaceVal = Teuchos::ScalarTraits<double>::zero();
870  paramName = "chebyshev: rowsumabs diagonal replacement value";
871  if (paramList.isParameter(paramName)) {
872  chebyReplaceVal = paramList.get<double>(paramName);
873  paramList.remove(paramName);
874  }
875  bool chebyReplaceSingleEntryRowWithZero = false;
876  paramName = "chebyshev: rowsumabs replace single entry row with zero";
877  if (paramList.isParameter(paramName)) {
878  chebyReplaceSingleEntryRowWithZero = paramList.get<bool>(paramName);
879  paramList.remove(paramName);
880  }
881  bool useAverageAbsDiagVal = false;
882  paramName = "chebyshev: rowsumabs use automatic diagonal tolerance";
883  if (paramList.isParameter(paramName)) {
884  useAverageAbsDiagVal = paramList.get<bool>(paramName);
885  paramList.remove(paramName);
886  }
887  if (doScale) {
888  TEUCHOS_ASSERT(!matA.is_null());
889  const bool doReciprocal = true;
890  RCP<Vector> lumpedDiagonal = Utilities::GetLumpedMatrixDiagonal(*matA, doReciprocal, chebyReplaceTol, chebyReplaceVal, chebyReplaceSingleEntryRowWithZero, useAverageAbsDiagVal);
892  paramList.set("chebyshev: operator inv diagonal", tmpVec.getTpetra_Vector());
893  }
894  }
895 
896  return lambdaMax;
897 }
898 
899 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
901  typedef Tpetra::RowMatrix<SC, LO, GO, NO> tRowMatrix;
902  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
903  if (!bA.is_null())
904  A_ = bA->Merge();
905 
907 
908  bool reusePreconditioner = false;
909  if (this->IsSetup() == true) {
910  // Reuse the constructed preconditioner
911  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
912 
914  if (!prec.is_null()) {
915  prec->setMatrix(tA);
916  reusePreconditioner = true;
917  } else {
918  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
919  "reverting to full construction"
920  << std::endl;
921  }
922  }
923 
924  if (!reusePreconditioner) {
925  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
926  SetPrecParameters();
927  prec_->initialize();
928  }
929 
930  prec_->compute();
931 }
932 
933 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
934 void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
935  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Apply(): Setup() has not been called");
936 
937  // Forward the InitialGuessIsZero option to Ifpack2
938  // TODO: It might be nice to switch back the internal
939  // "zero starting solution" option of the ifpack2 object prec_ to his
940  // initial value at the end but there is no way right now to get
941  // the current value of the "zero starting solution" in ifpack2.
942  // It's not really an issue, as prec_ can only be used by this method.
943  Teuchos::ParameterList paramList;
944  bool supportInitialGuess = false;
945  const Teuchos::ParameterList params = this->GetParameterList();
946 
947  if (prec_->supportsZeroStartingSolution()) {
948  prec_->setZeroStartingSolution(InitialGuessIsZero);
949  supportInitialGuess = true;
950  } else if (type_ == "SCHWARZ") {
951  paramList.set("schwarz: zero starting solution", InitialGuessIsZero);
952  // Because additive Schwarz has "delta" semantics, it's sufficient to
953  // toggle only the zero initial guess flag, and not pass in already
954  // set parameters. If we call SetPrecParameters, the subdomain solver
955  // will be destroyed.
956  prec_->setParameters(paramList);
957  supportInitialGuess = true;
958  }
959 
960  // TODO JJH 30Apr2014 Calling SetPrecParameters(paramList) when the smoother
961  // is Ifpack2::AdditiveSchwarz::setParameterList() will destroy the subdomain
962  //(aka inner) solver. This behavior is documented but a departure from what
963  // it previously did, and what other Ifpack2 solvers currently do. So I have
964  // moved SetPrecParameters(paramList) into the if-else block above.
965 
966  // Apply
967  if (InitialGuessIsZero || supportInitialGuess) {
970  prec_->apply(tpB, tpX);
971  } else {
972  typedef Teuchos::ScalarTraits<Scalar> TST;
973 
974  RCP<MultiVector> Residual;
975  {
976  std::string prefix = this->ShortClassName() + ": Apply: ";
977  RCP<TimeMonitor> tM = rcp(new TimeMonitor(*this, prefix + "residual calculation", Timings0));
978  Residual = Utilities::Residual(*A_, X, B);
979  }
980 
981  RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
982 
985 
986  prec_->apply(tpB, tpX);
987 
988  X.update(TST::one(), *Correction, TST::one());
989  }
990 }
991 
992 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
994  RCP<Ifpack2Smoother> smoother = rcp(new Ifpack2Smoother(*this));
995  smoother->SetParameterList(this->GetParameterList());
996  return smoother;
997 }
998 
999 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1001  std::ostringstream out;
1003  out << prec_->description();
1004  } else {
1006  out << "{type = " << type_ << "}";
1007  }
1008  return out.str();
1009 }
1010 
1011 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1014 
1015  if (verbLevel & Parameters0)
1016  out0 << "Prec. type: " << type_ << std::endl;
1017 
1018  if (verbLevel & Parameters1) {
1019  out0 << "Parameter list: " << std::endl;
1020  Teuchos::OSTab tab2(out);
1021  out << this->GetParameterList();
1022  out0 << "Overlap: " << overlap_ << std::endl;
1023  }
1024 
1025  if (verbLevel & External)
1026  if (prec_ != Teuchos::null) {
1027  Teuchos::OSTab tab2(out);
1028  out << *prec_ << std::endl;
1029  }
1030 
1031  if (verbLevel & Debug) {
1032  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
1033  << "-" << std::endl
1034  << "RCP<prec_>: " << prec_ << std::endl;
1035  }
1036 }
1037 
1038 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1040  typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
1041  // NOTE: Only works for a subset of Ifpack2's smoothers
1043  if (!pr.is_null()) return pr->getNodeSmootherComplexity();
1044 
1046  if (!pc.is_null()) return pc->getNodeSmootherComplexity();
1047 
1049  if (!pb.is_null()) return pb->getNodeSmootherComplexity();
1050 
1051  RCP<Ifpack2::ILUT<MatrixType>> pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType>>(prec_);
1052  if (!pi.is_null()) return pi->getNodeSmootherComplexity();
1053 
1054  RCP<Ifpack2::RILUK<MatrixType>> pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType>>(prec_);
1055  if (!pk.is_null()) return pk->getNodeSmootherComplexity();
1056 
1058 }
1059 
1060 } // namespace MueLu
1061 
1062 #endif // HAVE_MUELU_IFPACK2
1063 #endif // MUELU_IFPACK2SMOOTHER_DEF_HPP
void DeclareInput(Level &currentLevel) const
Input.
Important warning messages (one line)
void SetupGeneric(Level &currentLevel)
Exception indicating invalid cast attempted.
High level timing information (use Teuchos::TimeMonitor::summarize() to print)
MueLu::DefaultLocalOrdinal LocalOrdinal
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level-&gt;Get&lt; RCP&lt;Matrix&gt; &gt;(&quot;A&quot;, factory) if factory == NULL =&gt; use default factory.
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.
KOKKOS_INLINE_FUNCTION LO GetNumAggregates() const
static magnitudeType eps()
Print external lib objects.
GlobalOrdinal GO
T & get(const std::string &name, T def_value)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >> Op)
friend class Ifpack2Smoother
Constructor.
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Print more statistics.
Print additional debugging information.
LocalOrdinal LO
One-liner description of what is happening.
void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
size_type size() const
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
std::string type_
ifpack2-specific key phrase that denote smoother type
void SetupBlockRelaxation(Level &currentLevel)
void resize(const size_type n, const T &val=T())
bool isParameter(const std::string &name) const
void SetupSchwarz(Level &currentLevel)
void SetupLineSmoothing(Level &currentLevel)
Print statistics that do not involve significant additional computation.
bool remove(std::string const &name, bool throwIfNotExists=true)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
virtual void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
RCP< SmootherPrototype > Copy() const
void declareConstructionOutcome(bool fail, std::string msg)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
bool isSublist(const std::string &name) const
void Setup(Level &currentLevel)
Set up the smoother.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
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.
void SetupHiptmair(Level &currentLevel)
void SetupTopological(Level &currentLevel)
std::string description() const
Return a simple one-line description of this object.
bool IsSetup() const
Get the state of a smoother prototype.
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)
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
ParameterList & setParameters(const ParameterList &source)
Print all timing information.
Print class parameters.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< Level > & GetPreviousLevel()
Previous level.
Definition: MueLu_Level.cpp:60
MatrixType::scalar_type getLambdaMaxForApply() const
Scalar SC
void SetupAggregate(Level &currentLevel)
bool isType(const std::string &name) const
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Vector() const
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
Node NO
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > convertToBlockCrsMatrix(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize, bool use_local_ID=true)
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:51
Print class parameters (more parameters, more verbose)
Exception throws to report errors in the internal logical of the program.
void SetupChebyshev(Level &currentLevel)
#define TEUCHOS_ASSERT(assertion_test)
Description of what is happening (more verbose)
Class that encapsulates Ifpack2 smoothers.
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< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> vec)
virtual std::string description() const
Return a simple one-line description of this object.
std::string toString(const T &t)
bool is_null() const
Scalar SetupChebyshevEigenvalues(Level &currentLevel, const std::string &matrixName, const std::string &label, Teuchos::ParameterList &paramList) const