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