MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_IndefBlockedDiagonalSmoother_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 /*
11  * MueLu_IndefBlockedDiagonalSmoother_def.hpp
12  *
13  * Created on: 13 May 2014
14  * Author: wiesner
15  */
16 
17 #ifndef MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_
18 #define MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_
19 
20 #include "Teuchos_ArrayViewDecl.hpp"
21 #include "Teuchos_ScalarTraits.hpp"
22 
23 #include "MueLu_ConfigDefs.hpp"
24 
25 #include <Xpetra_Matrix.hpp>
26 #include <Xpetra_CrsMatrixWrap.hpp>
28 #include <Xpetra_MultiVectorFactory.hpp>
29 #include <Xpetra_VectorFactory.hpp>
31 
33 #include "MueLu_Level.hpp"
34 #include "MueLu_Monitor.hpp"
35 #include "MueLu_HierarchyUtils.hpp"
36 #include "MueLu_SmootherBase.hpp"
37 
38 // include files for default FactoryManager
39 #include "MueLu_SchurComplementFactory.hpp"
40 #include "MueLu_FactoryManager.hpp"
41 
42 namespace MueLu {
43 
44 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
46  : type_("IndefiniteBlockDiagonalSmoother")
47  , A_(Teuchos::null) {
48 }
49 
50 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
52 
53 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
55  RCP<ParameterList> validParamList = rcp(new ParameterList());
56 
57  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A (must be a 2x2 block matrix)");
58  validParamList->set<Scalar>("Damping factor", 1.0, "Damping/Scaling factor");
59  validParamList->set<LocalOrdinal>("Sweeps", 1, "Number of SIMPLE sweeps (default = 1)");
60  // validParamList->set< bool > ("UseSIMPLEC", false, "Use SIMPLEC instead of SIMPLE (default = false)");
61 
62  return validParamList;
63 }
64 
66 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68  TEUCHOS_TEST_FOR_EXCEPTION(pos < 0, Exceptions::RuntimeError, "MueLu::IndefBlockedDiagonalSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
69 
70  size_t myPos = Teuchos::as<size_t>(pos);
71 
72  if (myPos < FactManager_.size()) {
73  // replace existing entris in FactManager_ vector
74  FactManager_.at(myPos) = FactManager;
75  } else if (myPos == FactManager_.size()) {
76  // add new Factory manager in the end of the vector
77  FactManager_.push_back(FactManager);
78  } else { // if(myPos > FactManager_.size())
79  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
80  *out << "Warning: cannot add new FactoryManager at proper position " << pos << ". The FactoryManager is just appended to the end. Check this!" << std::endl;
81 
82  // add new Factory manager in the end of the vector
83  FactManager_.push_back(FactManager);
84  }
85 }
86 
87 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
89  currentLevel.DeclareInput("A", this->GetFactory("A").get());
90 
91  TEUCHOS_TEST_FOR_EXCEPTION(FactManager_.size() != 2, Exceptions::RuntimeError, "MueLu::IndefBlockedDiagonalSmoother::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\"!");
92 
93  // loop over all factory managers for the subblocks of blocked operator A
94  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
95  for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
96  SetFactoryManager currentSFM(rcpFromRef(currentLevel), *it);
97 
98  // request "Smoother" for current subblock row.
99  currentLevel.DeclareInput("PreSmoother", (*it)->GetFactory("Smoother").get());
100  }
101 }
102 
103 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
105  FactoryMonitor m(*this, "Setup for indefinite blocked diagonal smoother", currentLevel);
106 
107  if (SmootherPrototype::IsSetup() == true)
108  this->GetOStream(Warnings0) << "MueLu::IndefBlockedDiagonalSmoother::Setup(): Setup() has already been called";
109 
110  // extract blocked operator A from current level
111  A_ = Factory::Get<RCP<Matrix> >(currentLevel, "A");
112 
113  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
114  TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::IndefBlockedDiagonalSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
115 
116  // store map extractors
117  rangeMapExtractor_ = bA->getRangeMapExtractor();
118  domainMapExtractor_ = bA->getDomainMapExtractor();
119 
120  // Store the blocks in local member variables
121  Teuchos::RCP<Matrix> A00 = bA->getMatrix(0, 0);
122  Teuchos::RCP<Matrix> A01 = bA->getMatrix(0, 1);
123  Teuchos::RCP<Matrix> A10 = bA->getMatrix(1, 0);
124  Teuchos::RCP<Matrix> A11 = bA->getMatrix(1, 1);
125 
126  F_ = A00;
127  Z_ = A11;
128 
129  /*const ParameterList & pL = Factory::GetParameterList();
130  bool bSIMPLEC = pL.get<bool>("UseSIMPLEC");
131 
132  // Create the inverse of the diagonal of F
133  RCP<Vector> diagFVector = VectorFactory::Build(F_->getRowMap());
134  if(!bSIMPLEC) {
135  F_->getLocalDiagCopy(*diagFVector); // extract diagonal of F
136  diagFVector->reciprocal(*diagFVector); // build reciprocal
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  diagFVector->reciprocal(*diagFVector); // build reciprocal
152  }
153  diagFinv_ = diagFVector;*/
154 
155  // Set the Smoother
156  // carefully switch to the SubFactoryManagers (defined by the users)
157  {
158  RCP<const FactoryManagerBase> velpredictFactManager = FactManager_.at(0);
159  SetFactoryManager currentSFM(rcpFromRef(currentLevel), velpredictFactManager);
160  velPredictSmoo_ = currentLevel.Get<RCP<SmootherBase> >("PreSmoother", velpredictFactManager->GetFactory("Smoother").get());
161  }
162  {
163  RCP<const FactoryManagerBase> schurFactManager = FactManager_.at(1);
164  SetFactoryManager currentSFM(rcpFromRef(currentLevel), schurFactManager);
165  schurCompSmoo_ = currentLevel.Get<RCP<SmootherBase> >("PreSmoother", schurFactManager->GetFactory("Smoother").get());
166  }
168 }
169 
170 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
171 void IndefBlockedDiagonalSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector &X, const MultiVector &B, bool /* InitialGuessIsZero */) const {
172  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::IndefBlockedDiagonalSmoother::Apply(): Setup() has not been called");
173 
174  Teuchos::RCP<Teuchos::FancyOStream> fos = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
175 
177 
178  // extract parameters from internal parameter list
180  LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
181  Scalar omega = pL.get<Scalar>("Damping factor");
182 
183  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
184  bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
185 
186  // wrap current solution vector in RCP
187  RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
188  RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
189 
190  // make sure that both rcpX and rcpB are BlockedMultiVector objects
191  bool bCopyResultX = false;
192  bool bReorderX = false;
193  bool bReorderB = false;
194  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
195  MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
196  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
197  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
198 
199  // check the type of operator
202 
203  if (rbA.is_null() == false) {
204  // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
205  // map for the construction of the blocked multivectors
206 
207  // check type of X vector
208  if (bX.is_null() == true) {
209  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
210  rcpX.swap(vectorWithBlockedMap);
211  bCopyResultX = true;
212  bReorderX = true;
213  }
214 
215  // check type of B vector
216  if (bB.is_null() == true) {
217  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
218  rcpB.swap(vectorWithBlockedMap);
219  bReorderB = true;
220  }
221  } else {
222  // A is a BlockedCrsMatrix and uses a plain blocked map
223  if (bX.is_null() == true) {
224  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
225  rcpX.swap(vectorWithBlockedMap);
226  bCopyResultX = true;
227  }
228 
229  if (bB.is_null() == true) {
230  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(), rcpB));
231  rcpB.swap(vectorWithBlockedMap);
232  }
233  }
234 
235  // we now can guarantee that X and B are blocked multi vectors
236  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
237  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
238 
239  // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
240  if (rbA.is_null() == false) {
241  Teuchos::RCP<const Xpetra::BlockReorderManager> brm = rbA->getBlockReorderManager();
242 
243  // check type of X vector
244  if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
245  // X is a blocked multi vector but incompatible to the reordered blocked operator A
246  Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
247  rcpX.swap(reorderedBlockedVector);
248  }
249 
250  if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
251  // B is a blocked multi vector but incompatible to the reordered blocked operator A
252  Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
253  rcpB.swap(reorderedBlockedVector);
254  }
255  }
256 
257  // Throughout the rest of the algorithm rcpX and rcpB are used for solution vector and RHS
258 
259  // create residual vector
260  // contains current residual of current solution X with rhs B
261  RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
262  RCP<BlockedMultiVector> bresidual = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(residual);
263  Teuchos::RCP<MultiVector> r1 = bresidual->getMultiVector(0, bRangeThyraMode);
264  Teuchos::RCP<MultiVector> r2 = bresidual->getMultiVector(1, bRangeThyraMode);
265 
266  // helper vector 1
267  RCP<MultiVector> xtilde = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
268  RCP<BlockedMultiVector> bxtilde = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(xtilde);
269  RCP<MultiVector> xtilde1 = bxtilde->getMultiVector(0, bDomainThyraMode);
270  RCP<MultiVector> xtilde2 = bxtilde->getMultiVector(1, bDomainThyraMode);
271 
272  // incrementally improve solution vector X
273  for (LocalOrdinal run = 0; run < nSweeps; ++run) {
274  // 1) calculate current residual
275  residual->update(one, *rcpB, zero); // residual = B
276  A_->apply(*rcpX, *residual, Teuchos::NO_TRANS, -one, one);
277 
278  // split residual vector
279  // Teuchos::RCP<MultiVector> r1 = rangeMapExtractor_->ExtractVector(residual, 0, bRangeThyraMode);
280  // Teuchos::RCP<MultiVector> r2 = rangeMapExtractor_->ExtractVector(residual, 1, bRangeThyraMode);
281 
282  // 2) solve F * \Delta \tilde{x}_1 = r_1
283  // start with zero guess \Delta \tilde{x}_1
284  // RCP<MultiVector> xtilde1 = MultiVectorFactory::Build(r1->getMap(),rcpX->getNumVectors(),true);
285  // RCP<MultiVector> xtilde2 = MultiVectorFactory::Build(r2->getMap(),rcpX->getNumVectors(),true);
286  bxtilde->putScalar(zero);
287  velPredictSmoo_->Apply(*xtilde1, *r1);
288 
289  // 3) solve SchurComp equation
290  // start with zero guess \Delta \tilde{x}_2
291  schurCompSmoo_->Apply(*xtilde2, *r2);
292 #if 1
293  // 4) update solution vector
294  rcpX->update(omega, *bxtilde, one);
295 #else
296  Teuchos::RCP<MultiVector> x1 = domainMapExtractor_->ExtractVector(rcpX, 0, bDomainThyraMode);
297  Teuchos::RCP<MultiVector> x2 = domainMapExtractor_->ExtractVector(rcpX, 1, bDomainThyraMode);
298 
299  // 5) update solution vector with increments xhat1 and xhat2
300  // rescale increment for x2 with omega_
301  x1->update(omega, *xtilde1, one); // x1 = x1_old + omega xtilde1
302  x2->update(omega, *xtilde2, one); // x2 = x2_old + omega xtilde2
303 
304  // write back solution in global vector X
305  domainMapExtractor_->InsertVector(x1, 0, rcpX, bDomainThyraMode);
306  domainMapExtractor_->InsertVector(x2, 1, rcpX, bDomainThyraMode);
307 #endif
308  }
309 
310  if (bCopyResultX == true) {
311  RCP<MultiVector> Xmerged = bX->Merge();
312  X.update(one, *Xmerged, zero);
313  }
314 }
315 
316 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
319  return rcp(new IndefBlockedDiagonalSmoother(*this));
320 }
321 
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
324  std::ostringstream out;
326  out << "{type = " << type_ << "}";
327  return out.str();
328 }
329 
330 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
333 
334  if (verbLevel & Parameters0) {
335  out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
336  }
337 
338  if (verbLevel & Debug) {
339  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
340  }
341 }
342 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
344  // FIXME: This is a placeholder
346 }
347 
348 } // namespace MueLu
349 
350 #endif /* MUELU_INDEFBLOCKEDDIAGONALSMOOTHER_DEF_HPP_ */
Important warning messages (one line)
virtual const Teuchos::ParameterList & GetParameterList() const
Exception indicating invalid cast attempted.
std::string description() const
Return a simple one-line description of this object.
MueLu::DefaultLocalOrdinal LocalOrdinal
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
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.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void DeclareInput(Level &currentLevel) const
Input.
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
Cheap Blocked diagonal smoother for indefinite 2x2 block matrices.
RCP< const ParameterList > GetValidParameterList() const
Input.
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.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Scalar SC
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos=0)
Add a factory manager for BraessSarazin internal SchurComplement handling.
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