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 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_SIMPLESMOOTHER_DEF_HPP_
11 #define MUELU_SIMPLESMOOTHER_DEF_HPP_
12 
13 #include "Teuchos_ArrayViewDecl.hpp"
14 #include "Teuchos_ScalarTraits.hpp"
15 
16 #include "MueLu_ConfigDefs.hpp"
17 
18 #include <Xpetra_Matrix.hpp>
19 #include <Xpetra_CrsMatrixWrap.hpp>
21 #include <Xpetra_MultiVectorFactory.hpp>
22 #include <Xpetra_VectorFactory.hpp>
24 
26 #include "MueLu_Level.hpp"
27 #include "MueLu_Utilities.hpp"
28 #include "MueLu_Monitor.hpp"
29 #include "MueLu_HierarchyUtils.hpp"
30 #include "MueLu_SmootherBase.hpp"
31 
32 // include files for default FactoryManager
33 #include "MueLu_SchurComplementFactory.hpp"
34 #include "MueLu_FactoryManager.hpp"
35 
36 namespace MueLu {
37 
38 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40  : type_("SIMPLE")
41  , A_(Teuchos::null) {}
42 
43 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45 
46 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
48  RCP<ParameterList> validParamList = rcp(new ParameterList());
49 
50  validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory of the matrix A");
51  validParamList->set<Scalar>("Damping factor", 1.0, "Damping/Scaling factor in SIMPLE");
52  validParamList->set<LocalOrdinal>("Sweeps", 1, "Number of SIMPLE sweeps (default = 1)");
53  validParamList->set<bool>("UseSIMPLEC", false, "Use SIMPLEC instead of SIMPLE (default = false)");
54 
55  return validParamList;
56 }
57 
59 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61  TEUCHOS_TEST_FOR_EXCEPTION(pos < 0, Exceptions::RuntimeError, "MueLu::SimpleSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
62 
63  size_t myPos = Teuchos::as<size_t>(pos);
64 
65  if (myPos < FactManager_.size()) {
66  // replace existing entris in FactManager_ vector
67  FactManager_.at(myPos) = FactManager;
68  } else if (myPos == FactManager_.size()) {
69  // add new Factory manager in the end of the vector
70  FactManager_.push_back(FactManager);
71  } else { // if(myPos > FactManager_.size())
72  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
73  *out << "Warning: cannot add new FactoryManager at proper position " << pos << ". The FactoryManager is just appended to the end. Check this!" << std::endl;
74 
75  // add new Factory manager in the end of the vector
76  FactManager_.push_back(FactManager);
77  }
78 }
79 
80 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82  AddFactoryManager(FactManager, 0); // overwrite factory manager for predicting the primary variable
83 }
84 
85 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87  AddFactoryManager(FactManager, 1); // overwrite factory manager for SchurComplement
88 }
89 
90 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
92  currentLevel.DeclareInput("A", this->GetFactory("A").get());
93 
94  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!");
95 
96  // loop over all factory managers for the subblocks of blocked operator A
97  std::vector<Teuchos::RCP<const FactoryManagerBase>>::const_iterator it;
98  for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
99  SetFactoryManager currentSFM(rcpFromRef(currentLevel), *it);
100 
101  // request "Smoother" for current subblock row.
102  currentLevel.DeclareInput("PreSmoother", (*it)->GetFactory("Smoother").get());
103  }
104 }
105 
106 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
108  FactoryMonitor m(*this, "Setup blocked SIMPLE Smoother", currentLevel);
109 
110  if (SmootherPrototype::IsSetup() == true)
111  this->GetOStream(Warnings0) << "MueLu::SimpleSmoother::Setup(): Setup() has already been called";
112 
113  // extract blocked operator A from current level
114  A_ = Factory::Get<RCP<Matrix>>(currentLevel, "A");
115 
116  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
117  TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::SimpleSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
118 
119  // store map extractors
120  rangeMapExtractor_ = bA->getRangeMapExtractor();
121  domainMapExtractor_ = bA->getDomainMapExtractor();
122 
123  // Store the blocks in local member variables
124  F_ = bA->getMatrix(0, 0);
125  G_ = bA->getMatrix(0, 1);
126  D_ = bA->getMatrix(1, 0);
127  Z_ = bA->getMatrix(1, 1);
128 
130  bool bSIMPLEC = pL.get<bool>("UseSIMPLEC");
131 
132  // Create the inverse of the diagonal of F
133  // TODO add safety check for zeros on diagonal of F!
134  RCP<Vector> diagFVector = VectorFactory::Build(F_->getRowMap());
135  if (!bSIMPLEC) {
136  F_->getLocalDiagCopy(*diagFVector); // extract diagonal of F
137  } else {
138  /*const RCP<const Map> rowmap = F_->getRowMap();
139  size_t locSize = rowmap->getLocalNumElements();
140  Teuchos::ArrayRCP<SC> diag = diagFVector->getDataNonConst(0);
141  Teuchos::ArrayView<const LO> cols;
142  Teuchos::ArrayView<const SC> vals;
143  for (size_t i=0; i<locSize; ++i) { // loop over rows
144  F_->getLocalRowView(i,cols,vals);
145  Scalar absRowSum = Teuchos::ScalarTraits<Scalar>::zero();
146  for (LO j=0; j<cols.size(); ++j) { // loop over cols
147  absRowSum += Teuchos::ScalarTraits<Scalar>::magnitude(vals[j]);
148  }
149  diag[i] = absRowSum;
150  }*/
151  // TODO this does not work if F_ is nested!
152  diagFVector = Utilities::GetLumpedMatrixDiagonal(*F_);
153  }
154  diagFinv_ = Utilities::GetInverse(diagFVector);
155 
156  // check whether diagFinv_ is a blocked vector with only 1 block
157  RCP<BlockedVector> bdiagFinv = Teuchos::rcp_dynamic_cast<BlockedVector>(diagFinv_);
158  if (bdiagFinv.is_null() == false && bdiagFinv->getBlockedMap()->getNumMaps() == 1) {
159  RCP<Vector> nestedVec = bdiagFinv->getMultiVector(0, bdiagFinv->getBlockedMap()->getThyraMode())->getVectorNonConst(0);
160  diagFinv_.swap(nestedVec);
161  }
162 
163  // Set the Smoother
164  // carefully switch to the SubFactoryManagers (defined by the users)
165  {
166  RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
167  SetFactoryManager currentSFM(rcpFromRef(currentLevel), velpredictFactManager);
168  velPredictSmoo_ = currentLevel.Get<RCP<SmootherBase>>("PreSmoother", velpredictFactManager->GetFactory("Smoother").get());
169  }
170  {
171  RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
172  SetFactoryManager currentSFM(rcpFromRef(currentLevel), schurFactManager);
173  schurCompSmoo_ = currentLevel.Get<RCP<SmootherBase>>("PreSmoother", schurFactManager->GetFactory("Smoother").get());
174  }
175 
177 }
178 
179 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
180 void SimpleSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero) const {
182  "MueLu::SimpleSmoother::Apply(): Setup() has not been called");
183 #if 0
184  // TODO simplify this debug check
185  RCP<MultiVector> rcpDebugX = Teuchos::rcpFromRef(X);
186  RCP<const MultiVector> rcpDebugB = Teuchos::rcpFromRef(B);
187  RCP<BlockedMultiVector> rcpBDebugX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpDebugX);
188  RCP<const BlockedMultiVector> rcpBDebugB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpDebugB);
189  if(rcpBDebugB.is_null() == false) {
190  //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.");
191  } else {
192  //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.");
193  }
194  if(rcpBDebugX.is_null() == false) {
195  //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.");
196  } else {
197  //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.");
198  }
199 #endif
200 
201  const SC zero = Teuchos::ScalarTraits<SC>::zero();
202  const SC one = Teuchos::ScalarTraits<SC>::one();
203 
204  // extract parameters from internal parameter list
206  LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
207  Scalar omega = pL.get<Scalar>("Damping factor");
208 
209  // The boolean flags check whether we use Thyra or Xpetra style GIDs
210  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode(); // && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_) == Teuchos::null);
211  bool bDomainThyraMode = domainMapExtractor_->getThyraMode(); // && (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(F_) == Teuchos::null);
212 
213  // RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
214 
215  // wrap current solution vector in RCP
216  RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
217  RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
218 
219  // make sure that both rcpX and rcpB are BlockedMultiVector objects
220  bool bCopyResultX = false;
221  bool bReorderX = false;
222  bool bReorderB = false;
223  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
224  MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
225  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
226  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
227 
228  // check the type of operator
231 
232  if (rbA.is_null() == false) {
233  // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
234  // map for the construction of the blocked multivectors
235 
236  // check type of X vector
237  if (bX.is_null() == true) {
238  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
239  rcpX.swap(vectorWithBlockedMap);
240  bCopyResultX = true;
241  bReorderX = true;
242  }
243 
244  // check type of B vector
245  if (bB.is_null() == true) {
246  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
247  rcpB.swap(vectorWithBlockedMap);
248  bReorderB = true;
249  }
250  } else {
251  // A is a BlockedCrsMatrix and uses a plain blocked map
252  if (bX.is_null() == true) {
253  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
254  rcpX.swap(vectorWithBlockedMap);
255  bCopyResultX = true;
256  }
257 
258  if (bB.is_null() == true) {
259  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(), rcpB));
260  rcpB.swap(vectorWithBlockedMap);
261  }
262  }
263 
264  // we now can guarantee that X and B are blocked multi vectors
265  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
266  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
267 
268  // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
269  if (rbA.is_null() == false) {
270  Teuchos::RCP<const Xpetra::BlockReorderManager> brm = rbA->getBlockReorderManager();
271 
272  // check type of X vector
273  if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
274  // X is a blocked multi vector but incompatible to the reordered blocked operator A
275  Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
276  rcpX.swap(reorderedBlockedVector);
277  }
278 
279  if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
280  // B is a blocked multi vector but incompatible to the reordered blocked operator A
281  Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
282  rcpB.swap(reorderedBlockedVector);
283  }
284  }
285 
286  // Throughout the rest of the algorithm rcpX and rcpB are used for solution vector and RHS
287 
288  // create residual vector
289  // contains current residual of current solution X with rhs B
290  RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
291  RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
292  RCP<MultiVector> r1 = bresidual->getMultiVector(0, bRangeThyraMode);
293  RCP<MultiVector> r2 = bresidual->getMultiVector(1, bRangeThyraMode);
294 
295  // helper vector 1
296  RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
297  RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
298  RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0, bDomainThyraMode);
299  RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1, bDomainThyraMode);
300 
301  // helper vector 2
302  RCP<MultiVector> xhat = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
303  RCP<BlockedMultiVector> bxhat = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xhat);
304  RCP<MultiVector> xhat1 = bxhat->getMultiVector(0, bDomainThyraMode);
305  RCP<MultiVector> xhat2 = bxhat->getMultiVector(1, bDomainThyraMode);
306 
307  // incrementally improve solution vector X
308  for (LocalOrdinal run = 0; run < nSweeps; ++run) {
309  // 1) calculate current residual
310  residual->update(one, *rcpB, zero); // residual = B
311  if (InitialGuessIsZero == false || run > 0)
312  A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
313 
314  // 2) solve F * \Delta \tilde{x}_1 = r_1
315  // start with zero guess \Delta \tilde{x}_1
316  xtilde1->putScalar(zero);
317  xtilde2->putScalar(zero);
318  velPredictSmoo_->Apply(*xtilde1, *r1);
319 
320  // 3) calculate rhs for SchurComp equation
321  // r_2 - D \Delta \tilde{x}_1
322  RCP<MultiVector> schurCompRHS = rangeMapExtractor_->getVector(1, rcpB->getNumVectors(), bRangeThyraMode);
323  D_->apply(*xtilde1, *schurCompRHS);
324 
325  schurCompRHS->update(one, *r2, -one);
326 
327  // 4) solve SchurComp equation
328  // start with zero guess \Delta \tilde{x}_2
329  schurCompSmoo_->Apply(*xtilde2, *schurCompRHS);
330 
331  // 5) scale xtilde2 with omega
332  // store this in xhat2
333  xhat2->update(omega, *xtilde2, zero);
334 
335  // 6) calculate xhat1
336  RCP<MultiVector> xhat1_temp = domainMapExtractor_->getVector(0, rcpX->getNumVectors(), bDomainThyraMode);
337  G_->apply(*xhat2, *xhat1_temp); // store result temporarely in xtilde1_temp
338 
339  xhat1->elementWiseMultiply(one /*/omega*/, *diagFinv_, *xhat1_temp, zero);
340  xhat1->update(one, *xtilde1, -one);
341 
342  rcpX->update(one, *bxhat, one);
343  }
344 
345  if (bCopyResultX == true) {
346  RCP<const MultiVector> Xmerged = bX->Merge();
347  X.update(one, *Xmerged, zero);
348  }
349 }
350 
351 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354  return rcp(new SimpleSmoother(*this));
355 }
356 
357 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
359  std::ostringstream out;
361  out << "{type = " << type_ << "}";
362  return out.str();
363 }
364 
365 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
368 
369  if (verbLevel & Parameters0) {
370  out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
371  }
372 
373  if (verbLevel & Debug) {
374  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
375  }
376 }
377 
378 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
380  // FIXME: This is a placeholder
382 }
383 
384 } // namespace MueLu
385 
386 #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)
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.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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:63
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