MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_BraessSarazinSmoother_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_BRAESSSARAZINSMOOTHER_DEF_HPP_
11 #define MUELU_BRAESSSARAZINSMOOTHER_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_Utilities.hpp"
27 #include "MueLu_Monitor.hpp"
28 #include "MueLu_HierarchyUtils.hpp"
29 #include "MueLu_SmootherBase.hpp"
30 
31 #include "MueLu_FactoryManager.hpp"
32 
33 namespace MueLu {
34 
35 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37  TEUCHOS_TEST_FOR_EXCEPTION(pos != 0, Exceptions::RuntimeError, "MueLu::BraessSarazinSmoother::AddFactoryManager: parameter \'pos\' must be zero! error.");
38  FactManager_ = FactManager;
39 }
40 
41 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43  RCP<ParameterList> validParamList = rcp(new ParameterList());
44 
46 
47  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A");
48 
49  validParamList->set<bool>("lumping", false,
50  "Use lumping to construct diag(A(0,0)), "
51  "i.e. use row sum of the abs values on the diagonal as approximation of A00 (and A00^{-1})");
52  validParamList->set<SC>("Damping factor", one, "Damping/Scaling factor in BraessSarazin (usually has to be chosen > 1, default = 1.0)");
53  validParamList->set<LO>("Sweeps", 1, "Number of BraessSarazin sweeps (default = 1)");
54  validParamList->set<bool>("q2q1 mode", false, "Run in the mode matching Q2Q1 matlab code");
55 
56  return validParamList;
57 }
58 
59 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61  this->Input(currentLevel, "A");
62 
64  "MueLu::BraessSarazinSmoother::DeclareInput: FactManager_ must not be null! "
65  "Introduce a FactoryManager for the SchurComplement equation.");
66 
67  // carefully call DeclareInput after switching to internal FactoryManager
68  {
69  SetFactoryManager currentSFM(rcpFromRef(currentLevel), FactManager_);
70 
71  // request "Smoother" for current subblock row.
72  currentLevel.DeclareInput("PreSmoother", FactManager_->GetFactory("Smoother").get());
73 
74  // request Schur matrix just in case
75  currentLevel.DeclareInput("A", FactManager_->GetFactory("A").get());
76  }
77 }
78 
79 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81  FactoryMonitor m(*this, "Setup Smoother", currentLevel);
82 
83  if (SmootherPrototype::IsSetup() == true)
84  this->GetOStream(Warnings0) << "MueLu::BreaessSarazinSmoother::Setup(): Setup() has already been called";
85 
86  // Extract blocked operator A from current level
87  A_ = Factory::Get<RCP<Matrix> >(currentLevel, "A");
88  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
90  "MueLu::BraessSarazinSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
91 
92  // Store map extractors
93  rangeMapExtractor_ = bA->getRangeMapExtractor();
94  domainMapExtractor_ = bA->getDomainMapExtractor();
95 
96  // Store the blocks in local member variables
97  A00_ = bA->getMatrix(0, 0);
98  A01_ = bA->getMatrix(0, 1);
99  A10_ = bA->getMatrix(1, 0);
100  A11_ = bA->getMatrix(1, 1);
101 
103  SC omega = pL.get<SC>("Damping factor");
104 
105  // Create the inverse of the diagonal of F
106  // TODO add safety check for zeros on diagonal of F!
107  RCP<Vector> diagFVector = VectorFactory::Build(A00_->getRowMap());
108  if (pL.get<bool>("lumping") == false) {
109  A00_->getLocalDiagCopy(*diagFVector); // extract diagonal of F
110  } else {
111  diagFVector = Utilities::GetLumpedMatrixDiagonal(*A00_);
112  }
113  diagFVector->scale(omega);
114  D_ = Utilities::GetInverse(diagFVector);
115 
116  // check whether D_ is a blocked vector with only 1 block
117  RCP<BlockedVector> bD = Teuchos::rcp_dynamic_cast<BlockedVector>(D_);
118  if (bD.is_null() == false && bD->getBlockedMap()->getNumMaps() == 1) {
119  RCP<Vector> nestedVec = bD->getMultiVector(0, bD->getBlockedMap()->getThyraMode())->getVectorNonConst(0);
120  D_.swap(nestedVec);
121  }
122 
123  // Set the Smoother
124  // carefully switch to the SubFactoryManagers (defined by the users)
125  {
126  SetFactoryManager currentSFM(rcpFromRef(currentLevel), FactManager_);
127  smoo_ = currentLevel.Get<RCP<SmootherBase> >("PreSmoother", FactManager_->GetFactory("Smoother").get());
128  S_ = currentLevel.Get<RCP<Matrix> >("A", FactManager_->GetFactory("A").get());
129  }
130 
132 }
133 
134 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
135 void BraessSarazinSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
137  "MueLu::BraessSarazinSmoother::Apply(): Setup() has not been called");
138 
139  RCP<MultiVector> rcpX = rcpFromRef(X);
140  RCP<const MultiVector> rcpB = rcpFromRef(B);
141 
142  // make sure that both rcpX and rcpB are BlockedMultiVector objects
143  bool bCopyResultX = false;
144  bool bReorderX = false;
145  bool bReorderB = false;
146  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
147  MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
148  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
149  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
150 
151  // check the type of operator
154 
155  if (rbA.is_null() == false) {
156  // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
157  // map for the construction of the blocked multivectors
158 
159  // check type of X vector
160  if (bX.is_null() == true) {
161  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
162  rcpX.swap(vectorWithBlockedMap);
163  bCopyResultX = true;
164  bReorderX = true;
165  }
166 
167  // check type of B vector
168  if (bB.is_null() == true) {
169  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
170  rcpB.swap(vectorWithBlockedMap);
171  bReorderB = true;
172  }
173  } else {
174  // A is a BlockedCrsMatrix and uses a plain blocked map
175  if (bX.is_null() == true) {
176  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
177  rcpX.swap(vectorWithBlockedMap);
178  bCopyResultX = true;
179  }
180 
181  if (bB.is_null() == true) {
182  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(), rcpB));
183  rcpB.swap(vectorWithBlockedMap);
184  }
185  }
186 
187  // we now can guarantee that X and B are blocked multi vectors
188  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
189  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
190 
191  // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
192  if (rbA.is_null() == false) {
193  Teuchos::RCP<const Xpetra::BlockReorderManager> brm = rbA->getBlockReorderManager();
194 
195  // check type of X vector
196  if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
197  // X is a blocked multi vector but incompatible to the reordered blocked operator A
198  Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
199  rcpX.swap(reorderedBlockedVector);
200  }
201 
202  if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
203  // B is a blocked multi vector but incompatible to the reordered blocked operator A
204  Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
205  rcpB.swap(reorderedBlockedVector);
206  }
207  }
208 
209  // use the GIDs of the sub blocks
210  // This is valid as the subblocks actually represent the GIDs (either Thyra, Xpetra or pseudo Xpetra)
211 
212  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
213  bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
214 
215  RCP<MultiVector> deltaX = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
216  RCP<BlockedMultiVector> bdeltaX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(deltaX);
217  RCP<MultiVector> deltaX0 = bdeltaX->getMultiVector(0, bDomainThyraMode);
218  RCP<MultiVector> deltaX1 = bdeltaX->getMultiVector(1, bDomainThyraMode);
219 
220  RCP<MultiVector> Rtmp = rangeMapExtractor_->getVector(1, rcpB->getNumVectors(), bRangeThyraMode);
221 
222  typedef Teuchos::ScalarTraits<SC> STS;
223  SC one = STS::one(), zero = STS::zero();
224 
225  // extract parameters from internal parameter list
227  LO nSweeps = pL.get<LO>("Sweeps");
228 
229  RCP<MultiVector> R; // = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());;
230  if (InitialGuessIsZero) {
231  R = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
232  R->update(one, *rcpB, zero);
233  } else {
234  R = Utilities::Residual(*A_, *rcpX, *rcpB);
235  }
236 
237  // extract diagonal of Schur complement operator
238  RCP<Vector> diagSVector = VectorFactory::Build(S_->getRowMap());
239  S_->getLocalDiagCopy(*diagSVector);
240 
241  for (LO run = 0; run < nSweeps; ++run) {
242  // Extract corresponding subvectors from X and R
243  // Automatically detect whether we use Thyra or Xpetra GIDs
244  // The GIDs should be always compatible with the GIDs of A00, A01, etc...
245  RCP<MultiVector> R0 = rangeMapExtractor_->ExtractVector(R, 0, bRangeThyraMode);
246  RCP<MultiVector> R1 = rangeMapExtractor_->ExtractVector(R, 1, bRangeThyraMode);
247 
248  // Calculate Rtmp = R1 - D * deltaX0 (equation 8.14)
249  deltaX0->putScalar(zero);
250  deltaX0->elementWiseMultiply(one, *D_, *R0, zero); // deltaX0 = D * R0 (equation 8.13)
251  A10_->apply(*deltaX0, *Rtmp); // Rtmp = A10*D*deltaX0 (intermediate step)
252  Rtmp->update(one, *R1, -one); // Rtmp = R1 - A10*D*deltaX0
253 
254  if (!pL.get<bool>("q2q1 mode")) {
255  deltaX1->putScalar(zero);
256  } else {
257  // special code for q2q1
258  if (Teuchos::rcp_dynamic_cast<BlockedVector>(diagSVector) == Teuchos::null) {
259  ArrayRCP<SC> Sdiag = diagSVector->getDataNonConst(0);
260  ArrayRCP<SC> deltaX1data = deltaX1->getDataNonConst(0);
261  ArrayRCP<SC> Rtmpdata = Rtmp->getDataNonConst(0);
262  for (GO row = 0; row < deltaX1data.size(); row++)
263  deltaX1data[row] = Teuchos::as<SC>(1.1) * Rtmpdata[row] / Sdiag[row];
264  } else {
265  TEUCHOS_TEST_FOR_EXCEPTION(true, MueLu::Exceptions::RuntimeError, "MueLu::BraessSarazinSmoother: q2q1 mode only supported for non-blocked operators.")
266  }
267  }
268 
269  smoo_->Apply(*deltaX1, *Rtmp);
270 
271  // Compute deltaX0
272  deltaX0->putScalar(zero); // just for safety
273  A01_->apply(*deltaX1, *deltaX0); // deltaX0 = A01*deltaX1
274  deltaX0->update(one, *R0, -one); // deltaX0 = R0 - A01*deltaX1
275  R0.swap(deltaX0);
276  deltaX0->elementWiseMultiply(one, *D_, *R0, zero); // deltaX0 = D*(R0 - A01*deltaX1)
277 
278  RCP<MultiVector> X0 = domainMapExtractor_->ExtractVector(rcpX, 0, bDomainThyraMode);
279  RCP<MultiVector> X1 = domainMapExtractor_->ExtractVector(rcpX, 1, bDomainThyraMode);
280 
281  // Update solution
282  X0->update(one, *deltaX0, one);
283  X1->update(one, *deltaX1, one);
284 
285  domainMapExtractor_->InsertVector(X0, 0, rcpX, bDomainThyraMode);
286  domainMapExtractor_->InsertVector(X1, 1, rcpX, bDomainThyraMode);
287 
288  if (run < nSweeps - 1) {
289  R = Utilities::Residual(*A_, *rcpX, *rcpB);
290  }
291  }
292 
293  if (bCopyResultX == true) {
294  RCP<MultiVector> Xmerged = bX->Merge();
295  X.update(one, *Xmerged, zero);
296  }
297 }
298 
299 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
302  return rcp(new BraessSarazinSmoother(*this));
303 }
304 
305 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
307  description() const {
308  std::ostringstream out;
310  out << "{type = " << type_ << "}";
311  return out.str();
312 }
313 
314 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
317 
318  if (verbLevel & Parameters0) {
319  out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
320  }
321 
322  if (verbLevel & Debug) {
323  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
324  }
325 }
326 
327 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
329  // FIXME: This is a placeholder
331 }
332 
333 } // namespace MueLu
334 
335 #endif /* MUELU_BRAESSSARAZINSMOOTHER_DEF_HPP_ */
Important warning messages (one line)
virtual const Teuchos::ParameterList & GetParameterList() const
Exception indicating invalid cast attempted.
BraessSarazin smoother for 2x2 block matrices.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
GlobalOrdinal GO
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.
LocalOrdinal LO
RCP< SmootherPrototype > Copy() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
size_type size() const
void DeclareInput(Level &currentLevel) const
Input.
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
bool IsSetup() const
Get the state of a smoother prototype.
RCP< const ParameterList > GetValidParameterList() const
Input.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Print class parameters.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos=0)
Add a factory manager for BraessSarazin internal SchurComplement handling.
Scalar SC
void Setup(Level &currentLevel)
Setup routine.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the Braess Sarazin smoother.
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.
std::string description() const
Return a simple one-line description of this object.
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()
virtual std::string description() const
Return a simple one-line description of this object.
std::string toString(const T &t)
bool is_null() const