MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_BelosSmoother_decl.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_DECL_HPP
11#define MUELU_BELOSSMOOTHER_DECL_HPP
12
13#include <Teuchos_ParameterList.hpp>
14
15#include <Xpetra_Matrix_fwd.hpp>
16#include <Xpetra_MultiVectorFactory_fwd.hpp>
17
18#include "MueLu_ConfigDefs.hpp"
19
20#if defined(HAVE_MUELU_BELOS)
21
23
24#include "BelosTpetraAdapter.hpp"
25#include <Tpetra_Operator.hpp>
26#include <Tpetra_MultiVector.hpp>
27
28#include "MueLu_Level_fwd.hpp"
29#include "MueLu_SmootherPrototype.hpp"
31
32#include <BelosSolverManager.hpp>
33
34namespace MueLu {
35
46template <class Scalar = SmootherPrototype<>::scalar_type,
50class BelosSmoother : public SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
51#undef MUELU_BELOSSMOOTHER_SHORT
53
54 public:
56
57 // TODO: update doc for Belos.
67#ifndef _MSC_VER
68 // Avoid error C3772: invalid friend template declaration
69 template <class Scalar2, class LocalOrdinal2, class GlobalOrdinal2, class Node2>
70 friend class BelosSmoother;
71#endif
72
73 BelosSmoother(const std::string type, const Teuchos::ParameterList& paramList = Teuchos::ParameterList());
74
76 virtual ~BelosSmoother() = default;
77
79
80 void SetParameterList(const Teuchos::ParameterList& paramList);
81
83
84
85 void DeclareInput(Level& currentLevel) const;
86
88
90
91
97 void Setup(Level& currentLevel);
98
107 void Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero = false) const;
108
110
112
113
114 RCP<SmootherPrototype> Copy() const;
115
117
119
120
122 std::string description() const;
123
125 // using MueLu::Describable::describe; // overloading, not hiding
126 // void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel = Default) const
127 void print(Teuchos::FancyOStream& out, const VerbLevel verbLevel = Default) const;
128
130
132 // RCP<Ifpack2::Preconditioner<Scalar,LocalOrdinal,GlobalOrdinal,Node> > getPreconditioner(){return prec_;}
133
135 size_t getNodeSmootherComplexity() const;
136
137 private:
138 void SetupBelos(Level& currentLevel);
139
140 private:
141 std::string type_;
142
143 typedef Tpetra::MultiVector<SC, LO, GO, NO> tMV;
144 typedef Tpetra::Operator<SC, LO, GO, NO> tOP;
145 RCP<Belos::LinearProblem<Scalar, tMV, tOP> > tBelosProblem_;
146 RCP<Belos::SolverManager<Scalar, tMV, tOP> > tSolver_;
147
149 RCP<Matrix> A_;
150
151}; // class BelosSmoother
152
153} // namespace MueLu
154
155#define MUELU_BELOSSMOOTHER_SHORT
156#endif // HAVE_MUELU_BELOS
157#endif // MUELU_BELOSSMOOTHER_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Class that encapsulates Belos smoothers.
Tpetra::MultiVector< SC, LO, GO, NO > tMV
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
Tpetra::Operator< SC, LO, GO, NO > tOP
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)
virtual ~BelosSmoother()=default
Destructor.
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
Class that holds all level-specific information.
Base class for smoother prototypes.
Namespace for MueLu classes and methods.