10 #ifndef MUELU_BELOSSMOOTHER_DEF_HPP
11 #define MUELU_BELOSSMOOTHER_DEF_HPP
15 #if defined(HAVE_MUELU_BELOS)
21 #include <Xpetra_MultiVectorFactory.hpp>
22 #ifdef HAVE_XPETRA_TPETRA
23 #include <Xpetra_TpetraMultiVector.hpp>
28 #include "MueLu_Utilities.hpp"
31 #include <BelosLinearProblem.hpp>
32 #include <BelosSolverFactory.hpp>
36 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
39 bool solverSupported =
false;
40 Belos::SolverFactory<Scalar, tMV, tOP> tFactory;
41 solverSupported = solverSupported || tFactory.isSupported(type);
62 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
67 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 this->Input(currentLevel,
"A");
72 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 A_ = Factory::Get<RCP<Matrix> >(currentLevel,
"A");
77 SetupBelos(currentLevel);
79 this->GetOStream(
Statistics1) << description() << std::endl;
82 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
87 tBelosProblem_ =
rcp(
new Belos::LinearProblem<Scalar, tMV, tOP>());
89 tBelosProblem_->setOperator(tA);
91 Belos::SolverFactory<SC, tMV, tOP> solverFactory;
92 tSolver_ = solverFactory.create(type_, rcpFromRef(const_cast<ParameterList&>(this->GetParameterList())));
93 tSolver_->setProblem(tBelosProblem_);
99 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
104 if (InitialGuessIsZero) {
110 tBelosProblem_->setInitResVec(tpB);
111 tBelosProblem_->setProblem(tpX, tpB);
117 RCP<MultiVector> Correction = MultiVectorFactory::Build(A_->getDomainMap(), X.getNumVectors());
122 tBelosProblem_->setInitResVec(tpB);
123 tBelosProblem_->setProblem(tpX, tpB);
126 X.update(TST::one(), *Correction, TST::one());
133 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
142 std::ostringstream out;
145 out << tSolver_->description();
148 out <<
"BELOS {type = " << type_ <<
"}";
153 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
158 out0 <<
"Parameter list: " << std::endl;
160 out << this->GetParameterList();
164 if (tSolver_ != Teuchos::null) {
166 out << *tSolver_ << std::endl;
170 if (verbLevel &
Debug) {
174 <<
"RCP<solver_>: " << tSolver_ << std::endl;
179 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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 additional debugging information.
void Setup(Level ¤tLevel)
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 ¶mList)
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.
void SetupBelos(Level ¤tLevel)
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 ¶mList)
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 ¤tLevel) const
Input.
std::string toString(const T &t)
std::string description() const
Return a simple one-line description of this object.