MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_UzawaSmoother_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_UZAWASMOOTHER_DEF_HPP_
11 #define MUELU_UZAWASMOOTHER_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>
20 #include <Xpetra_MultiVectorFactory.hpp>
21 #include <Xpetra_VectorFactory.hpp>
23 
25 #include "MueLu_Level.hpp"
26 #include "MueLu_Monitor.hpp"
27 #include "MueLu_HierarchyUtils.hpp"
28 #include "MueLu_SmootherBase.hpp"
29 
30 #include "MueLu_FactoryManager.hpp"
31 
32 namespace MueLu {
33 
34 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36  RCP<ParameterList> validParamList = rcp(new ParameterList());
37 
38  validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory of the matrix A");
39  validParamList->set<Scalar>("Damping factor", 1.0, "Damping/Scaling factor in SIMPLE");
40  validParamList->set<LocalOrdinal>("Sweeps", 1, "Number of SIMPLE sweeps (default = 1)");
41 
42  return validParamList;
43 }
44 
45 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
47  TEUCHOS_TEST_FOR_EXCEPTION(pos < 0, Exceptions::RuntimeError, "MueLu::UzawaSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
48 
49  size_t myPos = Teuchos::as<size_t>(pos);
50 
51  if (myPos < FactManager_.size()) {
52  // replace existing entris in FactManager_ vector
53  FactManager_.at(myPos) = FactManager;
54  } else if (myPos == FactManager_.size()) {
55  // add new Factory manager in the end of the vector
56  FactManager_.push_back(FactManager);
57  } else { // if(myPos > FactManager_.size())
58  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
59  *out << "Warning: cannot add new FactoryManager at proper position " << pos << ". The FactoryManager is just appended to the end. Check this!" << std::endl;
60 
61  // add new Factory manager in the end of the vector
62  FactManager_.push_back(FactManager);
63  }
64 }
65 
66 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68  AddFactoryManager(FactManager, 0); // overwrite factory manager for predicting the primary variable
69 }
70 
71 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
73  AddFactoryManager(FactManager, 1); // overwrite factory manager for SchurComplement
74 }
75 
76 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78  currentLevel.DeclareInput("A", this->GetFactory("A").get());
79 
80  TEUCHOS_TEST_FOR_EXCEPTION(FactManager_.size() != 2, Exceptions::RuntimeError, "MueLu::UzawaSmoother::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!");
81 
82  // loop over all factory managers for the subblocks of blocked operator A
83  std::vector<Teuchos::RCP<const FactoryManagerBase>>::const_iterator it;
84  for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
85  SetFactoryManager currentSFM(rcpFromRef(currentLevel), *it);
86 
87  // request "Smoother" for current subblock row.
88  currentLevel.DeclareInput("PreSmoother", (*it)->GetFactory("Smoother").get());
89  }
90 }
91 
92 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
94  FactoryMonitor m(*this, "Setup blocked Uzawa Smoother", currentLevel);
95 
96  if (SmootherPrototype::IsSetup() == true)
97  this->GetOStream(Warnings0) << "MueLu::UzawaSmoother::Setup(): Setup() has already been called";
98 
99  // extract blocked operator A from current level
100  A_ = Factory::Get<RCP<Matrix>>(currentLevel, "A");
101 
102  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
103  TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::UzawaSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
104 
105  // store map extractors
106  rangeMapExtractor_ = bA->getRangeMapExtractor();
107  domainMapExtractor_ = bA->getDomainMapExtractor();
108 
109  // Store the blocks in local member variables
110  F_ = bA->getMatrix(0, 0);
111  G_ = bA->getMatrix(0, 1);
112  D_ = bA->getMatrix(1, 0);
113  Z_ = bA->getMatrix(1, 1);
114 
115  // Set the Smoother
116  // carefully switch to the SubFactoryManagers (defined by the users)
117  {
118  RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
119  SetFactoryManager currentSFM(rcpFromRef(currentLevel), velpredictFactManager);
120  velPredictSmoo_ = currentLevel.Get<RCP<SmootherBase>>("PreSmoother", velpredictFactManager->GetFactory("Smoother").get());
121  }
122  {
123  RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
124  SetFactoryManager currentSFM(rcpFromRef(currentLevel), schurFactManager);
125  schurCompSmoo_ = currentLevel.Get<RCP<SmootherBase>>("PreSmoother", schurFactManager->GetFactory("Smoother").get());
126  }
127 
129 }
130 
131 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
132 void UzawaSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector &X, const MultiVector &B, bool /* InitialGuessIsZero */) const {
133  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::UzawaSmoother::Apply(): Setup() has not been called");
134 
135  Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
136 
138 
139  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
140  bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
141 
142  // extract parameters from internal parameter list
144  LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
145  Scalar omega = pL.get<Scalar>("Damping factor");
146 
147  // wrap current solution vector in RCP
148  RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
149  RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
150 
151  // make sure that both rcpX and rcpB are BlockedMultiVector objects
152  bool bCopyResultX = false;
153  bool bReorderX = false;
154  bool bReorderB = false;
155  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
156  MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
157  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
158  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
159 
160  // check the type of operator
163 
164  if (rbA.is_null() == false) {
165  // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
166  // map for the construction of the blocked multivectors
167 
168  // check type of X vector
169  if (bX.is_null() == true) {
170  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
171  rcpX.swap(vectorWithBlockedMap);
172  bCopyResultX = true;
173  bReorderX = true;
174  }
175 
176  // check type of B vector
177  if (bB.is_null() == true) {
178  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
179  rcpB.swap(vectorWithBlockedMap);
180  bReorderB = true;
181  }
182  } else {
183  // A is a BlockedCrsMatrix and uses a plain blocked map
184  if (bX.is_null() == true) {
185  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
186  rcpX.swap(vectorWithBlockedMap);
187  bCopyResultX = true;
188  }
189 
190  if (bB.is_null() == true) {
191  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(), rcpB));
192  rcpB.swap(vectorWithBlockedMap);
193  }
194  }
195 
196  // we now can guarantee that X and B are blocked multi vectors
197  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
198  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
199 
200  // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
201  if (rbA.is_null() == false) {
202  Teuchos::RCP<const Xpetra::BlockReorderManager> brm = rbA->getBlockReorderManager();
203 
204  // check type of X vector
205  if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
206  // X is a blocked multi vector but incompatible to the reordered blocked operator A
207  Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
208  rcpX.swap(reorderedBlockedVector);
209  }
210 
211  if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
212  // B is a blocked multi vector but incompatible to the reordered blocked operator A
213  Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
214  rcpB.swap(reorderedBlockedVector);
215  }
216  }
217 
218  // Throughout the rest of the algorithm rcpX and rcpB are used for solution vector and RHS
219 
220  // create residual vector
221  // contains current residual of current solution X with rhs B
222  RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
223  RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
224  Teuchos::RCP<MultiVector> r1 = bresidual->getMultiVector(0, bRangeThyraMode);
225  Teuchos::RCP<MultiVector> r2 = bresidual->getMultiVector(1, bRangeThyraMode);
226 
227  // helper vector 1
228  RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
229  RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
230  RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0, bDomainThyraMode);
231  RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1, bDomainThyraMode);
232 
233  // incrementally improve solution vector X
234  for (LocalOrdinal run = 0; run < nSweeps; ++run) {
235  // 1) calculate current residual
236  residual->update(one, *rcpB, zero); // residual = B
237  A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
238 
239  // 2) solve F * \Delta \tilde{x}_1 = r_1
240  // start with zero guess \Delta \tilde{x}_1
241  bxtilde->putScalar(zero);
242  velPredictSmoo_->Apply(*xtilde1, *r1);
243 
244  // 3) calculate rhs for SchurComp equation
245  // r_2 - D \Delta \tilde{x}_1
246  RCP<MultiVector> schurCompRHS = MultiVectorFactory::Build(r2->getMap(), rcpB->getNumVectors());
247  D_->apply(*xtilde1, *schurCompRHS);
248  schurCompRHS->update(one, *r2, -one);
249 
250  // 4) solve SchurComp equation
251  // start with zero guess \Delta \tilde{x}_2
252  schurCompSmoo_->Apply(*xtilde2, *schurCompRHS);
253 
254  rcpX->update(omega, *bxtilde, one);
255  }
256 
257  if (bCopyResultX == true) {
258  RCP<MultiVector> Xmerged = bX->Merge();
259  X.update(one, *Xmerged, zero);
260  }
261 }
262 
263 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
266  return rcp(new UzawaSmoother(*this));
267 }
268 
269 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271  std::ostringstream out;
273  out << "{type = " << type_ << "}";
274  return out.str();
275 }
276 
277 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
280 
281  if (verbLevel & Parameters0) {
282  out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
283  }
284 
285  if (verbLevel & Debug) {
286  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
287  }
288 }
289 
290 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
292  // FIXME: This is a placeholder
294 }
295 
296 } // namespace MueLu
297 
298 #endif /* MUELU_UZAWASMOOTHER_DEF_HPP_ */
Important warning messages (one line)
virtual const Teuchos::ParameterList & GetParameterList() const
Exception indicating invalid cast attempted.
void Setup(Level &currentLevel)
Setup routine.
MueLu::DefaultLocalOrdinal LocalOrdinal
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
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 Apply(MultiVector &X, MultiVector const &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
RCP< const ParameterList > GetValidParameterList() const
Input.
void SetVelocityPredictionFactoryManager(RCP< FactoryManager > FactManager)
Set factory manager for internal velocity prediction.
virtual const RCP< const FactoryBase > GetFactory(const std::string &varName) const =0
Get.
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
Block triangle Uzawa smoother for 2x2 block matrices.
std::string description() const
Return a simple one-line description of this object.
RCP< SmootherPrototype > Copy() const
bool IsSetup() const
Get the state of a smoother prototype.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Print class parameters.
void SetSchurCompFactoryManager(RCP< FactoryManager > FactManager)
Set factory manager for internal SchurComplement handling.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Scalar SC
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos)
Add a factory manager at a specific position.
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)
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual std::string description() const
Return a simple one-line description of this object.
std::string toString(const T &t)
bool is_null() const