MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_BelosSmoother_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_BELOSSMOOTHER_DEF_HPP
11 #define MUELU_BELOSSMOOTHER_DEF_HPP
12 
13 #include "MueLu_ConfigDefs.hpp"
14 
15 #if defined(HAVE_MUELU_BELOS)
16 
18 
19 #include <Xpetra_CrsMatrix.hpp>
20 #include <Xpetra_Matrix.hpp>
21 #include <Xpetra_MultiVectorFactory.hpp>
22 #ifdef HAVE_XPETRA_TPETRA
23 #include <Xpetra_TpetraMultiVector.hpp>
24 #endif
25 
27 #include "MueLu_Level.hpp"
28 #include "MueLu_Utilities.hpp"
29 #include "MueLu_Monitor.hpp"
30 
31 #include <BelosLinearProblem.hpp>
32 #include <BelosSolverFactory.hpp>
33 
34 namespace MueLu {
35 
36 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38  : type_(type) {
39  bool solverSupported = false;
40  Belos::SolverFactory<Scalar, tMV, tOP> tFactory;
41  solverSupported = solverSupported || tFactory.isSupported(type);
42  // if (!solverSupported) {
43  // Teuchos::Array<std::string> supportedSolverNames = factory.supportedSolverNames();
44 
45  // std::ostringstream outString;
46  // outString << "[";
47  // for (auto iter = supportedSolverNames.begin(); iter != supportedSolverNames.end(); ++iter) {
48  // outString << "\"" << *iter << "\"";
49  // if (iter + 1 != supportedSolverNames.end()) {
50  // outString << ", ";
51  // }
52  // }
53  // outString << "]";
54 
55  // TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Belos does not provide the solver '" << type_ << "'.\nSupported Solvers: " << outString.str());
56  // }
57  this->declareConstructionOutcome(!solverSupported, "Belos does not provide the smoother '" + type_ + "'.");
58  if (solverSupported)
59  SetParameterList(paramList);
60 }
61 
62 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64  Factory::SetParameterList(paramList);
65 }
66 
67 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69  this->Input(currentLevel, "A");
70 }
71 
72 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74  FactoryMonitor m(*this, "Setup Smoother", currentLevel);
75 
76  A_ = Factory::Get<RCP<Matrix> >(currentLevel, "A");
77  SetupBelos(currentLevel);
79  this->GetOStream(Statistics1) << description() << std::endl;
80 }
81 
82 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
84  bool useTpetra = A_->getRowMap()->lib() == Xpetra::UseTpetra;
85 
86  if (useTpetra) {
87  tBelosProblem_ = rcp(new Belos::LinearProblem<Scalar, tMV, tOP>());
89  tBelosProblem_->setOperator(tA);
90 
91  Belos::SolverFactory<SC, tMV, tOP> solverFactory;
92  tSolver_ = solverFactory.create(type_, rcpFromRef(const_cast<ParameterList&>(this->GetParameterList())));
93  tSolver_->setProblem(tBelosProblem_);
94  } else {
95  TEUCHOS_ASSERT(false);
96  }
97 }
98 
99 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
100 void BelosSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
101  TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::BelosSmoother::Apply(): Setup() has not been called");
102 
103  if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
104  if (InitialGuessIsZero) {
105  X.putScalar(0.0);
106 
109 
110  tBelosProblem_->setInitResVec(tpB);
111  tBelosProblem_->setProblem(tpX, tpB);
112  tSolver_->solve();
113 
114  } else {
115  typedef Teuchos::ScalarTraits<Scalar> TST;
116  RCP<MultiVector> Residual = Utilities::Residual(*A_, X, B);
117  RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
118 
121 
122  tBelosProblem_->setInitResVec(tpB);
123  tBelosProblem_->setProblem(tpX, tpB);
124  tSolver_->solve();
125 
126  X.update(TST::one(), *Correction, TST::one());
127  }
128  } else {
129  TEUCHOS_ASSERT(false);
130  }
131 }
132 
133 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
135  RCP<BelosSmoother> smoother = rcp(new BelosSmoother(*this));
136  smoother->SetParameterList(this->GetParameterList());
137  return smoother;
138 }
139 
140 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
142  std::ostringstream out;
144  if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
145  out << tSolver_->description();
146  }
147  } else {
148  out << "BELOS {type = " << type_ << "}";
149  }
150  return out.str();
151 }
152 
153 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
156 
157  if (verbLevel & Parameters1) {
158  out0 << "Parameter list: " << std::endl;
159  Teuchos::OSTab tab2(out);
160  out << this->GetParameterList();
161  }
162 
163  if (verbLevel & External) {
164  if (tSolver_ != Teuchos::null) {
165  Teuchos::OSTab tab2(out);
166  out << *tSolver_ << std::endl;
167  }
168  }
169 
170  if (verbLevel & Debug) {
171  if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
172  out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
173  << "-" << std::endl
174  << "RCP<solver_>: " << tSolver_ << std::endl;
175  }
176  }
177 }
178 
179 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182 }
183 
184 } // namespace MueLu
185 
186 #endif // HAVE_MUELU_BELOS
187 #endif // MUELU_BELOSSMOOTHER_DEF_HPP
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
Print external lib objects.
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)
Class that encapsulates Belos smoothers.
Print more statistics.
Print additional debugging information.
void Setup(Level &currentLevel)
Set up the smoother.
size_t getNodeSmootherComplexity() const
For diagnostic purposes.
RCP< SmootherPrototype > Copy() const
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
virtual void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void declareConstructionOutcome(bool fail, std::string msg)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
void SetupBelos(Level &currentLevel)
friend class BelosSmoother
Constructor.
bool IsSetup() const
Get the state of a smoother prototype.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
Print class parameters (more parameters, more verbose)
Exception throws to report errors in the internal logical of the program.
#define TEUCHOS_ASSERT(assertion_test)
void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> Op)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> vec)
void DeclareInput(Level &currentLevel) const
Input.
std::string toString(const T &t)
std::string description() const
Return a simple one-line description of this object.