MueLu Version of the Day
Loading...
Searching...
No Matches
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
17#include <Teuchos_ParameterList.hpp>
18
19#include <Xpetra_CrsMatrix.hpp>
20#include <Xpetra_Matrix.hpp>
21#include <Xpetra_MultiVectorFactory.hpp>
22#include <Xpetra_TpetraMultiVector.hpp>
23
25#include "MueLu_Level.hpp"
26#include "MueLu_Utilities.hpp"
27#include "MueLu_Monitor.hpp"
28
29#include <BelosLinearProblem.hpp>
30#include <BelosSolverFactory.hpp>
31
32namespace MueLu {
33
34template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35BelosSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BelosSmoother(const std::string type, const Teuchos::ParameterList& paramList)
36 : type_(type) {
37 bool solverSupported = false;
38 Belos::SolverFactory<Scalar, tMV, tOP> tFactory;
39 solverSupported = solverSupported || tFactory.isSupported(type);
40 // if (!solverSupported) {
41 // Teuchos::Array<std::string> supportedSolverNames = factory.supportedSolverNames();
42
43 // std::ostringstream outString;
44 // outString << "[";
45 // for (auto iter = supportedSolverNames.begin(); iter != supportedSolverNames.end(); ++iter) {
46 // outString << "\"" << *iter << "\"";
47 // if (iter + 1 != supportedSolverNames.end()) {
48 // outString << ", ";
49 // }
50 // }
51 // outString << "]";
52
53 // TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Belos does not provide the solver '" << type_ << "'.\nSupported Solvers: " << outString.str());
54 // }
55 this->declareConstructionOutcome(!solverSupported, "Belos does not provide the smoother '" + type_ + "'.");
56 if (solverSupported)
57 SetParameterList(paramList);
58}
59
60template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63}
64
65template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67 this->Input(currentLevel, "A");
68}
69
70template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72 FactoryMonitor m(*this, "Setup Smoother", currentLevel);
73
74 A_ = Factory::Get<RCP<Matrix> >(currentLevel, "A");
75 SetupBelos(currentLevel);
77 this->GetOStream(Statistics1) << description() << std::endl;
78}
79
80template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82 bool useTpetra = A_->getRowMap()->lib() == Xpetra::UseTpetra;
83
84 if (useTpetra) {
85 tBelosProblem_ = rcp(new Belos::LinearProblem<Scalar, tMV, tOP>());
86 RCP<const tOP> tA = toTpetra(A_);
87 tBelosProblem_->setOperator(tA);
88
89 Belos::SolverFactory<SC, tMV, tOP> solverFactory;
90 tSolver_ = solverFactory.create(type_, rcpFromRef(const_cast<ParameterList&>(this->GetParameterList())));
91 tSolver_->setProblem(tBelosProblem_);
92 } else {
93 TEUCHOS_ASSERT(false);
94 }
95}
96
97template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98void BelosSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero) const {
99 TEUCHOS_TEST_FOR_EXCEPTION(SmootherPrototype::IsSetup() == false, Exceptions::RuntimeError, "MueLu::BelosSmoother::Apply(): Setup() has not been called");
100
101 if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
102 if (InitialGuessIsZero) {
103 X.putScalar(0.0);
104
105 RCP<Tpetra::MultiVector<SC, LO, GO, NO> > tpX = toTpetra(rcpFromRef(X));
106 RCP<const Tpetra::MultiVector<SC, LO, GO, NO> > tpB = toTpetra(rcpFromRef(B));
107
108 tBelosProblem_->setInitResVec(tpB);
109 tBelosProblem_->setProblem(tpX, tpB);
110 tSolver_->solve();
111
112 } else {
113 typedef Teuchos::ScalarTraits<Scalar> TST;
114 RCP<MultiVector> Residual = Utilities::Residual(*A_, X, B);
115 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
116
117 RCP<Tpetra::MultiVector<SC, LO, GO, NO> > tpX = toTpetra(Correction);
118 RCP<const Tpetra::MultiVector<SC, LO, GO, NO> > tpB = toTpetra(Residual);
119
120 tBelosProblem_->setInitResVec(tpB);
121 tBelosProblem_->setProblem(tpX, tpB);
122 tSolver_->solve();
123
124 X.update(TST::one(), *Correction, TST::one());
125 }
126 } else {
127 TEUCHOS_ASSERT(false);
128 }
129}
130
131template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
132RCP<MueLu::SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> > BelosSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Copy() const {
133 RCP<BelosSmoother> smoother = rcp(new BelosSmoother(*this));
134 smoother->SetParameterList(this->GetParameterList());
135 return smoother;
136}
137
138template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
140 std::ostringstream out;
142 if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
143 out << tSolver_->description();
144 }
145 } else {
146 out << "BELOS {type = " << type_ << "}";
147 }
148 return out.str();
149}
150
151template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152void BelosSmoother<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(Teuchos::FancyOStream& out, const VerbLevel verbLevel) const {
154
155 if (verbLevel & Parameters1) {
156 out0 << "Parameter list: " << std::endl;
157 Teuchos::OSTab tab2(out);
158 out << this->GetParameterList();
159 }
160
161 if (verbLevel & External) {
162 if (tSolver_ != Teuchos::null) {
163 Teuchos::OSTab tab2(out);
164 out << *tSolver_ << std::endl;
165 }
166 }
167
168 if (verbLevel & Debug) {
169 if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
170 out0 << "IsSetup: " << Teuchos::toString(SmootherPrototype::IsSetup()) << std::endl
171 << "-" << std::endl
172 << "RCP<solver_>: " << tSolver_ << std::endl;
173 }
174 }
175}
176
177template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
179 return Teuchos::OrdinalTraits<size_t>::invalid();
180}
181
182} // namespace MueLu
183
184#endif // HAVE_MUELU_BELOS
185#endif // MUELU_BELOSSMOOTHER_DEF_HPP
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Class that encapsulates Belos smoothers.
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
friend class BelosSmoother
Constructor.
void SetParameterList(const Teuchos::ParameterList &paramList)
Set parameters from a parameter list and return with default values.
void DeclareInput(Level &currentLevel) const
Input.
void SetupBelos(Level &currentLevel)
std::string description() const
Return a simple one-line description of this object.
size_t getNodeSmootherComplexity() const
For diagnostic purposes.
void Setup(Level &currentLevel)
Set up the smoother.
RCP< SmootherPrototype > Copy() const
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
virtual void SetParameterList(const Teuchos::ParameterList &paramList)=0
Set parameters from a parameter list and return with default values.
void declareConstructionOutcome(bool fail, std::string msg)
bool IsSetup() const
Get the state of a smoother prototype.
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
Namespace for MueLu classes and methods.
@ Debug
Print additional debugging information.
@ Statistics1
Print more statistics.
@ External
Print external lib objects.
@ Parameters1
Print class parameters (more parameters, more verbose)