46 #ifndef MUELU_XPETRAOPERATOR_DECL_HPP 47 #define MUELU_XPETRAOPERATOR_DECL_HPP 51 #include <Xpetra_Operator.hpp> 52 #include <Xpetra_MultiVector.hpp> 60 template <class Scalar = Xpetra::Operator<>::scalar_type,
61 class LocalOrdinal =
typename Xpetra::Operator<Scalar>::local_ordinal_type,
62 class GlobalOrdinal =
typename Xpetra::Operator<Scalar, LocalOrdinal>::global_ordinal_type,
63 class Node =
typename Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
64 class XpetraOperator :
public Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
81 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
getDomainMap()
const {
82 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Matrix;
84 RCP<Matrix> A =
Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
85 return A->getDomainMap();
89 Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
getRangeMap()
const {
90 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Matrix;
92 RCP<Matrix> A =
Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
93 return A->getRangeMap();
101 void apply(
const Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& X,
102 Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Y,
103 Teuchos::ETransp mode = Teuchos::NO_TRANS,
104 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
105 Scalar beta = Teuchos::ScalarTraits<Scalar>::one())
const{
107 #ifdef HAVE_MUELU_DEBUG 108 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Matrix;
109 RCP<Matrix> A =
Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
112 RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Xop =
113 Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getRangeMap(),X.getNumVectors());
114 RCP<Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > Yop =
115 Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(A->getDomainMap(),Y.getNumVectors());
116 TEUCHOS_TEST_FOR_EXCEPTION(A->getRangeMap()->isSameAs(*(Xop->getMap())) ==
false, std::logic_error,
117 "MueLu::XpetraOperator::apply: map of X is incompatible with range map of A");
118 TEUCHOS_TEST_FOR_EXCEPTION(A->getDomainMap()->isSameAs(*(Yop->getMap())) ==
false, std::logic_error,
119 "MueLu::XpetraOperator::apply: map of Y is incompatible with domain map of A");
122 Y.putScalar(Teuchos::ScalarTraits<Scalar>::zero());
124 }
catch (std::exception& e) {
126 std::cerr <<
"Caught an exception in MueLu::XpetraOperator::apply():" << std::endl
127 << e.what() << std::endl;
134 template <
class NewNode>
135 Teuchos::RCP< XpetraOperator<Scalar, LocalOrdinal, GlobalOrdinal, NewNode> >
136 clone(
const RCP<NewNode>& new_node)
const {
149 RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
Hierarchy_;
154 #endif // MUELU_XPETRAOPERATOR_DECL_HPP
Namespace for MueLu classes and methods.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetHierarchy() const
Direct access to the underlying MueLu::Hierarchy.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
XpetraOperator(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H)
Constructor.
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Teuchos::RCP< XpetraOperator< Scalar, LocalOrdinal, GlobalOrdinal, NewNode > > clone(const RCP< NewNode > &new_node) const
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::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 Xpetra::Operator applied to a Xpetra::MultiVector X. ...
virtual ~XpetraOperator()
Destructor.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Wraps an existing MueLu::Hierarchy as a Xpetra::Operator.