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