46 #ifndef MUELU_IFPACK2SMOOTHER_DECL_HPP
47 #define MUELU_IFPACK2SMOOTHER_DECL_HPP
51 #include <Xpetra_BlockedCrsMatrix_fwd.hpp>
52 #include <Xpetra_CrsMatrixWrap.hpp>
53 #include <Xpetra_Matrix_fwd.hpp>
54 #include <Xpetra_Matrix.hpp>
55 #include <Xpetra_MultiVectorFactory_fwd.hpp>
56 #ifdef HAVE_XPETRA_TPETRA // needed for clone()
57 #include <Xpetra_TpetraCrsMatrix.hpp>
63 #if defined(HAVE_MUELU_TPETRA) && defined(HAVE_MUELU_IFPACK2)
65 #include <Tpetra_CrsMatrix.hpp>
67 #include <Ifpack2_Factory_decl.hpp>
68 #include <Ifpack2_Factory_def.hpp>
74 #include "MueLu_SmootherPrototype.hpp"
89 template <class Scalar = SmootherPrototype<>::scalar_type,
90 class LocalOrdinal =
typename SmootherPrototype<Scalar>::local_ordinal_type,
91 class GlobalOrdinal =
typename SmootherPrototype<Scalar, LocalOrdinal>::global_ordinal_type,
92 class Node =
typename SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
94 #undef MUELU_IFPACK2SMOOTHER_SHORT
139 template<
class Scalar2,
class LocalOrdinal2,
class GlobalOrdinal2,
class Node2>
179 void Apply(MultiVector& X,
const MultiVector& B,
bool InitialGuessIsZero =
false)
const;
238 #ifdef HAVE_MUELU_EPETRA
240 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
241 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
249 #undef MUELU_AMESOS2SMOOTHER_SHORT
255 template<
class Scalar2,
class LocalOrdinal2,
class GlobalOrdinal2,
class Node2>
268 void Apply(MultiVector& X,
const MultiVector& B,
bool InitialGuessIsZero =
false)
const {}
282 #endif // HAVE_MUELU_EPETRA
286 #define MUELU_IFPACK2SMOOTHER_SHORT
287 #endif // HAVE_MUELU_TPETRA && HAVE_MUELU_IFPACK2
288 #endif // MUELU_IFPACK2SMOOTHER_DECL_HPP
Ifpack2Smoother(const std::string &type, const Teuchos::ParameterList ¶mList=Teuchos::ParameterList(), const LocalOrdinal &overlap=0)
void DeclareInput(Level ¤tLevel) const
Input.
void Setup(Level ¤tLevel)
void SetupGeneric(Level ¤tLevel)
RCP< Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getPreconditioner()
For diagnostic purposes.
MueLu::DefaultLocalOrdinal LocalOrdinal
friend class Ifpack2Smoother
Constructor.
RCP< Matrix > A_
matrix, used in apply if solving residual equation
Base class for smoother prototypes.
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
RCP< SmootherPrototype > Copy() const
virtual ~Ifpack2Smoother()
LO overlap_
overlap when using the smoother in additive Schwarz mode
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
std::string type_
ifpack2-specific key phrase that denote smoother type
void SetupBlockRelaxation(Level ¤tLevel)
void SetupSchwarz(Level ¤tLevel)
void SetupLineSmoothing(Level ¤tLevel)
std::string description() const
Return a simple one-line description of this object.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
RCP< SmootherPrototype > Copy() const
Class that holds all level-specific information.
void Setup(Level ¤tLevel)
Set up the smoother.
void DeclareInput(Level ¤tLevel) const
Input.
void SetupTopological(Level ¤tLevel)
#define MUELU_TPETRA_ETI_EXCEPTION(cl, obj, go)
std::string description() const
Return a simple one-line description of this object.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply smoother.
void SetPrecParameters(const Teuchos::ParameterList &list=Teuchos::ParameterList()) const
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
virtual ~Ifpack2Smoother()
Destructor.
bool constructionSuccessful_
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the preconditioner.
void SetupChebyshev(Level ¤tLevel)
Kokkos::Compat::KokkosSerialWrapperNode EpetraNode
Class that encapsulates Ifpack2 smoothers.
RCP< Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node > > prec_
pointer to Ifpack2 preconditioner object