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