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 // ***********************************************************************
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_BRAESSSARAZINSMOOTHER_DEF_HPP_
48 #define MUELU_BRAESSSARAZINSMOOTHER_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>
57 #include <Xpetra_MultiVectorFactory.hpp>
58 #include <Xpetra_VectorFactory.hpp>
60 
62 #include "MueLu_Level.hpp"
63 #include "MueLu_Utilities.hpp"
64 #include "MueLu_Monitor.hpp"
65 #include "MueLu_HierarchyUtils.hpp"
66 #include "MueLu_SmootherBase.hpp"
67 
68 #include "MueLu_FactoryManager.hpp"
69 
70 namespace MueLu {
71 
72 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74  TEUCHOS_TEST_FOR_EXCEPTION(pos != 0, Exceptions::RuntimeError, "MueLu::BraessSarazinSmoother::AddFactoryManager: parameter \'pos\' must be zero! error.");
75  FactManager_ = FactManager;
76 }
77 
78 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80  RCP<ParameterList> validParamList = rcp(new ParameterList());
81 
83 
84  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A");
85 
86  validParamList->set<bool>("lumping", false,
87  "Use lumping to construct diag(A(0,0)), "
88  "i.e. use row sum of the abs values on the diagonal as approximation of A00 (and A00^{-1})");
89  validParamList->set<SC>("Damping factor", one, "Damping/Scaling factor in BraessSarazin (usually has to be chosen > 1, default = 1.0)");
90  validParamList->set<LO>("Sweeps", 1, "Number of BraessSarazin sweeps (default = 1)");
91  validParamList->set<bool>("q2q1 mode", false, "Run in the mode matching Q2Q1 matlab code");
92 
93  return validParamList;
94 }
95 
96 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98  this->Input(currentLevel, "A");
99 
101  "MueLu::BraessSarazinSmoother::DeclareInput: FactManager_ must not be null! "
102  "Introduce a FactoryManager for the SchurComplement equation.");
103 
104  // carefully call DeclareInput after switching to internal FactoryManager
105  {
106  SetFactoryManager currentSFM(rcpFromRef(currentLevel), FactManager_);
107 
108  // request "Smoother" for current subblock row.
109  currentLevel.DeclareInput("PreSmoother", FactManager_->GetFactory("Smoother").get());
110 
111  // request Schur matrix just in case
112  currentLevel.DeclareInput("A", FactManager_->GetFactory("A").get());
113  }
114 }
115 
116 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
118  FactoryMonitor m(*this, "Setup Smoother", currentLevel);
119 
120  if (SmootherPrototype::IsSetup() == true)
121  this->GetOStream(Warnings0) << "MueLu::BreaessSarazinSmoother::Setup(): Setup() has already been called";
122 
123  // Extract blocked operator A from current level
124  A_ = Factory::Get<RCP<Matrix> >(currentLevel, "A");
125  RCP<BlockedCrsMatrix> bA = rcp_dynamic_cast<BlockedCrsMatrix>(A_);
127  "MueLu::BraessSarazinSmoother::Setup: input matrix A is not of type BlockedCrsMatrix! error.");
128 
129  // Store map extractors
130  rangeMapExtractor_ = bA->getRangeMapExtractor();
131  domainMapExtractor_ = bA->getDomainMapExtractor();
132 
133  // Store the blocks in local member variables
134  A00_ = bA->getMatrix(0, 0);
135  A01_ = bA->getMatrix(0, 1);
136  A10_ = bA->getMatrix(1, 0);
137  A11_ = bA->getMatrix(1, 1);
138 
140  SC omega = pL.get<SC>("Damping factor");
141 
142  // Create the inverse of the diagonal of F
143  // TODO add safety check for zeros on diagonal of F!
144  RCP<Vector> diagFVector = VectorFactory::Build(A00_->getRowMap());
145  if (pL.get<bool>("lumping") == false) {
146  A00_->getLocalDiagCopy(*diagFVector); // extract diagonal of F
147  } else {
148  diagFVector = Utilities::GetLumpedMatrixDiagonal(*A00_);
149  }
150  diagFVector->scale(omega);
151  D_ = Utilities::GetInverse(diagFVector);
152 
153  // check whether D_ is a blocked vector with only 1 block
154  RCP<BlockedVector> bD = Teuchos::rcp_dynamic_cast<BlockedVector>(D_);
155  if (bD.is_null() == false && bD->getBlockedMap()->getNumMaps() == 1) {
156  RCP<Vector> nestedVec = bD->getMultiVector(0, bD->getBlockedMap()->getThyraMode())->getVectorNonConst(0);
157  D_.swap(nestedVec);
158  }
159 
160  // Set the Smoother
161  // carefully switch to the SubFactoryManagers (defined by the users)
162  {
163  SetFactoryManager currentSFM(rcpFromRef(currentLevel), FactManager_);
164  smoo_ = currentLevel.Get<RCP<SmootherBase> >("PreSmoother", FactManager_->GetFactory("Smoother").get());
165  S_ = currentLevel.Get<RCP<Matrix> >("A", FactManager_->GetFactory("A").get());
166  }
167 
169 }
170 
171 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
172 void BraessSarazinSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
174  "MueLu::BraessSarazinSmoother::Apply(): Setup() has not been called");
175 
176  RCP<MultiVector> rcpX = rcpFromRef(X);
177  RCP<const MultiVector> rcpB = rcpFromRef(B);
178 
179  // make sure that both rcpX and rcpB are BlockedMultiVector objects
180  bool bCopyResultX = false;
181  bool bReorderX = false;
182  bool bReorderB = false;
183  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
184  MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
185  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
186  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
187 
188  // check the type of operator
191 
192  if (rbA.is_null() == false) {
193  // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
194  // map for the construction of the blocked multivectors
195 
196  // check type of X vector
197  if (bX.is_null() == true) {
198  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
199  rcpX.swap(vectorWithBlockedMap);
200  bCopyResultX = true;
201  bReorderX = true;
202  }
203 
204  // check type of B vector
205  if (bB.is_null() == true) {
206  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
207  rcpB.swap(vectorWithBlockedMap);
208  bReorderB = true;
209  }
210  } else {
211  // A is a BlockedCrsMatrix and uses a plain blocked map
212  if (bX.is_null() == true) {
213  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
214  rcpX.swap(vectorWithBlockedMap);
215  bCopyResultX = true;
216  }
217 
218  if (bB.is_null() == true) {
219  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(), rcpB));
220  rcpB.swap(vectorWithBlockedMap);
221  }
222  }
223 
224  // we now can guarantee that X and B are blocked multi vectors
225  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
226  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
227 
228  // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
229  if (rbA.is_null() == false) {
230  Teuchos::RCP<const Xpetra::BlockReorderManager> brm = rbA->getBlockReorderManager();
231 
232  // check type of X vector
233  if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
234  // X is a blocked multi vector but incompatible to the reordered blocked operator A
235  Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
236  rcpX.swap(reorderedBlockedVector);
237  }
238 
239  if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
240  // B is a blocked multi vector but incompatible to the reordered blocked operator A
241  Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
242  rcpB.swap(reorderedBlockedVector);
243  }
244  }
245 
246  // use the GIDs of the sub blocks
247  // This is valid as the subblocks actually represent the GIDs (either Thyra, Xpetra or pseudo Xpetra)
248 
249  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
250  bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
251 
252  RCP<MultiVector> deltaX = MultiVectorFactory::Build(rcpX->getMap(), rcpX->getNumVectors());
253  RCP<BlockedMultiVector> bdeltaX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(deltaX);
254  RCP<MultiVector> deltaX0 = bdeltaX->getMultiVector(0, bDomainThyraMode);
255  RCP<MultiVector> deltaX1 = bdeltaX->getMultiVector(1, bDomainThyraMode);
256 
257  RCP<MultiVector> Rtmp = rangeMapExtractor_->getVector(1, rcpB->getNumVectors(), bRangeThyraMode);
258 
259  typedef Teuchos::ScalarTraits<SC> STS;
260  SC one = STS::one(), zero = STS::zero();
261 
262  // extract parameters from internal parameter list
264  LO nSweeps = pL.get<LO>("Sweeps");
265 
266  RCP<MultiVector> R; // = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());;
267  if (InitialGuessIsZero) {
268  R = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
269  R->update(one, *rcpB, zero);
270  } else {
271  R = Utilities::Residual(*A_, *rcpX, *rcpB);
272  }
273 
274  // extract diagonal of Schur complement operator
275  RCP<Vector> diagSVector = VectorFactory::Build(S_->getRowMap());
276  S_->getLocalDiagCopy(*diagSVector);
277 
278  for (LO run = 0; run < nSweeps; ++run) {
279  // Extract corresponding subvectors from X and R
280  // Automatically detect whether we use Thyra or Xpetra GIDs
281  // The GIDs should be always compatible with the GIDs of A00, A01, etc...
282  RCP<MultiVector> R0 = rangeMapExtractor_->ExtractVector(R, 0, bRangeThyraMode);
283  RCP<MultiVector> R1 = rangeMapExtractor_->ExtractVector(R, 1, bRangeThyraMode);
284 
285  // Calculate Rtmp = R1 - D * deltaX0 (equation 8.14)
286  deltaX0->putScalar(zero);
287  deltaX0->elementWiseMultiply(one, *D_, *R0, zero); // deltaX0 = D * R0 (equation 8.13)
288  A10_->apply(*deltaX0, *Rtmp); // Rtmp = A10*D*deltaX0 (intermediate step)
289  Rtmp->update(one, *R1, -one); // Rtmp = R1 - A10*D*deltaX0
290 
291  if (!pL.get<bool>("q2q1 mode")) {
292  deltaX1->putScalar(zero);
293  } else {
294  // special code for q2q1
295  if (Teuchos::rcp_dynamic_cast<BlockedVector>(diagSVector) == Teuchos::null) {
296  ArrayRCP<SC> Sdiag = diagSVector->getDataNonConst(0);
297  ArrayRCP<SC> deltaX1data = deltaX1->getDataNonConst(0);
298  ArrayRCP<SC> Rtmpdata = Rtmp->getDataNonConst(0);
299  for (GO row = 0; row < deltaX1data.size(); row++)
300  deltaX1data[row] = Teuchos::as<SC>(1.1) * Rtmpdata[row] / Sdiag[row];
301  } else {
302  TEUCHOS_TEST_FOR_EXCEPTION(true, MueLu::Exceptions::RuntimeError, "MueLu::BraessSarazinSmoother: q2q1 mode only supported for non-blocked operators.")
303  }
304  }
305 
306  smoo_->Apply(*deltaX1, *Rtmp);
307 
308  // Compute deltaX0
309  deltaX0->putScalar(zero); // just for safety
310  A01_->apply(*deltaX1, *deltaX0); // deltaX0 = A01*deltaX1
311  deltaX0->update(one, *R0, -one); // deltaX0 = R0 - A01*deltaX1
312  R0.swap(deltaX0);
313  deltaX0->elementWiseMultiply(one, *D_, *R0, zero); // deltaX0 = D*(R0 - A01*deltaX1)
314 
315  RCP<MultiVector> X0 = domainMapExtractor_->ExtractVector(rcpX, 0, bDomainThyraMode);
316  RCP<MultiVector> X1 = domainMapExtractor_->ExtractVector(rcpX, 1, bDomainThyraMode);
317 
318  // Update solution
319  X0->update(one, *deltaX0, one);
320  X1->update(one, *deltaX1, one);
321 
322  domainMapExtractor_->InsertVector(X0, 0, rcpX, bDomainThyraMode);
323  domainMapExtractor_->InsertVector(X1, 1, rcpX, bDomainThyraMode);
324 
325  if (run < nSweeps - 1) {
326  R = Utilities::Residual(*A_, *rcpX, *rcpB);
327  }
328  }
329 
330  if (bCopyResultX == true) {
331  RCP<MultiVector> Xmerged = bX->Merge();
332  X.update(one, *Xmerged, zero);
333  }
334 }
335 
336 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
339  return rcp(new BraessSarazinSmoother(*this));
340 }
341 
342 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
344  description() const {
345  std::ostringstream out;
347  out << "{type = " << type_ << "}";
348  return out.str();
349 }
350 
351 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354 
355  if (verbLevel & Parameters0) {
356  out0 << "Prec. type: " << type_ << /*" Sweeps: " << nSweeps_ << " damping: " << omega_ <<*/ std::endl;
357  }
358 
359  if (verbLevel & Debug) {
360  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
361  }
362 }
363 
364 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
366  // FIXME: This is a placeholder
368 }
369 
370 } // namespace MueLu
371 
372 #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)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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
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:99
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