MueLu  Version of the Day
MueLu_Ifpack2Smoother_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_IFPACK2SMOOTHER_DEF_HPP
47 #define MUELU_IFPACK2SMOOTHER_DEF_HPP
48 
49 #include "MueLu_ConfigDefs.hpp"
50 
51 #if defined(HAVE_MUELU_TPETRA) && defined(HAVE_MUELU_IFPACK2)
52 
53 #include <Teuchos_ParameterList.hpp>
54 
55 #include <Tpetra_RowMatrix.hpp>
56 
57 #include <Ifpack2_Chebyshev.hpp>
58 #include <Ifpack2_Relaxation.hpp>
59 #include <Ifpack2_ILUT.hpp>
60 #include <Ifpack2_BlockRelaxation.hpp>
61 #include <Ifpack2_Factory.hpp>
62 #include <Ifpack2_Parameters.hpp>
63 
64 #include <Xpetra_BlockedCrsMatrix.hpp>
65 #include <Xpetra_CrsMatrix.hpp>
66 #include <Xpetra_CrsMatrixWrap.hpp>
67 #include <Xpetra_Matrix.hpp>
68 #include <Xpetra_MultiVectorFactory.hpp>
69 
71 #include "MueLu_Level.hpp"
73 #include "MueLu_Utilities.hpp"
74 #include "MueLu_Monitor.hpp"
75 
76 #ifdef HAVE_MUELU_INTREPID2
79  #include "Intrepid2_Basis.hpp"
80 
81  #ifdef HAVE_MUELU_INTREPID2_REFACTOR
82  #include "Kokkos_DynRankView.hpp"
83  #else
84  #include "Intrepid2_FieldContainer.hpp"
85  #endif
86 #endif
87 
88 // #define IFPACK2_HAS_PROPER_REUSE
89 
90 namespace MueLu {
91 
92  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
93  Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Ifpack2Smoother(const std::string& type, const Teuchos::ParameterList& paramList, const LO& overlap)
94  : type_(type), overlap_(overlap)
95  {
96  SetParameterList(paramList);
97  }
98 
99  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
101  Factory::SetParameterList(paramList);
102 
104  // It might be invalid to change parameters after the setup, but it depends entirely on Ifpack implementation.
105  // TODO: I don't know if Ifpack returns an error code or exception or ignore parameters modification in this case...
107  }
108  }
109 
110  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
112  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
113  paramList.setParameters(list);
114 
115  RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
116 
117  prec_->setParameters(*precList);
118 
119  paramList.setParameters(*precList); // what about that??
120  }
121 
122  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
124  this->Input(currentLevel, "A");
125 
126  if (type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
127  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
128  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
129  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
130  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
131  type_ == "LINESMOOTHING_BLOCKRELAXATION") {
132  this->Input(currentLevel, "CoarseNumZLayers"); // necessary for fallback criterion
133  this->Input(currentLevel, "LineDetection_VertLineIds"); // necessary to feed block smoother
134  }
135  else if (type_ == "TOPOLOGICAL")
136  {
137  // for the topological smoother, we require an element to node map:
138  this->Input(currentLevel, "pcoarsen: element to node map");
139  }
140  }
141 
142  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
144  FactoryMonitor m(*this, "Setup Smoother", currentLevel);
145 
146  A_ = Factory::Get< RCP<Matrix> >(currentLevel, "A");
147 
148  if (type_ == "SCHWARZ")
149  SetupSchwarz(currentLevel);
150 
151  else if (type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
152  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
153  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
154  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
155  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
156  type_ == "LINESMOOTHING_BLOCKRELAXATION")
157  SetupLineSmoothing(currentLevel);
158 
159  else if (type_ == "CHEBYSHEV")
160  SetupChebyshev(currentLevel);
161 
162  else if (type_ == "TOPOLOGICAL")
163  {
164 #ifdef HAVE_MUELU_INTREPID2
165  SetupTopological(currentLevel);
166 #else
167  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, "'TOPOLOGICAL' smoother choice requires Intrepid2");
168 #endif
169  }
170  else
171  {
172  SetupGeneric(currentLevel);
173  }
174 
176 
177  this->GetOStream(Statistics1) << description() << std::endl;
178  }
179 
180  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
182  typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
183 
184  bool reusePreconditioner = false;
185  if (this->IsSetup() == true) {
186  // Reuse the constructed preconditioner
187  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
188 
189  bool isTRowMatrix = true;
190  RCP<const tRowMatrix> tA;
191  try {
193  } catch (Exceptions::BadCast) {
194  isTRowMatrix = false;
195  }
196 
197  RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
198  if (!prec.is_null() && isTRowMatrix) {
199 #ifdef IFPACK2_HAS_PROPER_REUSE
200  prec->resetMatrix(tA);
201  reusePreconditioner = true;
202 #else
203  this->GetOStream(Errors) << "Ifpack2 does not have proper reuse yet." << std::endl;
204 #endif
205 
206  } else {
207  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
208  "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
209  }
210  }
211 
212  if (!reusePreconditioner) {
213  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
214 
215  bool isBlockedMatrix = false;
216  RCP<Matrix> merged2Mat;
217 
218  std::string sublistName = "subdomain solver parameters";
219  if (paramList.isSublist(sublistName)) {
220  // If we are doing "user" partitioning, we assume that what the user
221  // really wants to do is make tiny little subdomains with one row
222  // assigned to each subdomain. The rows used for these little
223  // subdomains correspond to those in the 2nd block row. Then,
224  // if we overlap these mini-subdomains, we will do something that
225  // looks like Vanka (grabbing all velocities associated with each
226  // each pressure unknown). In addition, we put all Dirichlet points
227  // as a little mini-domain.
228  ParameterList& subList = paramList.sublist(sublistName);
229 
230  std::string partName = "partitioner: type";
231  if (subList.isParameter(partName) && subList.get<std::string>(partName) == "user") {
232  isBlockedMatrix = true;
233 
234  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
235  TEUCHOS_TEST_FOR_EXCEPTION(bA.is_null(), Exceptions::BadCast,
236  "Matrix A must be of type BlockedCrsMatrix.");
237 
238  size_t numVels = bA->getMatrix(0,0)->getNodeNumRows();
239  size_t numPres = bA->getMatrix(1,0)->getNodeNumRows();
240  size_t numRows = A_->getNodeNumRows();
241 
242  ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
243 
244  size_t numBlocks = 0;
245  for (size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
246  blockSeeds[rowOfB] = numBlocks++;
247 
248  RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
249  TEUCHOS_TEST_FOR_EXCEPTION(bA2.is_null(), Exceptions::BadCast,
250  "Matrix A must be of type BlockedCrsMatrix.");
251 
252  merged2Mat = bA2->Merge();
253 
254  // Add Dirichlet rows to the list of seeds
255  ArrayRCP<const bool> boundaryNodes = Utilities::DetectDirichletRows(*merged2Mat, 0.0);
256  bool haveBoundary = false;
257  for (LO i = 0; i < boundaryNodes.size(); i++)
258  if (boundaryNodes[i]) {
259  // FIXME:
260  // 1. would not this [] overlap with some in the previos blockSeed loop?
261  // 2. do we need to distinguish between pressure and velocity Dirichlet b.c.
262  blockSeeds[i] = numBlocks;
263  haveBoundary = true;
264  }
265  if (haveBoundary)
266  numBlocks++;
267 
268  subList.set("partitioner: map", blockSeeds);
269  subList.set("partitioner: local parts", as<int>(numBlocks));
270 
271  } else {
272  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
273  if (!bA.is_null()) {
274  isBlockedMatrix = true;
275  merged2Mat = bA->Merge();
276  }
277  }
278  }
279 
280  RCP<const tRowMatrix> tA;
281  if (isBlockedMatrix == true) tA = Utilities::Op2NonConstTpetraRow(merged2Mat);
283 
286 
287  prec_->initialize();
288  }
289 
290  prec_->compute();
291  }
292 
293 #ifdef HAVE_MUELU_INTREPID2
294  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
296  /*
297 
298  basic notion:
299 
300  Look for user input indicating topo dimension, something like "topological domain type: {node|edge|face}"
301  Call something like what you can find in Poisson example line 1180 to set seeds for a smoother
302 
303  */
304  if (this->IsSetup() == true) {
305  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
306  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
307  }
308 
309  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
310 
311  typedef typename Node::device_type::execution_space ES;
312 
313 #ifdef HAVE_MUELU_INTREPID2_REFACTOR
314  typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO; // "Field Container" for ordinals
315 #else
316  typedef Intrepid2::FieldContainer<LocalOrdinal> FCO;
317 #endif
318 
319  LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
320 
321  using namespace std;
322 
323  const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO> >(currentLevel,"pcoarsen: element to node map");
324 
325  string basisString = paramList.get<string>("pcoarsen: hi basis");
326  int degree;
327  // NOTE: To make sure Stokhos works we only instantiate these guys with double. There's a lot
328  // of stuff in the guts of Intrepid2 that doesn't play well with Stokhos as of yet. Here, we only
329  // care about the assignment of basis ordinals to topological entities, so this code is actually
330  // independent of the Scalar type--hard-coding double here won't hurt us.
331 #ifdef HAVE_MUELU_INTREPID2_REFACTOR
332  auto basis = MueLuIntrepid::BasisFactory<double,ES>(basisString, degree);
333 #else
334  auto basis = MueLuIntrepid::BasisFactory<double>(basisString, degree);
335 #endif
336 
337  string topologyTypeString = paramList.get<string>("smoother: neighborhood type");
338  int dimension;
339  if (topologyTypeString == "node")
340  dimension = 0;
341  else if (topologyTypeString == "edge")
342  dimension = 1;
343  else if (topologyTypeString == "face")
344  dimension = 2;
345  else if (topologyTypeString == "cell")
346  dimension = basis->getBaseCellTopology().getDimension();
347  else
348  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
349  vector<vector<LocalOrdinal>> seeds;
350  MueLuIntrepid::FindGeometricSeedOrdinals(basis, *elemToNode, seeds, *A_->getRowMap(), *A_->getColMap());
351 
352  // Ifpack2 wants the seeds in an array of the same length as the number of local elements,
353  // with local partition #s marked for the ones that are seeds, and invalid for the rest
354  int myNodeCount = A_->getRowMap()->getNodeNumElements();
355  ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
356  int localPartitionNumber = 0;
357  for (LocalOrdinal seed : seeds[dimension])
358  {
359  nodeSeeds[seed] = localPartitionNumber++;
360  }
361 
362  paramList.remove("smoother: neighborhood type");
363  paramList.remove("pcoarsen: hi basis");
364 
365  paramList.set("partitioner: map", nodeSeeds);
366  paramList.set("partitioner: type", "user");
367  paramList.set("partitioner: overlap", 1);
368  paramList.set("partitioner: local parts", int(seeds[dimension].size()));
369 
370  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
371 
372  type_ = "BLOCKRELAXATION";
375  prec_->initialize();
376  prec_->compute();
377  }
378 #endif
379 
380  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
382  if (this->IsSetup() == true) {
383  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
384  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
385  }
386 
387  ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
388 
389  LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,"CoarseNumZLayers");
390  if (CoarseNumZLayers > 0) {
391  Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get< Teuchos::ArrayRCP<LO> >(currentLevel, "LineDetection_VertLineIds");
392 
393  // determine number of local parts
394  LO maxPart = 0;
395  for(size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
396  if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
397  }
398 
399  size_t numLocalRows = A_->getNodeNumRows();
400  TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % TVertLineIdSmoo.size() != 0, Exceptions::RuntimeError,
401  "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
402 
403  if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
404  myparamList.set("partitioner: type","user");
405  myparamList.set("partitioner: map",TVertLineIdSmoo);
406  myparamList.set("partitioner: local parts",maxPart+1);
407  } else {
408  // we assume a constant number of DOFs per node
409  size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
410 
411  // Create a new Teuchos::ArrayRCP<LO> of size numLocalRows and fill it with the corresponding information
412  Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
413  for (size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
414  for (size_t dof = 0; dof < numDofsPerNode; dof++)
415  partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
416  myparamList.set("partitioner: type","user");
417  myparamList.set("partitioner: map",partitionerMap);
418  myparamList.set("partitioner: local parts",maxPart + 1);
419  }
420 
421  if (type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
422  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
423  type_ == "LINESMOOTHING_BANDEDRELAXATION")
424  type_ = "BANDEDRELAXATION";
425  else
426  type_ = "BLOCKRELAXATION";
427  } else {
428  // line detection failed -> fallback to point-wise relaxation
429  this->GetOStream(Runtime0) << "Line detection failed: fall back to point-wise relaxation" << std::endl;
430  myparamList.remove("partitioner: type",false);
431  myparamList.remove("partitioner: map", false);
432  myparamList.remove("partitioner: local parts",false);
433  type_ = "RELAXATION";
434  }
435 
436  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
437 
440  prec_->initialize();
441  prec_->compute();
442  }
443 
444  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
446  if (this->IsSetup() == true) {
447  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): SetupChebyshev() has already been called" << std::endl;
448  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available, reverting to full construction" << std::endl;
449  }
450 
451  typedef Teuchos::ScalarTraits<SC> STS;
452  SC negone = -STS::one();
453 
454  SC lambdaMax = negone;
455  {
456  std::string maxEigString = "chebyshev: max eigenvalue";
457  std::string eigRatioString = "chebyshev: ratio eigenvalue";
458 
459  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
460 
461  // Get/calculate the maximum eigenvalue
462  if (paramList.isParameter(maxEigString)) {
463  if (paramList.isType<double>(maxEigString))
464  lambdaMax = paramList.get<double>(maxEigString);
465  else
466  lambdaMax = paramList.get<SC>(maxEigString);
467  this->GetOStream(Statistics1) << maxEigString << " (cached with smoother parameter list) = " << lambdaMax << std::endl;
468 
469  } else {
470  lambdaMax = A_->GetMaxEigenvalueEstimate();
471  if (lambdaMax != negone) {
472  this->GetOStream(Statistics1) << maxEigString << " (cached with matrix) = " << lambdaMax << std::endl;
473  paramList.set(maxEigString, lambdaMax);
474  }
475  }
476 
477  // Calculate the eigenvalue ratio
478  const SC defaultEigRatio = 20;
479 
480  SC ratio = defaultEigRatio;
481  if (paramList.isParameter(eigRatioString)) {
482  if (paramList.isType<double>(eigRatioString))
483  ratio = paramList.get<double>(eigRatioString);
484  else
485  ratio = paramList.get<SC>(eigRatioString);
486  }
487  if (currentLevel.GetLevelID()) {
488  // Update ratio to be
489  // ratio = max(number of fine DOFs / number of coarse DOFs, defaultValue)
490  //
491  // NOTE: We don't need to request previous level matrix as we know for sure it was constructed
492  RCP<const Matrix> fineA = currentLevel.GetPreviousLevel()->Get<RCP<Matrix> >("A");
493  size_t nRowsFine = fineA->getGlobalNumRows();
494  size_t nRowsCoarse = A_->getGlobalNumRows();
495 
496  SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
497  if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
498  ratio = levelRatio;
499  }
500 
501  this->GetOStream(Statistics1) << eigRatioString << " (computed) = " << ratio << std::endl;
502  paramList.set(eigRatioString, ratio);
503  }
504 
505  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
506 
509  prec_->initialize();
510  prec_->compute();
511 
512  if (lambdaMax == negone) {
513  typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
514 
515  Teuchos::RCP<Ifpack2::Chebyshev<MatrixType> > chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
516  if (chebyPrec != Teuchos::null) {
517  lambdaMax = chebyPrec->getLambdaMaxForApply();
518  A_->SetMaxEigenvalueEstimate(lambdaMax);
519  this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)" << " = " << lambdaMax << std::endl;
520  }
521  TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
522  }
523  }
524 
525  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
527  typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
528 
529  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
530  if (!bA.is_null())
531  A_ = bA->Merge();
532 
533  RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
534 
535  bool reusePreconditioner = false;
536  if (this->IsSetup() == true) {
537  // Reuse the constructed preconditioner
538  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
539 
540  RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
541  if (!prec.is_null()) {
542 #ifdef IFPACK2_HAS_PROPER_REUSE
543  prec->resetMatrix(tA);
544  reusePreconditioner = true;
545 #else
546  this->GetOStream(Errors) << "Ifpack2 does not have proper reuse yet." << std::endl;
547 #endif
548 
549  } else {
550  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available (failed cast to CanChangeMatrix), "
551  "reverting to full construction" << std::endl;
552  }
553  }
554 
555  if (!reusePreconditioner) {
558  prec_->initialize();
559  }
560 
561  prec_->compute();
562  }
563 
564  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
565  void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
566  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Apply(): Setup() has not been called");
567 
568  // Forward the InitialGuessIsZero option to Ifpack2
569  // TODO: It might be nice to switch back the internal
570  // "zero starting solution" option of the ifpack2 object prec_ to his
571  // initial value at the end but there is no way right now to get
572  // the current value of the "zero starting solution" in ifpack2.
573  // It's not really an issue, as prec_ can only be used by this method.
574  // TODO: When https://software.sandia.gov/bugzilla/show_bug.cgi?id=5283#c2 is done
575  // we should remove the if/else/elseif and just test if this
576  // option is supported by current ifpack2 preconditioner
577  Teuchos::ParameterList paramList;
578  bool supportInitialGuess = false;
579  if (type_ == "CHEBYSHEV") {
580  paramList.set("chebyshev: zero starting solution", InitialGuessIsZero);
581  SetPrecParameters(paramList);
582  supportInitialGuess = true;
583 
584  } else if (type_ == "RELAXATION") {
585  paramList.set("relaxation: zero starting solution", InitialGuessIsZero);
586  SetPrecParameters(paramList);
587  supportInitialGuess = true;
588 
589  } else if (type_ == "KRYLOV") {
590  paramList.set("krylov: zero starting solution", InitialGuessIsZero);
591  SetPrecParameters(paramList);
592  supportInitialGuess = true;
593 
594  } else if (type_ == "SCHWARZ") {
595  paramList.set("schwarz: zero starting solution", InitialGuessIsZero);
596  //Because additive Schwarz has "delta" semantics, it's sufficient to
597  //toggle only the zero initial guess flag, and not pass in already
598  //set parameters. If we call SetPrecParameters, the subdomain solver
599  //will be destroyed.
600  prec_->setParameters(paramList);
601  supportInitialGuess = true;
602  }
603 
604  //TODO JJH 30Apr2014 Calling SetPrecParameters(paramList) when the smoother
605  //is Ifpack2::AdditiveSchwarz::setParameterList() will destroy the subdomain
606  //(aka inner) solver. This behavior is documented but a departure from what
607  //it previously did, and what other Ifpack2 solvers currently do. So I have
608  //moved SetPrecParameters(paramList) into the if-else block above.
609 
610  // Apply
611  if (InitialGuessIsZero || supportInitialGuess) {
612  Tpetra::MultiVector<SC,LO,GO,NO>& tpX = Utilities::MV2NonConstTpetraMV(X);
613  const Tpetra::MultiVector<SC,LO,GO,NO>& tpB = Utilities::MV2TpetraMV(B);
614  prec_->apply(tpB, tpX);
615  } else {
616  typedef Teuchos::ScalarTraits<Scalar> TST;
617  RCP<MultiVector> Residual = Utilities::Residual(*A_, X, B);
618  RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
619 
620  Tpetra::MultiVector<SC,LO,GO,NO>& tpX = Utilities::MV2NonConstTpetraMV(*Correction);
621  const Tpetra::MultiVector<SC,LO,GO,NO>& tpB = Utilities::MV2TpetraMV(*Residual);
622 
623  prec_->apply(tpB, tpX);
624 
625  X.update(TST::one(), *Correction, TST::one());
626  }
627  }
628 
629  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
630  RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
631  RCP<Ifpack2Smoother> smoother = rcp(new Ifpack2Smoother(*this) );
632  smoother->SetParameterList(this->GetParameterList());
633  return smoother;
634  }
635 
636  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
638  std::ostringstream out;
640  out << prec_->description();
641  } else {
643  out << "{type = " << type_ << "}";
644  }
645  return out.str();
646  }
647 
648  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
649  void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream &out, const VerbLevel verbLevel) const {
651 
652  if (verbLevel & Parameters0)
653  out0 << "Prec. type: " << type_ << std::endl;
654 
655  if (verbLevel & Parameters1) {
656  out0 << "Parameter list: " << std::endl;
657  Teuchos::OSTab tab2(out);
658  out << this->GetParameterList();
659  out0 << "Overlap: " << overlap_ << std::endl;
660  }
661 
662  if (verbLevel & External)
663  if (prec_ != Teuchos::null) {
664  Teuchos::OSTab tab2(out);
665  out << *prec_ << std::endl;
666  }
667 
668  if (verbLevel & Debug) {
669  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
670  << "-" << std::endl
671  << "RCP<prec_>: " << prec_ << std::endl;
672  }
673  }
674 
675  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
677  typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
678  // NOTE: Only works for a subset of Ifpack2's smoothers
679  RCP<Ifpack2::Relaxation<MatrixType> > pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType> >(prec_);
680  if(!pr.is_null()) return pr->getNodeSmootherComplexity();
681 
682  RCP<Ifpack2::Chebyshev<MatrixType> > pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
683  if(!pc.is_null()) return pc->getNodeSmootherComplexity();
684 
685  RCP<Ifpack2::BlockRelaxation<MatrixType> > pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType> >(prec_);
686  if(!pb.is_null()) return pb->getNodeSmootherComplexity();
687 
688  RCP<Ifpack2::ILUT<MatrixType> > pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType> >(prec_);
689  if(!pi.is_null()) return pi->getNodeSmootherComplexity();
690 
691 
692  return Teuchos::OrdinalTraits<size_t>::invalid();
693  }
694 
695 
696 } // namespace MueLu
697 
698 #endif // HAVE_MUELU_TPETRA && HAVE_MUELU_IFPACK2
699 #endif // MUELU_IFPACK2SMOOTHER_DEF_HPP
Important warning messages (one line)
void SetupGeneric(Level &currentLevel)
RCP< SmootherPrototype > Copy() const
Exception indicating invalid cast attempted.
RCP< Level > & GetPreviousLevel()
Previous level.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
size_t getNodeSmootherComplexity() const
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
std::string description() const
Return a simple one-line description of this object.
Print external lib objects.
bool IsSetup() const
Get the state of a smoother prototype.
friend class Ifpack2Smoother
Constructor.
Timer to be used in factories. Similar to Monitor but with additional timers.
RCP< Matrix > A_
matrix, used in apply if solving residual equation
RCP< ParameterList > RemoveFactoriesFromList(const ParameterList &list) const
Print more statistics.
Print additional debugging information.
One-liner description of what is happening.
void SetParameterList(const Teuchos::ParameterList &paramList)
Namespace for MueLu classes and methods.
LO overlap_
overlap when using the smoother in additive Schwarz mode
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
std::string type_
ifpack2-specific key phrase that denote smoother type
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetupSchwarz(Level &currentLevel)
void SetupLineSmoothing(Level &currentLevel)
MatrixType::scalar_type getLambdaMaxForApply() const
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS)
virtual void SetParameterList(const ParameterList &paramList)
Set parameters from a parameter list and return with default values.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void Setup(Level &currentLevel)
Set up the smoother.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
void SetupTopological(Level &currentLevel)
size_t getNodeSmootherComplexity() const
static Teuchos::RCP< Preconditioner< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > create(const std::string &precType, const Teuchos::RCP< const MatrixType > &matrix)
void DeclareInput(Level &currentLevel) const
Input.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
virtual const Teuchos::ParameterList & GetParameterList() const
Print class parameters.
void Input(Level &level, const std::string &varName) const
size_t getNodeSmootherComplexity() const
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::magnitude(0.))
Print class parameters (more parameters, more verbose)
size_t getNodeSmootherComplexity() const
Exception throws to report errors in the internal logical of the program.
void SetupChebyshev(Level &currentLevel)
Description of what is happening (more verbose)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
virtual std::string description() const
Return a simple one-line description of this object.
RCP< Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node > > prec_
pointer to Ifpack2 preconditioner object