46 #ifndef MUELU_TPETRAOPERATOR_DECL_HPP
47 #define MUELU_TPETRAOPERATOR_DECL_HPP
51 #ifdef HAVE_MUELU_TPETRA
52 #include <Tpetra_Operator.hpp>
53 #include <Tpetra_MultiVector_decl.hpp>
61 template <class Scalar = Tpetra::Operator<>::scalar_type,
62 class LocalOrdinal =
typename Tpetra::Operator<Scalar>::local_ordinal_type,
63 class GlobalOrdinal =
typename Tpetra::Operator<Scalar, LocalOrdinal>::global_ordinal_type,
64 class Node =
typename Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
65 class TpetraOperator :
public Tpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
95 void apply(
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
96 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Y,
104 #ifdef HAVE_MUELU_DEPRECATED_CODE
105 template <
class NewNode>
130 #endif //ifdef HAVE_MUELU_TPETRA
132 #endif // MUELU_TPETRAOPERATOR_DECL_HPP
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Operator_
TpetraOperator(const RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Op)
Constructor.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetOperator() const
Direct access to the underlying MueLu::Operator.
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
virtual ~TpetraOperator()
Destructor.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetHierarchy() const
Direct access to the underlying MueLu::Hierarchy.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
TpetraOperator(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H)
Constructor.
void apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::one()) const
Returns in Y the result of a Tpetra::Operator applied to a Tpetra::MultiVector X. ...