MueLu  Version of the Day
MueLu_SimpleSmoother_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 
47 #ifndef MUELU_SIMPLESMOOTHER_DEF_HPP_
48 #define MUELU_SIMPLESMOOTHER_DEF_HPP_
49 
50 #include "Teuchos_ArrayViewDecl.hpp"
51 #include "Teuchos_ScalarTraits.hpp"
52 
53 #include "MueLu_ConfigDefs.hpp"
54 
55 #include <Xpetra_Matrix.hpp>
56 #include <Xpetra_CrsMatrixWrap.hpp>
57 #include <Xpetra_BlockedCrsMatrix.hpp>
58 #include <Xpetra_MultiVectorFactory.hpp>
59 #include <Xpetra_VectorFactory.hpp>
60 #include <Xpetra_ReorderedBlockedCrsMatrix.hpp>
61 
63 #include "MueLu_Level.hpp"
64 #include "MueLu_Utilities.hpp"
65 #include "MueLu_Monitor.hpp"
66 #include "MueLu_HierarchyUtils.hpp"
67 #include "MueLu_SmootherBase.hpp"
68 #include "MueLu_SubBlockAFactory.hpp"
69 
70 // include files for default FactoryManager
71 #include "MueLu_SchurComplementFactory.hpp"
72 #include "MueLu_DirectSolver.hpp"
73 #include "MueLu_SmootherFactory.hpp"
74 #include "MueLu_FactoryManager.hpp"
75 
76 namespace MueLu {
77 
78  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
80  : type_("SIMPLE"), A_(Teuchos::null) {}
81 
82  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
84 
85  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87  RCP<ParameterList> validParamList = rcp(new ParameterList());
88 
89  validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory of the matrix A");
90  validParamList->set<Scalar>("Damping factor", 1.0, "Damping/Scaling factor in SIMPLE");
91  validParamList->set<LocalOrdinal>("Sweeps", 1, "Number of SIMPLE sweeps (default = 1)");
92  validParamList->set<bool>("UseSIMPLEC", false, "Use SIMPLEC instead of SIMPLE (default = false)");
93 
94  return validParamList;
95  }
96 
98  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
99  void SimpleSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::AddFactoryManager(RCP<const FactoryManagerBase> FactManager, int pos) {
100  TEUCHOS_TEST_FOR_EXCEPTION(pos < 0, Exceptions::RuntimeError, "MueLu::SimpleSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
101 
102  size_t myPos = Teuchos::as<size_t>(pos);
103 
104  if (myPos < FactManager_.size()) {
105  // replace existing entris in FactManager_ vector
106  FactManager_.at(myPos) = FactManager;
107  } else if( myPos == FactManager_.size()) {
108  // add new Factory manager in the end of the vector
109  FactManager_.push_back(FactManager);
110  } else { // if(myPos > FactManager_.size())
111  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
112  *out << "Warning: cannot add new FactoryManager at proper position " << pos << ". The FactoryManager is just appended to the end. Check this!" << std::endl;
113 
114  // add new Factory manager in the end of the vector
115  FactManager_.push_back(FactManager);
116  }
117 
118  }
119 
120 
121  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
123  AddFactoryManager(FactManager, 0); // overwrite factory manager for predicting the primary variable
124  }
125 
126  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
128  AddFactoryManager(FactManager, 1); // overwrite factory manager for SchurComplement
129  }
130 
131  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
133  currentLevel.DeclareInput("A",this->GetFactory("A").get());
134 
135  TEUCHOS_TEST_FOR_EXCEPTION(FactManager_.size() != 2, Exceptions::RuntimeError,"MueLu::SimpleSmoother::DeclareInput: You have to declare two FactoryManagers with a \"Smoother\" object: One for predicting the primary variable and one for the SchurComplement system. The smoother for the SchurComplement system needs a SchurComplementFactory as input for variable \"A\". make sure that you use the same proper damping factors for omega both in the SchurComplementFactory and in the SIMPLE smoother!");
136 
137  // loop over all factory managers for the subblocks of blocked operator A
138  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
139  for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
140  SetFactoryManager currentSFM (rcpFromRef(currentLevel), *it);
141 
142  // request "Smoother" for current subblock row.
143  currentLevel.DeclareInput("PreSmoother",(*it)->GetFactory("Smoother").get());
144  }
145  }
146 
147  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
149 
150  FactoryMonitor m(*this, "Setup blocked SIMPLE Smoother", currentLevel);
151 
152  if (SmootherPrototype::IsSetup() == true)
153  this->GetOStream(Warnings0) << "MueLu::SimpleSmoother::Setup(): Setup() has already been called";
154 
155  // extract blocked operator A from current level
156  A_ = Factory::Get<RCP<Matrix> > (currentLevel, "A");
157 
158  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
159  TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::SimpleSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
160 
161  // store map extractors
162  rangeMapExtractor_ = bA->getRangeMapExtractor();
163  domainMapExtractor_ = bA->getDomainMapExtractor();
164 
165  // Store the blocks in local member variables
166  F_ = bA->getMatrix(0, 0);
167  G_ = bA->getMatrix(0, 1);
168  D_ = bA->getMatrix(1, 0);
169  Z_ = bA->getMatrix(1, 1);
170 
171  const ParameterList & pL = Factory::GetParameterList();
172  bool bSIMPLEC = pL.get<bool>("UseSIMPLEC");
173 
174  // Create the inverse of the diagonal of F
175  // TODO add safety check for zeros on diagonal of F!
176  RCP<Vector> diagFVector = VectorFactory::Build(F_->getRowMap());
177  if(!bSIMPLEC) {
178  F_->getLocalDiagCopy(*diagFVector); // extract diagonal of F
179  } else {
180  /*const RCP<const Map> rowmap = F_->getRowMap();
181  size_t locSize = rowmap->getNodeNumElements();
182  Teuchos::ArrayRCP<SC> diag = diagFVector->getDataNonConst(0);
183  Teuchos::ArrayView<const LO> cols;
184  Teuchos::ArrayView<const SC> vals;
185  for (size_t i=0; i<locSize; ++i) { // loop over rows
186  F_->getLocalRowView(i,cols,vals);
187  Scalar absRowSum = Teuchos::ScalarTraits<Scalar>::zero();
188  for (LO j=0; j<cols.size(); ++j) { // loop over cols
189  absRowSum += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
190  }
191  diag[i] = absRowSum;
192  }*/
193  // TODO this does not work if F_ is nested!
194  diagFVector = Utilities::GetLumpedMatrixDiagonal(*F_);
195  }
196  diagFinv_ = Utilities::GetInverse(diagFVector);
197 
198  // check whether diagFinv_ is a blocked vector with only 1 block
199  RCP<BlockedVector> bdiagFinv = Teuchos::rcp_dynamic_cast<BlockedVector>(diagFinv_);
200  if(bdiagFinv.is_null() == false && bdiagFinv->getBlockedMap()->getNumMaps() == 1) {
201  RCP<Vector> nestedVec = bdiagFinv->getMultiVector(0,bdiagFinv->getBlockedMap()->getThyraMode())->getVectorNonConst(0);
202  diagFinv_.swap(nestedVec);
203  }
204 
205  // Set the Smoother
206  // carefully switch to the SubFactoryManagers (defined by the users)
207  {
208  RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
209  SetFactoryManager currentSFM (rcpFromRef(currentLevel), velpredictFactManager);
210  velPredictSmoo_ = currentLevel.Get< RCP<SmootherBase> >("PreSmoother",velpredictFactManager->GetFactory("Smoother").get());
211  }
212  {
213  RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
214  SetFactoryManager currentSFM (rcpFromRef(currentLevel), schurFactManager);
215  schurCompSmoo_ = currentLevel.Get< RCP<SmootherBase> >("PreSmoother", schurFactManager->GetFactory("Smoother").get());
216  }
217 
219  }
220 
221  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
222  void SimpleSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const
223  {
224  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError,
225  "MueLu::SimpleSmoother::Apply(): Setup() has not been called");
226 #if 0
227  // TODO simplify this debug check
228  RCP<MultiVector> rcpDebugX = Teuchos::rcpFromRef(X);
229  RCP<const MultiVector> rcpDebugB = Teuchos::rcpFromRef(B);
230  RCP<BlockedMultiVector> rcpBDebugX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpDebugX);
231  RCP<const BlockedMultiVector> rcpBDebugB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpDebugB);
232  if(rcpBDebugB.is_null() == false) {
233  //TEUCHOS_TEST_FOR_EXCEPTION(A_->getRangeMap()->isSameAs(*(B.getMap())) == false, Exceptions::RuntimeError, "MueLu::SimpleSmoother::Apply(): The map of RHS vector B is not the same as range map of the blocked operator A. Please check the map of B and A.");
234  } else {
235  //TEUCHOS_TEST_FOR_EXCEPTION(bA->getFullRangeMap()->isSameAs(*(B.getMap())) == false, Exceptions::RuntimeError, "MueLu::SimpleSmoother::Apply(): The map of RHS vector B is not the same as range map of the blocked operator A. Please check the map of B and A.");
236  }
237  if(rcpBDebugX.is_null() == false) {
238  //TEUCHOS_TEST_FOR_EXCEPTION(A_->getDomainMap()->isSameAs(*(X.getMap())) == false, Exceptions::RuntimeError, "MueLu::SimpleSmoother::Apply(): The map of the solution vector X is not the same as domain map of the blocked operator A. Please check the map of X and A.");
239  } else {
240  //TEUCHOS_TEST_FOR_EXCEPTION(bA->getFullDomainMap()->isSameAs(*(X.getMap())) == false, Exceptions::RuntimeError, "MueLu::SimpleSmoother::Apply(): The map of the solution vector X is not the same as domain map of the blocked operator A. Please check the map of X and A.");
241  }
242 #endif
243 
244  const SC zero = Teuchos::ScalarTraits<SC>::zero();
245  const SC one = Teuchos::ScalarTraits<SC>::one();
246 
247  // extract parameters from internal parameter list
248  const ParameterList & pL = Factory::GetParameterList();
249  LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
250  Scalar omega = pL.get<Scalar>("Damping factor");
251 
252  // The boolean flags check whether we use Thyra or Xpetra style GIDs
253  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode(); // && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_) == Teuchos::null);
254  bool bDomainThyraMode = domainMapExtractor_->getThyraMode(); // && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_) == Teuchos::null);
255 
256  //RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
257 
258  // wrap current solution vector in RCP
259  RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
260  RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
261 
262  // make sure that both rcpX and rcpB are BlockedMultiVector objects
263  bool bCopyResultX = false;
264  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
266  "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
267  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
268  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
269 
270  if(bX.is_null() == true) {
271  RCP<MultiVector> test = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(),rcpX));
272  rcpX.swap(test);
273  bCopyResultX = true;
274  }
275 
276  if(bB.is_null() == true) {
277  RCP<const MultiVector> test = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(),rcpB));
278  rcpB.swap(test);
279  }
280 
281  // we now can guarantee that X and B are blocked multi vectors
282  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
283  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
284 
285  // check the type of operator
286  RCP<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> > rbA =
287  Teuchos::rcp_dynamic_cast<Xpetra::ReorderedBlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(bA);
288  if(rbA.is_null() == false) {
289  // A is a ReorderedBlockedCrsMatrix
290  Teuchos::RCP<const Xpetra::BlockReorderManager > brm = rbA->getBlockReorderManager();
291 
292  // check vector types
293  if(bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps()) {
294  // X is a blocked multi vector but incompatible to the reordered blocked operator A
295  Teuchos::RCP<MultiVector> test = buildReorderedBlockedMultiVector(brm, bX);
296  rcpX.swap(test);
297  }
298  if(bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps()) {
299  // B is a blocked multi vector but incompatible to the reordered blocked operator A
300  Teuchos::RCP<const MultiVector> test = buildReorderedBlockedMultiVector(brm, bB);
301  rcpB.swap(test);
302  }
303  }
304 
305  // Throughout the rest of the algorithm rcpX and rcpB are used for solution vector and RHS
306 
307  // create residual vector
308  // contains current residual of current solution X with rhs B
309  RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
310  RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
311  RCP<MultiVector> r1 = bresidual->getMultiVector(0,bRangeThyraMode);
312  RCP<MultiVector> r2 = bresidual->getMultiVector(1,bRangeThyraMode);
313 
314  // helper vector 1
315  RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
316  RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
317  RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0, bDomainThyraMode);
318  RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1, bDomainThyraMode);
319 
320  // helper vector 2
321  RCP<MultiVector> xhat = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
322  RCP<BlockedMultiVector> bxhat = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xhat);
323  RCP<MultiVector> xhat1 = bxhat->getMultiVector(0, bDomainThyraMode);
324  RCP<MultiVector> xhat2 = bxhat->getMultiVector(1, bDomainThyraMode);
325 
326 
327  // incrementally improve solution vector X
328  for (LocalOrdinal run = 0; run < nSweeps; ++run) {
329  // 1) calculate current residual
330  residual->update(one, *rcpB, zero); // residual = B
331  if(InitialGuessIsZero == false || run > 0)
332  A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
333 
334  // 2) solve F * \Delta \tilde{x}_1 = r_1
335  // start with zero guess \Delta \tilde{x}_1
336  xtilde1->putScalar(zero);
337  xtilde2->putScalar(zero);
338  velPredictSmoo_->Apply(*xtilde1, *r1);
339 
340  // 3) calculate rhs for SchurComp equation
341  // r_2 - D \Delta \tilde{x}_1
342  RCP<MultiVector> schurCompRHS = rangeMapExtractor_->getVector(1, rcpB->getNumVectors(), bRangeThyraMode);
343  D_->apply(*xtilde1, *schurCompRHS);
344 
345  schurCompRHS->update(one, *r2, -one);
346 
347  // 4) solve SchurComp equation
348  // start with zero guess \Delta \tilde{x}_2
349  schurCompSmoo_->Apply(*xtilde2, *schurCompRHS);
350 
351  // 5) scale xtilde2 with omega
352  // store this in xhat2
353  xhat2->update(omega, *xtilde2, zero);
354 
355  // 6) calculate xhat1
356  RCP<MultiVector> xhat1_temp = domainMapExtractor_->getVector(0, rcpX->getNumVectors(), bDomainThyraMode);
357  G_->apply(*xhat2, *xhat1_temp); // store result temporarely in xtilde1_temp
358 
359  xhat1->elementWiseMultiply(one/*/omega*/, *diagFinv_, *xhat1_temp, zero);
360  xhat1->update(one, *xtilde1, -one);
361 
362  rcpX->update(one, *bxhat, one);
363  }
364 
365  if (bCopyResultX == true) {
366  RCP<const MultiVector> Xmerged = bX->Merge();
367  X.update(one, *Xmerged, zero);
368  }
369 
370  }
371 
372  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
373  RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
375  return rcp( new SimpleSmoother(*this) );
376  }
377 
378  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
380  std::ostringstream out;
382  out << "{type = " << type_ << "}";
383  return out.str();
384  }
385 
386  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
387  void SimpleSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream &out, const VerbLevel verbLevel) const {
389 
390  if (verbLevel & Parameters0) {
391  out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
392  }
393 
394  if (verbLevel & Debug) {
395  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
396  }
397  }
398 
399  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
401  // FIXME: This is a placeholder
402  return Teuchos::OrdinalTraits<size_t>::invalid();
403  }
404 
405 
406 } // namespace MueLu
407 
408 
409 #endif /* MUELU_SIMPLESMOOTHER_DEF_HPP_ */
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
#define MUELU_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
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 holds all level-specific information.
Definition: MueLu_Level.hpp:99
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual const Teuchos::ParameterList & GetParameterList() const
An exception safe way to call the method 'Level::SetFactoryManager()'.
SIMPLE smoother for 2x2 block matrices.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
void Setup(Level &currentLevel)
Setup routine.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
std::string description() const
Return a simple one-line description of this object.
void DeclareInput(Level &currentLevel) const
Input.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
RCP< const ParameterList > GetValidParameterList() const
Input.
virtual ~SimpleSmoother()
Destructor.
RCP< SmootherPrototype > Copy() const
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos)
Add a factory manager at a specific position.
void Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
bool IsSetup() const
Get the state of a smoother prototype.
static Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetInverse(Teuchos::RCP< const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero())
static Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetLumpedMatrixDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false)
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Debug
Print additional debugging information.
@ Parameters0
Print class parameters.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.