10#ifndef MUELU_BELOSSMOOTHER_DEF_HPP
11#define MUELU_BELOSSMOOTHER_DEF_HPP
15#if defined(HAVE_MUELU_BELOS)
17#include <Teuchos_ParameterList.hpp>
19#include <Xpetra_CrsMatrix.hpp>
20#include <Xpetra_Matrix.hpp>
21#include <Xpetra_MultiVectorFactory.hpp>
22#include <Xpetra_TpetraMultiVector.hpp>
26#include "MueLu_Utilities.hpp"
29#include <BelosLinearProblem.hpp>
30#include <BelosSolverFactory.hpp>
34template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
37 bool solverSupported =
false;
38 Belos::SolverFactory<Scalar, tMV, tOP> tFactory;
39 solverSupported = solverSupported || tFactory.isSupported(type);
60template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
65template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
67 this->Input(currentLevel,
"A");
70template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
74 A_ = Factory::Get<RCP<Matrix> >(currentLevel,
"A");
75 SetupBelos(currentLevel);
77 this->GetOStream(
Statistics1) << description() << std::endl;
80template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
82 bool useTpetra = A_->getRowMap()->lib() == Xpetra::UseTpetra;
85 tBelosProblem_ = rcp(
new Belos::LinearProblem<Scalar, tMV, tOP>());
86 RCP<const tOP> tA = toTpetra(A_);
87 tBelosProblem_->setOperator(tA);
89 Belos::SolverFactory<SC, tMV, tOP> solverFactory;
90 tSolver_ = solverFactory.create(type_, rcpFromRef(
const_cast<ParameterList&
>(this->GetParameterList())));
91 tSolver_->setProblem(tBelosProblem_);
93 TEUCHOS_ASSERT(
false);
97template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
101 if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
102 if (InitialGuessIsZero) {
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));
108 tBelosProblem_->setInitResVec(tpB);
109 tBelosProblem_->setProblem(tpX, tpB);
113 typedef Teuchos::ScalarTraits<Scalar> TST;
115 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
117 RCP<Tpetra::MultiVector<SC, LO, GO, NO> > tpX = toTpetra(Correction);
118 RCP<const Tpetra::MultiVector<SC, LO, GO, NO> > tpB = toTpetra(Residual);
120 tBelosProblem_->setInitResVec(tpB);
121 tBelosProblem_->setProblem(tpX, tpB);
124 X.update(TST::one(), *Correction, TST::one());
127 TEUCHOS_ASSERT(
false);
131template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 smoother->SetParameterList(this->GetParameterList());
138template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 std::ostringstream out;
142 if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
143 out << tSolver_->description();
146 out <<
"BELOS {type = " << type_ <<
"}";
151template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
156 out0 <<
"Parameter list: " << std::endl;
157 Teuchos::OSTab tab2(out);
158 out << this->GetParameterList();
162 if (tSolver_ != Teuchos::null) {
163 Teuchos::OSTab tab2(out);
164 out << *tSolver_ << std::endl;
168 if (verbLevel &
Debug) {
169 if (A_->getRowMap()->lib() == Xpetra::UseTpetra) {
172 <<
"RCP<solver_>: " << tSolver_ << std::endl;
177template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
179 return Teuchos::OrdinalTraits<size_t>::invalid();
#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 ¶mList)
Set parameters from a parameter list and return with default values.
void DeclareInput(Level ¤tLevel) const
Input.
void SetupBelos(Level ¤tLevel)
std::string description() const
Return a simple one-line description of this object.
size_t getNodeSmootherComplexity() const
For diagnostic purposes.
void Setup(Level ¤tLevel)
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 ¶mList)=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)