10#ifndef MUELU_XPETRAOPERATOR_DEF_HPP
11#define MUELU_XPETRAOPERATOR_DEF_HPP
17template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
22template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
26template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
27const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> >
30 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Matrix;
32 RCP<Matrix> A = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
33 return A->getDomainMap();
36template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
37const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> >
40 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Matrix;
42 RCP<Matrix> A = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
43 return A->getRangeMap();
46template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
48 apply(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
49 Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Y,
50 Teuchos::ETransp mode,
53 TEUCHOS_TEST_FOR_EXCEPTION(mode != Teuchos::NO_TRANS, std::logic_error,
"MueLu::XpetraOperator does not support applying the adjoint operator");
55#ifdef HAVE_MUELU_DEBUG
56 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Matrix;
57 RCP<Matrix> A = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
60 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Xop =
61 Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A->getRangeMap(), X.getNumVectors());
62 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Yop =
63 Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(A->getDomainMap(), Y.getNumVectors());
64 TEUCHOS_TEST_FOR_EXCEPTION(A->getRangeMap()->isSameAs(*(Xop->getMap())) ==
false, std::logic_error,
65 "MueLu::XpetraOperator::apply: map of X is incompatible with range map of A");
66 TEUCHOS_TEST_FOR_EXCEPTION(A->getDomainMap()->isSameAs(*(Yop->getMap())) ==
false, std::logic_error,
67 "MueLu::XpetraOperator::apply: map of Y is incompatible with domain map of A");
69 Hierarchy_->Iterate(X, Y, 1,
true);
70 }
catch (std::exception& e) {
72 std::cerr <<
"Caught an exception in MueLu::XpetraOperator::apply():" << std::endl
73 << e.what() << std::endl;
77template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
83template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
85 residual(
const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
86 const Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& B,
87 Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& R)
const {
88 using STS = Teuchos::ScalarTraits<Scalar>;
89 R.update(STS::one(), B, STS::zero());
90 this->apply(X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
93template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
94RCP<MueLu::Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
MueLu::DefaultScalar Scalar
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
virtual ~XpetraOperator()
Destructor.
void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar=Teuchos::ScalarTraits< Scalar >::one(), Scalar=Teuchos::ScalarTraits< Scalar >::one()) const
Returns in Y the result of a Xpetra::Operator applied to a Xpetra::MultiVector X.
void residual(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetHierarchy() const
Direct access to the underlying MueLu::Hierarchy.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Namespace for MueLu classes and methods.