MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_TpetraOperator_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_TPETRAOPERATOR_DECL_HPP
11 #define MUELU_TPETRAOPERATOR_DECL_HPP
12 
13 #include "MueLu_ConfigDefs.hpp"
14 
15 #include <Tpetra_Operator.hpp>
17 #include "MueLu_Level.hpp"
18 #include "MueLu_Hierarchy_decl.hpp"
19 
20 namespace MueLu {
21 
24 template <class Scalar = Tpetra::Operator<>::scalar_type,
28 class TpetraOperator : public Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
29  protected:
30  TpetraOperator() = delete;
31 
32  public:
34 
35 
38  : Operator_(Op) {}
39 
42  : Hierarchy_(H) {}
43 
45  virtual ~TpetraOperator();
46 
48 
51 
54 
56 
65 
67  bool hasTransposeApply() const;
68 
70 
71 
74 
77 
79 
80  private:
83 };
84 
85 } // namespace MueLu
86 
87 #endif // MUELU_TPETRAOPERATOR_DECL_HPP
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
MueLu::DefaultLocalOrdinal LocalOrdinal
virtual ~TpetraOperator()
Destructor.
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.
MueLu::DefaultNode Node
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
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.
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. ...