MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_BlockedJacobiSmoother_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_BLOCKEDJACOBISMOOTHER_DEF_HPP_
11 #define MUELU_BLOCKEDJACOBISMOOTHER_DEF_HPP_
12 
13 #include "Teuchos_ArrayViewDecl.hpp"
14 #include "Teuchos_ScalarTraits.hpp"
15 
16 #include "MueLu_ConfigDefs.hpp"
17 
19 #include <Xpetra_Matrix.hpp>
23 #include <Xpetra_MultiVectorFactory.hpp>
24 
26 #include "MueLu_Level.hpp"
27 #include "MueLu_Monitor.hpp"
28 #include "MueLu_HierarchyUtils.hpp"
29 #include "MueLu_SmootherBase.hpp"
30 
31 namespace MueLu {
32 
33 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35  : type_("blocked Jacobi")
36  , A_(Teuchos::null) {
37  FactManager_.reserve(10); // TODO fix me!
38 }
39 
40 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42 
43 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
45  RCP<ParameterList> validParamList = rcp(new ParameterList());
46 
47  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A");
48  validParamList->set<Scalar>("Damping factor", 1.0, "Damping/Scaling factor in Jacobi");
49  validParamList->set<LocalOrdinal>("Sweeps", 1, "Number of sweeps (default = 1)");
50 
51  return validParamList;
52 }
53 
54 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56  TEUCHOS_TEST_FOR_EXCEPTION(pos < 0, Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::AddFactoryManager: parameter \'pos\' must not be negative! error.");
57 
58  size_t myPos = Teuchos::as<size_t>(pos);
59 
60  if (myPos < FactManager_.size()) {
61  // replace existing entris in FactManager_ vector
62  FactManager_.at(myPos) = FactManager;
63  } else if (myPos == FactManager_.size()) {
64  // add new Factory manager in the end of the vector
65  FactManager_.push_back(FactManager);
66  } else { // if(myPos > FactManager_.size())
67  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
68  *out << "Warning: cannot add new FactoryManager at proper position " << pos << ". The FactoryManager is just appended to the end. Check this!" << std::endl;
69 
70  // add new Factory manager in the end of the vector
71  FactManager_.push_back(FactManager);
72  }
73 }
74 
75 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
77  // this->Input(currentLevel, "A");
78  // TODO: check me: why is this->Input not freeing properly A in release mode?
79  currentLevel.DeclareInput("A", this->GetFactory("A").get());
80 
81  // loop over all factory managers for the subblocks of blocked operator A
82  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
83  for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
84  SetFactoryManager currentSFM(rcpFromRef(currentLevel), *it);
85 
86  // request "Smoother" for current subblock row.
87  currentLevel.DeclareInput("PreSmoother", (*it)->GetFactory("Smoother").get());
88 
89  // request "A" for current subblock row (only needed for Thyra mode)
90  currentLevel.DeclareInput("A", (*it)->GetFactory("A").get());
91  }
92 }
93 
94 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
96  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
97 
98  FactoryMonitor m(*this, "Setup blocked Jacobi Smoother", currentLevel);
99  if (SmootherPrototype::IsSetup() == true) this->GetOStream(Warnings0) << "MueLu::BlockedJacobiSmoother::Setup(): Setup() has already been called";
100 
101  // extract blocked operator A from current level
102  A_ = Factory::Get<RCP<Matrix> >(currentLevel, "A"); // A needed for extracting map extractors
103  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
104  TEUCHOS_TEST_FOR_EXCEPTION(bA == Teuchos::null, Exceptions::BadCast, "MueLu::BlockedPFactory::Build: input matrix A is not of type BlockedCrsMatrix! error.");
105 
106  // plausibility check
107  TEUCHOS_TEST_FOR_EXCEPTION(bA->Rows() != FactManager_.size(), Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::Setup: number of block rows of A is " << bA->Rows() << " and does not match number of SubFactoryManagers " << FactManager_.size() << ". error.");
108  TEUCHOS_TEST_FOR_EXCEPTION(bA->Cols() != FactManager_.size(), Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::Setup: number of block cols of A is " << bA->Cols() << " and does not match number of SubFactoryManagers " << FactManager_.size() << ". error.");
109 
110  // store map extractors
111  rangeMapExtractor_ = bA->getRangeMapExtractor();
112  domainMapExtractor_ = bA->getDomainMapExtractor();
113 
114  // loop over all factory managers for the subblocks of blocked operator A
115  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
116  for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
117  SetFactoryManager currentSFM(rcpFromRef(currentLevel), *it);
118 
119  // extract Smoother for current block row
120  RCP<const SmootherBase> Smoo = currentLevel.Get<RCP<SmootherBase> >("PreSmoother", (*it)->GetFactory("Smoother").get());
121  Inverse_.push_back(Smoo);
122 
123  // store whether subblock matrix is blocked or not!
124  RCP<Matrix> Aii = currentLevel.Get<RCP<Matrix> >("A", (*it)->GetFactory("A").get());
125  bIsBlockedOperator_.push_back(Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(Aii) != Teuchos::null);
126  }
127 
129 }
130 
131 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
132 void BlockedJacobiSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero) const {
133  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::Apply(): Setup() has not been called");
134 
135 #if 0 // def HAVE_MUELU_DEBUG
136  // TODO simplify this debug check
137  RCP<MultiVector> rcpDebugX = Teuchos::rcpFromRef(X);
138  RCP<const MultiVector> rcpDebugB = Teuchos::rcpFromRef(B);
139  RCP<BlockedMultiVector> rcpBDebugX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpDebugX);
140  RCP<const BlockedMultiVector> rcpBDebugB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpDebugB);
141  //RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
142  if(rcpBDebugB.is_null() == false) {
143  //this->GetOStream(Runtime1) << "BlockedJacobi: B is a BlockedMultiVector of size " << B.getMap()->getGlobalNumElements() << " with " << rcpBDebugB->getBlockedMap()->getNumMaps() << " blocks." << std::endl;
144  //TEUCHOS_TEST_FOR_EXCEPTION(A_->getRangeMap()->isSameAs(*(B.getMap())) == false, Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::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.");
145  } else {
146  //this->GetOStream(Runtime1) << "BlockedJacobi: B is a MultiVector of size " << B.getMap()->getGlobalNumElements() << std::endl;
147  //TEUCHOS_TEST_FOR_EXCEPTION(bA->getFullRangeMap()->isSameAs(*(B.getMap())) == false, Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::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.");
148  }
149  if(rcpBDebugX.is_null() == false) {
150  //this->GetOStream(Runtime1) << "BlockedJacobi: X is a BlockedMultiVector of size " << X.getMap()->getGlobalNumElements() << " with " << rcpBDebugX->getBlockedMap()->getNumMaps() << " blocks." << std::endl;
151  //TEUCHOS_TEST_FOR_EXCEPTION(A_->getDomainMap()->isSameAs(*(X.getMap())) == false, Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::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.");
152  } else {
153  //this->GetOStream(Runtime1) << "BlockedJacobi: X is a MultiVector of size " << X.getMap()->getGlobalNumElements() << std::endl;
154  //TEUCHOS_TEST_FOR_EXCEPTION(bA->getFullDomainMap()->isSameAs(*(X.getMap())) == false, Exceptions::RuntimeError, "MueLu::BlockedJacobiSmoother::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.");
155  }
156 
157 #endif
159 
160  // Input variables used for the rest of the algorithm
161  RCP<MultiVector> rcpX = Teuchos::rcpFromRef(X);
162  RCP<const MultiVector> rcpB = Teuchos::rcpFromRef(B);
163 
164  // make sure that both rcpX and rcpB are BlockedMultiVector objects
165  bool bCopyResultX = false;
166  bool bReorderX = false;
167  bool bReorderB = false;
168  RCP<BlockedCrsMatrix> bA = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A_);
169  MUELU_TEST_FOR_EXCEPTION(bA.is_null() == true, Exceptions::RuntimeError, "MueLu::BlockedGaussSeidelSmoother::Apply(): A_ must be a BlockedCrsMatrix");
170  RCP<BlockedMultiVector> bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
171  RCP<const BlockedMultiVector> bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
172 
173  // check the type of operator
176 
177  if (rbA.is_null() == false) {
178  // A is a ReorderedBlockedCrsMatrix and thus uses nested maps, retrieve BlockedCrsMatrix and use plain blocked
179  // map for the construction of the blocked multivectors
180 
181  // check type of X vector
182  if (bX.is_null() == true) {
183  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedDomainMap(), rcpX));
184  rcpX.swap(vectorWithBlockedMap);
185  bCopyResultX = true;
186  bReorderX = true;
187  }
188 
189  // check type of B vector
190  if (bB.is_null() == true) {
191  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(rbA->getBlockedCrsMatrix()->getBlockedRangeMap(), rcpB));
192  rcpB.swap(vectorWithBlockedMap);
193  bReorderB = true;
194  }
195  } else {
196  // A is a BlockedCrsMatrix and uses a plain blocked map
197  if (bX.is_null() == true) {
198  RCP<MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedDomainMap(), rcpX));
199  rcpX.swap(vectorWithBlockedMap);
200  bCopyResultX = true;
201  }
202 
203  if (bB.is_null() == true) {
204  RCP<const MultiVector> vectorWithBlockedMap = Teuchos::rcp(new BlockedMultiVector(bA->getBlockedRangeMap(), rcpB));
205  rcpB.swap(vectorWithBlockedMap);
206  }
207  }
208 
209  // we now can guarantee that X and B are blocked multi vectors
210  bX = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(rcpX);
211  bB = Teuchos::rcp_dynamic_cast<const BlockedMultiVector>(rcpB);
212 
213  // Finally we can do a reordering of the blocked multivectors if A is a ReorderedBlockedCrsMatrix
214  if (rbA.is_null() == false) {
215  Teuchos::RCP<const Xpetra::BlockReorderManager> brm = rbA->getBlockReorderManager();
216 
217  // check type of X vector
218  if (bX->getBlockedMap()->getNumMaps() != bA->getDomainMapExtractor()->NumMaps() || bReorderX) {
219  // X is a blocked multi vector but incompatible to the reordered blocked operator A
220  Teuchos::RCP<MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bX);
221  rcpX.swap(reorderedBlockedVector);
222  }
223 
224  if (bB->getBlockedMap()->getNumMaps() != bA->getRangeMapExtractor()->NumMaps() || bReorderB) {
225  // B is a blocked multi vector but incompatible to the reordered blocked operator A
226  Teuchos::RCP<const MultiVector> reorderedBlockedVector = buildReorderedBlockedMultiVector(brm, bB);
227  rcpB.swap(reorderedBlockedVector);
228  }
229  }
230 
231  // Throughout the rest of the algorithm rcpX and rcpB are used for solution vector and RHS
232 
233  RCP<MultiVector> residual = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
234  RCP<MultiVector> tempres = MultiVectorFactory::Build(rcpB->getMap(), rcpB->getNumVectors());
235 
236  // extract parameters from internal parameter list
238  LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
239  Scalar omega = pL.get<Scalar>("Damping factor");
240 
241  // outer Richardson loop
242  for (LocalOrdinal run = 0; run < nSweeps; ++run) {
243  residual->update(1.0, *rcpB, 0.0); // r = B
244  if (InitialGuessIsZero == false || run > 0) {
245  bA->apply(*rcpX, *residual, Teuchos::NO_TRANS, -1.0, 1.0);
246  }
247  // one sweep of Jacobi: loop over all block rows
248  for (size_t i = 0; i < Inverse_.size(); i++) {
249  // extract corresponding subvectors from X and residual
250  bool bRangeThyraMode = rangeMapExtractor_->getThyraMode();
251  bool bDomainThyraMode = domainMapExtractor_->getThyraMode();
252  Teuchos::RCP<MultiVector> Xi = domainMapExtractor_->ExtractVector(rcpX, i, bDomainThyraMode);
253  Teuchos::RCP<MultiVector> ri = rangeMapExtractor_->ExtractVector(residual, i, bRangeThyraMode);
254  Teuchos::RCP<MultiVector> tXi = domainMapExtractor_->getVector(i, X.getNumVectors(), bDomainThyraMode);
255 
256  // apply solver/smoother
257  Inverse_.at(i)->Apply(*tXi, *ri, false);
258 
259  // update vector
260  if (InitialGuessIsZero && run == 0) {
261  Xi->update(omega, *tXi, 0.0); // X_{i+1} = X_i + omega \Delta X_i
262  } else {
263  Xi->update(omega, *tXi, 1.0); // X_{i+1} = X_i + omega \Delta X_i
264  }
265  }
266  }
267 
268  if (bCopyResultX == true) {
269  RCP<MultiVector> Xmerged = bX->Merge();
270  X.update(one, *Xmerged, zero);
271  }
272 }
273 
274 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
276  return rcp(new BlockedJacobiSmoother(*this));
277 }
278 
279 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281  std::ostringstream out;
283  out << "{type = " << type_ << "}";
284  return out.str();
285 }
286 
287 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
290 
291  // extract parameters from internal parameter list
293  LocalOrdinal nSweeps = pL.get<LocalOrdinal>("Sweeps");
294  Scalar omega = pL.get<Scalar>("Damping factor");
295 
296  if (verbLevel & Parameters0) {
297  out0 << "Prec. type: " << type_ << " Sweeps: " << nSweeps << " damping: " << omega << std::endl;
298  }
299 
300  if (verbLevel & Debug) {
301  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl;
302  }
303 }
304 
305 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
307  // FIXME: This is a placeholder
309 }
310 
311 } // namespace MueLu
312 
313 #endif /* MUELU_BLOCKEDJACOBISMOOTHER_DEF_HPP_ */
Important warning messages (one line)
virtual const Teuchos::ParameterList & GetParameterList() const
Exception indicating invalid cast attempted.
MueLu::DefaultLocalOrdinal LocalOrdinal
RCP< const ParameterList > GetValidParameterList() const
Input.
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)
std::string description() const
Return a simple one-line description of this object.
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 Setup(Level &currentLevel)
Setup routine In the Setup method the Inverse_ vector is filled with the corresponding SmootherBase o...
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
void DeclareInput(Level &currentLevel) const
Input.
bool IsSetup() const
Get the state of a smoother prototype.
std::vector< Teuchos::RCP< const FactoryManagerBase > > FactManager_
vector of factory managers
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Print class parameters.
block Jacobi method for blocked matrices
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the direct solver.
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)
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
RCP< SmootherPrototype > Copy() const
virtual std::string description() const
Return a simple one-line description of this object.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager, int pos)
Add a factory manager.
std::string toString(const T &t)
bool is_null() const