46#ifndef MUELU_BELOSSMOOTHER_DEF_HPP
47#define MUELU_BELOSSMOOTHER_DEF_HPP
51#if defined(HAVE_MUELU_BELOS)
53#include <Teuchos_ParameterList.hpp>
55#include <Xpetra_CrsMatrix.hpp>
56#include <Xpetra_Matrix.hpp>
57#include <Xpetra_MultiVectorFactory.hpp>
58#ifdef HAVE_XPETRA_TPETRA
59#include <Xpetra_TpetraMultiVector.hpp>
64#include "MueLu_Utilities.hpp"
67#include <BelosLinearProblem.hpp>
68#include <BelosSolverFactory.hpp>
74 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
78 bool solverSupported =
false;
79 Belos::SolverFactory<Scalar,tMV,tOP> tFactory;
80 solverSupported = solverSupported || tFactory.isSupported(type);
101 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
106 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
109 this->
Input(currentLevel,
"A");
113 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
124 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
127 bool useTpetra =
A_->getRowMap()->lib() == Xpetra::UseTpetra;
130 tBelosProblem_ = rcp(
new Belos::LinearProblem<Scalar, tMV, tOP>());
134 Belos::SolverFactory<SC, tMV, tOP> solverFactory;
138 TEUCHOS_ASSERT(
false);
142 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
146 if (
A_->getRowMap()->lib() == Xpetra::UseTpetra) {
147 if (InitialGuessIsZero) {
158 typedef Teuchos::ScalarTraits<Scalar> TST;
160 RCP<MultiVector> Correction = MultiVectorFactory::Build(
A_->getDomainMap(), X.getNumVectors());
169 X.update(TST::one(), *Correction, TST::one());
172 TEUCHOS_ASSERT(
false);
177 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
179 RCP<BelosSmoother> smoother = rcp(
new BelosSmoother(*
this) );
184 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
186 std::ostringstream out;
188 if (
A_->getRowMap()->lib() == Xpetra::UseTpetra) {
192 out <<
"BELOS {type = " <<
type_ <<
"}";
197 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
202 out0 <<
"Parameter list: " << std::endl;
203 Teuchos::OSTab tab2(out);
209 Teuchos::OSTab tab2(out);
214 if (verbLevel &
Debug) {
215 if (
A_->getRowMap()->lib() == Xpetra::UseTpetra) {
218 <<
"RCP<solver_>: " <<
tSolver_ << std::endl;
223 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
225 return Teuchos::OrdinalTraits<size_t>::invalid();
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
RCP< Belos::SolverManager< Scalar, tMV, tOP > > tSolver_
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.
RCP< Belos::LinearProblem< Scalar, tMV, tOP > > tBelosProblem_
RCP< Matrix > A_
matrix, used in apply if solving residual equation
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.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
Class that holds all level-specific information.
virtual const Teuchos::ParameterList & GetParameterList() const
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, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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.
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
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).