46 #ifndef MUELU_IFPACK2SMOOTHER_DECL_HPP
47 #define MUELU_IFPACK2SMOOTHER_DECL_HPP
49 #include <Teuchos_ParameterList.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,
94 #undef MUELU_IFPACK2SMOOTHER_SHORT
139 template<
class Scalar2,
class LocalOrdinal2,
class GlobalOrdinal2,
class Node2>
143 Ifpack2Smoother(
const std::string& type,
const Teuchos::ParameterList& paramList = Teuchos::ParameterList(),
const LO& overlap = 0);
179 void Apply(MultiVector& X,
const MultiVector& B,
bool InitialGuessIsZero =
false)
const;
186 RCP<SmootherPrototype>
Copy()
const;
191 template<
typename Node2>
192 RCP<MueLu::Ifpack2Smoother<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >
193 clone(
const RCP<Node2>& node2,
const Teuchos::RCP<
const Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >& A_newnode)
const;
215 void SetPrecParameters(
const Teuchos::ParameterList& list = Teuchos::ParameterList())
const;
233 RCP<Ifpack2::Preconditioner<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
prec_;
240 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 template<
typename Node2>
242 RCP<MueLu::Ifpack2Smoother<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >
244 #ifdef HAVE_XPETRA_TPETRA
245 const ParameterList& paramList = this->GetParameterList();
246 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> Matrix1;
247 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> Matrix2;
248 RCP<Ifpack2Smoother<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > cloneSmoother =
252 RCP<const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > crsOp =
253 rcp_dynamic_cast<
const Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node2> >(A_newnode);
254 const RCP<const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node2> >& tmp =
255 rcp_dynamic_cast<
const Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node2> >(crsOp->getCrsMatrix());
258 cloneSmoother->prec_ = factory.
clone<Matrix1, Matrix2>(prec_, tmp->getTpetra_CrsMatrix(), paramList);
259 cloneSmoother->type_ = type_;
260 cloneSmoother->SetParameterList(paramList);
261 cloneSmoother->IsSetup(this->IsSetup());
262 return cloneSmoother;
265 "MueLu::Ifpack2Smoother::clone(): clone only available with Tpetra.");
269 #ifdef HAVE_MUELU_EPETRA
271 # if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
272 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
280 #undef MUELU_AMESOS2SMOOTHER_SHORT
286 template<
class Scalar2,
class LocalOrdinal2,
class GlobalOrdinal2,
class Node2>
299 void Apply(MultiVector& X,
const MultiVector& B,
bool InitialGuessIsZero =
false)
const {}
300 RCP<SmootherPrototype>
Copy()
const {
return Teuchos::null;}
302 template<
typename Node2>
303 RCP<MueLu::Ifpack2Smoother<Scalar,LocalOrdinal,GlobalOrdinal,Node2> >
304 clone(
const RCP<Node2>& node2,
const Teuchos::RCP<const Matrix >& A_newnode)
const {
return Teuchos::null; }
316 #endif // HAVE_MUELU_EPETRA
320 #define MUELU_IFPACK2SMOOTHER_SHORT
321 #endif // HAVE_MUELU_TPETRA && HAVE_MUELU_IFPACK2
322 #endif // MUELU_IFPACK2SMOOTHER_DECL_HPP