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 #include <Xpetra_TpetraMultiVector.hpp>
70 
72 #include "MueLu_Level.hpp"
74 #include "MueLu_Utilities.hpp"
75 #include "MueLu_Monitor.hpp"
76 
77 #ifdef HAVE_MUELU_INTREPID2
80 #include "Intrepid2_Basis.hpp"
81 #include "Kokkos_DynRankView.hpp"
82 #endif
83 
84 // #define IFPACK2_HAS_PROPER_REUSE
85 
86 namespace MueLu {
87 
88  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
89  Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Ifpack2Smoother(const std::string& type, const Teuchos::ParameterList& paramList, const LO& overlap)
90  : type_(type), overlap_(overlap)
91  {
92  SetParameterList(paramList);
93  }
94 
95  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
97  Factory::SetParameterList(paramList);
98 
100  // It might be invalid to change parameters after the setup, but it depends entirely on Ifpack implementation.
101  // TODO: I don't know if Ifpack returns an error code or exception or ignore parameters modification in this case...
102  SetPrecParameters();
103  }
104  }
105 
106  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
108  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
109  paramList.setParameters(list);
110 
111  RCP<ParameterList> precList = this->RemoveFactoriesFromList(this->GetParameterList());
112 
113  prec_->setParameters(*precList);
114 
115  paramList.setParameters(*precList); // what about that??
116  }
117 
118  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
120  this->Input(currentLevel, "A");
121 
122  if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
123  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
124  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
125  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
126  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
127  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
128  type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
129  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
130  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
131  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
132  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
133  type_ == "LINESMOOTHING_BLOCKRELAXATION") {
134  this->Input(currentLevel, "CoarseNumZLayers"); // necessary for fallback criterion
135  this->Input(currentLevel, "LineDetection_VertLineIds"); // necessary to feed block smoother
136  }
137  else if (type_ == "BLOCK RELAXATION" ||
138  type_ == "BLOCK_RELAXATION" ||
139  type_ == "BLOCKRELAXATION" ||
140  // Banded
141  type_ == "BANDED_RELAXATION" ||
142  type_ == "BANDED RELAXATION" ||
143  type_ == "BANDEDRELAXATION" ||
144  // Tridiagonal
145  type_ == "TRIDI_RELAXATION" ||
146  type_ == "TRIDI RELAXATION" ||
147  type_ == "TRIDIRELAXATION" ||
148  type_ == "TRIDIAGONAL_RELAXATION" ||
149  type_ == "TRIDIAGONAL RELAXATION" ||
150  type_ == "TRIDIAGONALRELAXATION")
151  {
152  //We need to check for the "partitioner type" = "line"
153  ParameterList precList = this->GetParameterList();
154  if(precList.isParameter("partitioner: type") &&
155  precList.get<std::string>("partitioner: type") == "line") {
156  this->Input(currentLevel, "Coordinates");
157  }
158  }
159  else if (type_ == "TOPOLOGICAL")
160  {
161  // for the topological smoother, we require an element to node map:
162  this->Input(currentLevel, "pcoarsen: element to node map");
163  }
164  }
165 
166  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
168  FactoryMonitor m(*this, "Setup Smoother", currentLevel);
169 
170  A_ = Factory::Get< RCP<Matrix> >(currentLevel, "A");
171 
172  if (type_ == "SCHWARZ")
173  SetupSchwarz(currentLevel);
174 
175  else if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
176  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
177  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
178  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
179  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
180  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION" ||
181  type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
182  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
183  type_ == "LINESMOOTHING_BANDEDRELAXATION" ||
184  type_ == "LINESMOOTHING_BLOCK_RELAXATION" ||
185  type_ == "LINESMOOTHING_BLOCK RELAXATION" ||
186  type_ == "LINESMOOTHING_BLOCKRELAXATION")
187  SetupLineSmoothing(currentLevel);
188 
189  else if (type_ == "BLOCK_RELAXATION" ||
190  type_ == "BLOCK RELAXATION" ||
191  type_ == "BLOCKRELAXATION" ||
192  // Banded
193  type_ == "BANDED_RELAXATION" ||
194  type_ == "BANDED RELAXATION" ||
195  type_ == "BANDEDRELAXATION" ||
196  // Tridiagonal
197  type_ == "TRIDI_RELAXATION" ||
198  type_ == "TRIDI RELAXATION" ||
199  type_ == "TRIDIRELAXATION" ||
200  type_ == "TRIDIAGONAL_RELAXATION" ||
201  type_ == "TRIDIAGONAL RELAXATION" ||
202  type_ == "TRIDIAGONALRELAXATION")
203  SetupBlockRelaxation(currentLevel);
204 
205  else if (type_ == "CHEBYSHEV")
206  SetupChebyshev(currentLevel);
207 
208  else if (type_ == "TOPOLOGICAL")
209  {
210 #ifdef HAVE_MUELU_INTREPID2
211  SetupTopological(currentLevel);
212 #else
213  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, "'TOPOLOGICAL' smoother choice requires Intrepid2");
214 #endif
215  }
216  else
217  {
218  SetupGeneric(currentLevel);
219  }
220 
222 
223  this->GetOStream(Statistics1) << description() << std::endl;
224  }
225 
226  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
228  typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
229 
230  bool reusePreconditioner = false;
231  if (this->IsSetup() == true) {
232  // Reuse the constructed preconditioner
233  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupSchwarz(): Setup() has already been called, assuming reuse" << std::endl;
234 
235  bool isTRowMatrix = true;
236  RCP<const tRowMatrix> tA;
237  try {
239  } catch (Exceptions::BadCast&) {
240  isTRowMatrix = false;
241  }
242 
243  RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
244  if (!prec.is_null() && isTRowMatrix) {
245 #ifdef IFPACK2_HAS_PROPER_REUSE
246  prec->resetMatrix(tA);
247  reusePreconditioner = true;
248 #else
249  this->GetOStream(Errors) << "Ifpack2 does not have proper reuse yet." << std::endl;
250 #endif
251 
252  } else {
253  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupSchwarz(): reuse of this type is not available "
254  "(either failed cast to CanChangeMatrix, or to Tpetra Row Matrix), reverting to full construction" << std::endl;
255  }
256  }
257 
258  if (!reusePreconditioner) {
259  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
260 
261  bool isBlockedMatrix = false;
262  RCP<Matrix> merged2Mat;
263 
264  std::string sublistName = "subdomain solver parameters";
265  if (paramList.isSublist(sublistName)) {
266  // If we are doing "user" partitioning, we assume that what the user
267  // really wants to do is make tiny little subdomains with one row
268  // assigned to each subdomain. The rows used for these little
269  // subdomains correspond to those in the 2nd block row. Then,
270  // if we overlap these mini-subdomains, we will do something that
271  // looks like Vanka (grabbing all velocities associated with each
272  // each pressure unknown). In addition, we put all Dirichlet points
273  // as a little mini-domain.
274  ParameterList& subList = paramList.sublist(sublistName);
275 
276  std::string partName = "partitioner: type";
277  if (subList.isParameter(partName) && subList.get<std::string>(partName) == "user") {
278  isBlockedMatrix = true;
279 
280  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
281  TEUCHOS_TEST_FOR_EXCEPTION(bA.is_null(), Exceptions::BadCast,
282  "Matrix A must be of type BlockedCrsMatrix.");
283 
284  size_t numVels = bA->getMatrix(0,0)->getNodeNumRows();
285  size_t numPres = bA->getMatrix(1,0)->getNodeNumRows();
286  size_t numRows = A_->getNodeNumRows();
287 
288  ArrayRCP<LocalOrdinal> blockSeeds(numRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
289 
290  size_t numBlocks = 0;
291  for (size_t rowOfB = numVels; rowOfB < numVels+numPres; ++rowOfB)
292  blockSeeds[rowOfB] = numBlocks++;
293 
294  RCP<BlockedCrsMatrix> bA2 = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
295  TEUCHOS_TEST_FOR_EXCEPTION(bA2.is_null(), Exceptions::BadCast,
296  "Matrix A must be of type BlockedCrsMatrix.");
297 
298  merged2Mat = bA2->Merge();
299 
300  // Add Dirichlet rows to the list of seeds
301  ArrayRCP<const bool> boundaryNodes = Utilities::DetectDirichletRows(*merged2Mat, 0.0);
302  bool haveBoundary = false;
303  for (LO i = 0; i < boundaryNodes.size(); i++)
304  if (boundaryNodes[i]) {
305  // FIXME:
306  // 1. would not this [] overlap with some in the previos blockSeed loop?
307  // 2. do we need to distinguish between pressure and velocity Dirichlet b.c.
308  blockSeeds[i] = numBlocks;
309  haveBoundary = true;
310  }
311  if (haveBoundary)
312  numBlocks++;
313 
314  subList.set("partitioner: map", blockSeeds);
315  subList.set("partitioner: local parts", as<int>(numBlocks));
316 
317  } else {
318  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
319  if (!bA.is_null()) {
320  isBlockedMatrix = true;
321  merged2Mat = bA->Merge();
322  }
323  }
324  }
325 
326  RCP<const tRowMatrix> tA;
327  if (isBlockedMatrix == true) tA = Utilities::Op2NonConstTpetraRow(merged2Mat);
328  else tA = Utilities::Op2NonConstTpetraRow(A_);
329 
330  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
331  SetPrecParameters();
332 
333  prec_->initialize();
334  }
335 
336  prec_->compute();
337  }
338 
339 #ifdef HAVE_MUELU_INTREPID2
340  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
342  /*
343 
344  basic notion:
345 
346  Look for user input indicating topo dimension, something like "topological domain type: {node|edge|face}"
347  Call something like what you can find in Poisson example line 1180 to set seeds for a smoother
348 
349  */
350  if (this->IsSetup() == true) {
351  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): Setup() has already been called" << std::endl;
352  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupTopological(): reuse of this type is not available, reverting to full construction" << std::endl;
353  }
354 
355  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
356 
357  typedef typename Node::device_type::execution_space ES;
358 
359  typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCO; //
360 
361  LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
362 
363  using namespace std;
364 
365  const Teuchos::RCP<FCO> elemToNode = Factory::Get<Teuchos::RCP<FCO> >(currentLevel,"pcoarsen: element to node map");
366 
367  string basisString = paramList.get<string>("pcoarsen: hi basis");
368  int degree;
369  // NOTE: To make sure Stokhos works we only instantiate these guys with double. There's a lot
370  // of stuff in the guts of Intrepid2 that doesn't play well with Stokhos as of yet. Here, we only
371  // care about the assignment of basis ordinals to topological entities, so this code is actually
372  // independent of the Scalar type--hard-coding double here won't hurt us.
373  auto basis = MueLuIntrepid::BasisFactory<double,ES>(basisString, degree);
374 
375  string topologyTypeString = paramList.get<string>("smoother: neighborhood type");
376  int dimension;
377  if (topologyTypeString == "node")
378  dimension = 0;
379  else if (topologyTypeString == "edge")
380  dimension = 1;
381  else if (topologyTypeString == "face")
382  dimension = 2;
383  else if (topologyTypeString == "cell")
384  dimension = basis->getBaseCellTopology().getDimension();
385  else
386  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,"Unrecognized smoother neighborhood type. Supported types are node, edge, face.");
387  vector<vector<LocalOrdinal>> seeds;
388  MueLuIntrepid::FindGeometricSeedOrdinals(basis, *elemToNode, seeds, *A_->getRowMap(), *A_->getColMap());
389 
390  // Ifpack2 wants the seeds in an array of the same length as the number of local elements,
391  // with local partition #s marked for the ones that are seeds, and invalid for the rest
392  int myNodeCount = A_->getRowMap()->getNodeNumElements();
393  ArrayRCP<LocalOrdinal> nodeSeeds(myNodeCount,lo_invalid);
394  int localPartitionNumber = 0;
395  for (LocalOrdinal seed : seeds[dimension])
396  {
397  nodeSeeds[seed] = localPartitionNumber++;
398  }
399 
400  paramList.remove("smoother: neighborhood type");
401  paramList.remove("pcoarsen: hi basis");
402 
403  paramList.set("partitioner: map", nodeSeeds);
404  paramList.set("partitioner: type", "user");
405  paramList.set("partitioner: overlap", 1);
406  paramList.set("partitioner: local parts", int(seeds[dimension].size()));
407 
408  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
409 
410  type_ = "BLOCKRELAXATION";
411  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
412  SetPrecParameters();
413  prec_->initialize();
414  prec_->compute();
415  }
416 #endif
417 
418  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
420  if (this->IsSetup() == true) {
421  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): Setup() has already been called" << std::endl;
422  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupLineSmoothing(): reuse of this type is not available, reverting to full construction" << std::endl;
423  }
424 
425  ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
426 
427  LO CoarseNumZLayers = Factory::Get<LO>(currentLevel,"CoarseNumZLayers");
428  if (CoarseNumZLayers > 0) {
429  Teuchos::ArrayRCP<LO> TVertLineIdSmoo = Factory::Get< Teuchos::ArrayRCP<LO> >(currentLevel, "LineDetection_VertLineIds");
430 
431  // determine number of local parts
432  LO maxPart = 0;
433  for(size_t k = 0; k < Teuchos::as<size_t>(TVertLineIdSmoo.size()); k++) {
434  if(maxPart < TVertLineIdSmoo[k]) maxPart = TVertLineIdSmoo[k];
435  }
436  size_t numLocalRows = A_->getNodeNumRows();
437 
438  TEUCHOS_TEST_FOR_EXCEPTION(numLocalRows % TVertLineIdSmoo.size() != 0, Exceptions::RuntimeError,
439  "MueLu::Ifpack2Smoother::Setup(): the number of local nodes is incompatible with the TVertLineIdsSmoo.");
440 
441  //actualDofsPerNode is the actual number of matrix rows per mesh element.
442  //It is encoded in either the MueLu Level, or in the Xpetra matrix block size.
443  //This value is needed by Ifpack2 to do decoupled block relaxation.
444  int actualDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
445  LO matrixBlockSize = A_->GetFixedBlockSize();
446  if(matrixBlockSize > 1 && actualDofsPerNode > 1)
447  {
448  TEUCHOS_TEST_FOR_EXCEPTION(actualDofsPerNode != A_->GetFixedBlockSize(), Exceptions::RuntimeError,
449  "MueLu::Ifpack2Smoother::Setup(): A is a block matrix but its block size and DOFs/node from partitioner disagree");
450  }
451  else if(matrixBlockSize > 1)
452  {
453  actualDofsPerNode = A_->GetFixedBlockSize();
454  }
455  myparamList.set("partitioner: PDE equations", actualDofsPerNode);
456 
457  if (numLocalRows == Teuchos::as<size_t>(TVertLineIdSmoo.size())) {
458  myparamList.set("partitioner: type","user");
459  myparamList.set("partitioner: map",TVertLineIdSmoo);
460  myparamList.set("partitioner: local parts",maxPart+1);
461  } else {
462  // we assume a constant number of DOFs per node
463  size_t numDofsPerNode = numLocalRows / TVertLineIdSmoo.size();
464 
465  // Create a new Teuchos::ArrayRCP<LO> of size numLocalRows and fill it with the corresponding information
466  Teuchos::ArrayRCP<LO> partitionerMap(numLocalRows, Teuchos::OrdinalTraits<LocalOrdinal>::invalid());
467  for (size_t blockRow = 0; blockRow < Teuchos::as<size_t>(TVertLineIdSmoo.size()); ++blockRow)
468  for (size_t dof = 0; dof < numDofsPerNode; dof++)
469  partitionerMap[blockRow * numDofsPerNode + dof] = TVertLineIdSmoo[blockRow];
470  myparamList.set("partitioner: type","user");
471  myparamList.set("partitioner: map",partitionerMap);
472  myparamList.set("partitioner: local parts",maxPart + 1);
473  }
474 
475  if (type_ == "LINESMOOTHING_BANDED_RELAXATION" ||
476  type_ == "LINESMOOTHING_BANDED RELAXATION" ||
477  type_ == "LINESMOOTHING_BANDEDRELAXATION")
478  type_ = "BANDEDRELAXATION";
479  else if (type_ == "LINESMOOTHING_TRIDI_RELAXATION" ||
480  type_ == "LINESMOOTHING_TRIDI RELAXATION" ||
481  type_ == "LINESMOOTHING_TRIDIRELAXATION" ||
482  type_ == "LINESMOOTHING_TRIDIAGONAL_RELAXATION" ||
483  type_ == "LINESMOOTHING_TRIDIAGONAL RELAXATION" ||
484  type_ == "LINESMOOTHING_TRIDIAGONALRELAXATION")
485  type_ = "TRIDIAGONALRELAXATION";
486  else
487  type_ = "BLOCKRELAXATION";
488  } else {
489  // line detection failed -> fallback to point-wise relaxation
490  this->GetOStream(Runtime0) << "Line detection failed: fall back to point-wise relaxation" << std::endl;
491  myparamList.remove("partitioner: type",false);
492  myparamList.remove("partitioner: map", false);
493  myparamList.remove("partitioner: local parts",false);
494  type_ = "RELAXATION";
495  }
496 
497  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
498 
499  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
500  SetPrecParameters();
501  prec_->initialize();
502  prec_->compute();
503  }
504 
505  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
507  typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
508 
509  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
510  if (!bA.is_null())
511  A_ = bA->Merge();
512 
513  RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
514 
515  bool reusePreconditioner = false;
516  if (this->IsSetup() == true) {
517  // Reuse the constructed preconditioner
518  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): Setup() has already been called, assuming reuse" << std::endl;
519 
520  RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
521  if (!prec.is_null()) {
522 #ifdef IFPACK2_HAS_PROPER_REUSE
523  prec->resetMatrix(tA);
524  reusePreconditioner = true;
525 #else
526  this->GetOStream(Errors) << "Ifpack2 does not have proper reuse yet." << std::endl;
527 #endif
528 
529  } else {
530  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupBlockRelaxation(): reuse of this type is not available (failed cast to CanChangeMatrix), "
531  "reverting to full construction" << std::endl;
532  }
533  }
534 
535  if (!reusePreconditioner) {
536  ParameterList& myparamList = const_cast<ParameterList&>(this->GetParameterList());
537  myparamList.print();
538  if(myparamList.isParameter("partitioner: type") &&
539  myparamList.get<std::string>("partitioner: type") == "line") {
540  Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > xCoordinates =
541  Factory::Get<Teuchos::RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > >(currentLevel, "Coordinates");
542  Teuchos::RCP<Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > coordinates = Teuchos::rcpFromRef(Xpetra::toTpetra<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO>(*xCoordinates));
543 
544  size_t numDofsPerNode = A_->getNodeNumRows() / xCoordinates->getMap()->getNodeNumElements();
545  myparamList.set("partitioner: coordinates", coordinates);
546  myparamList.set("partitioner: PDE equations", (int) numDofsPerNode);
547  }
548 
549  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
550  SetPrecParameters();
551  prec_->initialize();
552  }
553 
554  prec_->compute();
555  }
556 
557  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
559  if (this->IsSetup() == true) {
560  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): SetupChebyshev() has already been called" << std::endl;
561  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupChebyshev(): reuse of this type is not available, reverting to full construction" << std::endl;
562  }
563 
564  typedef Teuchos::ScalarTraits<SC> STS;
565  SC negone = -STS::one();
566 
567  SC lambdaMax = negone;
568  {
569  std::string maxEigString = "chebyshev: max eigenvalue";
570  std::string eigRatioString = "chebyshev: ratio eigenvalue";
571 
572  ParameterList& paramList = const_cast<ParameterList&>(this->GetParameterList());
573 
574  // Get/calculate the maximum eigenvalue
575  if (paramList.isParameter(maxEigString)) {
576  if (paramList.isType<double>(maxEigString))
577  lambdaMax = paramList.get<double>(maxEigString);
578  else
579  lambdaMax = paramList.get<SC>(maxEigString);
580  this->GetOStream(Statistics1) << maxEigString << " (cached with smoother parameter list) = " << lambdaMax << std::endl;
581 
582  } else {
583  lambdaMax = A_->GetMaxEigenvalueEstimate();
584  if (lambdaMax != negone) {
585  this->GetOStream(Statistics1) << maxEigString << " (cached with matrix) = " << lambdaMax << std::endl;
586  paramList.set(maxEigString, lambdaMax);
587  }
588  }
589 
590  // Calculate the eigenvalue ratio
591  const SC defaultEigRatio = 20;
592 
593  SC ratio = defaultEigRatio;
594  if (paramList.isParameter(eigRatioString)) {
595  if (paramList.isType<double>(eigRatioString))
596  ratio = paramList.get<double>(eigRatioString);
597  else
598  ratio = paramList.get<SC>(eigRatioString);
599  }
600  if (currentLevel.GetLevelID()) {
601  // Update ratio to be
602  // ratio = max(number of fine DOFs / number of coarse DOFs, defaultValue)
603  //
604  // NOTE: We don't need to request previous level matrix as we know for sure it was constructed
605  RCP<const Matrix> fineA = currentLevel.GetPreviousLevel()->Get<RCP<Matrix> >("A");
606  size_t nRowsFine = fineA->getGlobalNumRows();
607  size_t nRowsCoarse = A_->getGlobalNumRows();
608 
609  SC levelRatio = as<SC>(as<float>(nRowsFine)/nRowsCoarse);
610  if (STS::magnitude(levelRatio) > STS::magnitude(ratio))
611  ratio = levelRatio;
612  }
613 
614  this->GetOStream(Statistics1) << eigRatioString << " (computed) = " << ratio << std::endl;
615  paramList.set(eigRatioString, ratio);
616  }
617 
618  RCP<const Tpetra::RowMatrix<SC, LO, GO, NO> > tA = Utilities::Op2NonConstTpetraRow(A_);
619 
620  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
621  SetPrecParameters();
622  {
623  SubFactoryMonitor(*this, "Preconditioner init", currentLevel);
624  prec_->initialize();
625  }
626  {
627  SubFactoryMonitor(*this, "Preconditioner compute", currentLevel);
628  prec_->compute();
629  }
630 
631  if (lambdaMax == negone) {
632  typedef Tpetra::RowMatrix<SC, LO, GO, NO> MatrixType;
633 
634  Teuchos::RCP<Ifpack2::Chebyshev<MatrixType> > chebyPrec = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
635  if (chebyPrec != Teuchos::null) {
636  lambdaMax = chebyPrec->getLambdaMaxForApply();
637  A_->SetMaxEigenvalueEstimate(lambdaMax);
638  this->GetOStream(Statistics1) << "chebyshev: max eigenvalue (calculated by Ifpack2)" << " = " << lambdaMax << std::endl;
639  }
640  TEUCHOS_TEST_FOR_EXCEPTION(lambdaMax == negone, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Setup(): no maximum eigenvalue estimate");
641  }
642  }
643 
644  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
646  typedef Tpetra::RowMatrix<SC,LO,GO,NO> tRowMatrix;
647 
648  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
649  if (!bA.is_null())
650  A_ = bA->Merge();
651 
652  RCP<const tRowMatrix> tA = Utilities::Op2NonConstTpetraRow(A_);
653 
654  bool reusePreconditioner = false;
655  if (this->IsSetup() == true) {
656  // Reuse the constructed preconditioner
657  this->GetOStream(Runtime1) << "MueLu::Ifpack2Smoother::SetupGeneric(): Setup() has already been called, assuming reuse" << std::endl;
658 
659  RCP<Ifpack2::Details::CanChangeMatrix<tRowMatrix> > prec = rcp_dynamic_cast<Ifpack2::Details::CanChangeMatrix<tRowMatrix> >(prec_);
660  if (!prec.is_null()) {
661 #ifdef IFPACK2_HAS_PROPER_REUSE
662  prec->resetMatrix(tA);
663  reusePreconditioner = true;
664 #else
665  this->GetOStream(Errors) << "Ifpack2 does not have proper reuse yet." << std::endl;
666 #endif
667 
668  } else {
669  this->GetOStream(Warnings0) << "MueLu::Ifpack2Smoother::SetupGeneric(): reuse of this type is not available (failed cast to CanChangeMatrix), "
670  "reverting to full construction" << std::endl;
671  }
672  }
673 
674  if (!reusePreconditioner) {
675  prec_ = Ifpack2::Factory::create(type_, tA, overlap_);
676  SetPrecParameters();
677  prec_->initialize();
678  }
679 
680  prec_->compute();
681  }
682 
683  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
684  void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
685  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::Ifpack2Smoother::Apply(): Setup() has not been called");
686 
687  // Forward the InitialGuessIsZero option to Ifpack2
688  // TODO: It might be nice to switch back the internal
689  // "zero starting solution" option of the ifpack2 object prec_ to his
690  // initial value at the end but there is no way right now to get
691  // the current value of the "zero starting solution" in ifpack2.
692  // It's not really an issue, as prec_ can only be used by this method.
693  // TODO: When https://software.sandia.gov/bugzilla/show_bug.cgi?id=5283#c2 is done
694  // we should remove the if/else/elseif and just test if this
695  // option is supported by current ifpack2 preconditioner
696  Teuchos::ParameterList paramList;
697  bool supportInitialGuess = false;
698  if (type_ == "CHEBYSHEV") {
699  paramList.set("chebyshev: zero starting solution", InitialGuessIsZero);
700  SetPrecParameters(paramList);
701  supportInitialGuess = true;
702 
703  } else if (type_ == "RELAXATION") {
704  paramList.set("relaxation: zero starting solution", InitialGuessIsZero);
705  SetPrecParameters(paramList);
706  supportInitialGuess = true;
707 
708  } else if (type_ == "KRYLOV") {
709  paramList.set("krylov: zero starting solution", InitialGuessIsZero);
710  SetPrecParameters(paramList);
711  supportInitialGuess = true;
712 
713  } else if (type_ == "SCHWARZ") {
714  paramList.set("schwarz: zero starting solution", InitialGuessIsZero);
715  //Because additive Schwarz has "delta" semantics, it's sufficient to
716  //toggle only the zero initial guess flag, and not pass in already
717  //set parameters. If we call SetPrecParameters, the subdomain solver
718  //will be destroyed.
719  prec_->setParameters(paramList);
720  supportInitialGuess = true;
721  }
722 
723  //TODO JJH 30Apr2014 Calling SetPrecParameters(paramList) when the smoother
724  //is Ifpack2::AdditiveSchwarz::setParameterList() will destroy the subdomain
725  //(aka inner) solver. This behavior is documented but a departure from what
726  //it previously did, and what other Ifpack2 solvers currently do. So I have
727  //moved SetPrecParameters(paramList) into the if-else block above.
728 
729  // Apply
730  if (InitialGuessIsZero || supportInitialGuess) {
731  Tpetra::MultiVector<SC,LO,GO,NO>& tpX = Utilities::MV2NonConstTpetraMV(X);
732  const Tpetra::MultiVector<SC,LO,GO,NO>& tpB = Utilities::MV2TpetraMV(B);
733  prec_->apply(tpB, tpX);
734  } else {
735  typedef Teuchos::ScalarTraits<Scalar> TST;
736  RCP<MultiVector> Residual = Utilities::Residual(*A_, X, B);
737  RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
738 
739  Tpetra::MultiVector<SC,LO,GO,NO>& tpX = Utilities::MV2NonConstTpetraMV(*Correction);
740  const Tpetra::MultiVector<SC,LO,GO,NO>& tpB = Utilities::MV2TpetraMV(*Residual);
741 
742  prec_->apply(tpB, tpX);
743 
744  X.update(TST::one(), *Correction, TST::one());
745  }
746  }
747 
748  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
749  RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
750  RCP<Ifpack2Smoother> smoother = rcp(new Ifpack2Smoother(*this) );
751  smoother->SetParameterList(this->GetParameterList());
752  return smoother;
753  }
754 
755  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
757  std::ostringstream out;
759  out << prec_->description();
760  } else {
762  out << "{type = " << type_ << "}";
763  }
764  return out.str();
765  }
766 
767  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
768  void Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream &out, const VerbLevel verbLevel) const {
770 
771  if (verbLevel & Parameters0)
772  out0 << "Prec. type: " << type_ << std::endl;
773 
774  if (verbLevel & Parameters1) {
775  out0 << "Parameter list: " << std::endl;
776  Teuchos::OSTab tab2(out);
777  out << this->GetParameterList();
778  out0 << "Overlap: " << overlap_ << std::endl;
779  }
780 
781  if (verbLevel & External)
782  if (prec_ != Teuchos::null) {
783  Teuchos::OSTab tab2(out);
784  out << *prec_ << std::endl;
785  }
786 
787  if (verbLevel & Debug) {
788  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
789  << "-" << std::endl
790  << "RCP<prec_>: " << prec_ << std::endl;
791  }
792  }
793 
794  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
796  typedef Tpetra::RowMatrix<SC,LO,GO,NO> MatrixType;
797  // NOTE: Only works for a subset of Ifpack2's smoothers
798  RCP<Ifpack2::Relaxation<MatrixType> > pr = rcp_dynamic_cast<Ifpack2::Relaxation<MatrixType> >(prec_);
799  if(!pr.is_null()) return pr->getNodeSmootherComplexity();
800 
801  RCP<Ifpack2::Chebyshev<MatrixType> > pc = rcp_dynamic_cast<Ifpack2::Chebyshev<MatrixType> >(prec_);
802  if(!pc.is_null()) return pc->getNodeSmootherComplexity();
803 
804  RCP<Ifpack2::BlockRelaxation<MatrixType> > pb = rcp_dynamic_cast<Ifpack2::BlockRelaxation<MatrixType> >(prec_);
805  if(!pb.is_null()) return pb->getNodeSmootherComplexity();
806 
807  RCP<Ifpack2::ILUT<MatrixType> > pi = rcp_dynamic_cast<Ifpack2::ILUT<MatrixType> >(prec_);
808  if(!pi.is_null()) return pi->getNodeSmootherComplexity();
809 
810  RCP<Ifpack2::RILUK<MatrixType> > pk = rcp_dynamic_cast<Ifpack2::RILUK<MatrixType> >(prec_);
811  if(!pk.is_null()) return pk->getNodeSmootherComplexity();
812 
813 
814  return Teuchos::OrdinalTraits<size_t>::invalid();
815  }
816 
817 
818 } // namespace MueLu
819 
820 #endif // HAVE_MUELU_TPETRA && HAVE_MUELU_IFPACK2
821 #endif // MUELU_IFPACK2SMOOTHER_DEF_HPP
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
MueLu::DefaultLocalOrdinal LocalOrdinal
static Teuchos::RCP< Preconditioner< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > create(const std::string &precType, const Teuchos::RCP< const MatrixType > &matrix)
virtual std::string description() const
Return a simple one-line description of this object.
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that encapsulates Ifpack2 smoothers.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
void SetupSchwarz(Level &currentLevel)
void Setup(Level &currentLevel)
Set up the smoother.
void SetupBlockRelaxation(Level &currentLevel)
void SetupChebyshev(Level &currentLevel)
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
RCP< SmootherPrototype > Copy() const
std::string description() const
Return a simple one-line description of this object.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
void SetupTopological(Level &currentLevel)
friend class Ifpack2Smoother
Constructor.
void DeclareInput(Level &currentLevel) const
Input.
void SetupLineSmoothing(Level &currentLevel)
void SetupGeneric(Level &currentLevel)
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
RCP< Level > & GetPreviousLevel()
Previous level.
virtual void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
bool IsSetup() const
Get the state of a smoother prototype.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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)
static Teuchos::ArrayRCP< const bool > DetectDirichletRows(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::magnitude(0.), const bool count_twos_as_dirichlet=false)
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.
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Debug
Print additional debugging information.
@ Errors
Errors.
@ Statistics1
Print more statistics.
@ External
Print external lib objects.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
@ Parameters0
Print class parameters.
@ Parameters1
Print class parameters (more parameters, more verbose)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.