10 #ifndef MUELU_XPETRAOPERATOR_DEF_HPP
11 #define MUELU_XPETRAOPERATOR_DEF_HPP
17 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
22 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
26 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
32 RCP<Matrix> A = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
33 return A->getDomainMap();
36 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
42 RCP<Matrix> A = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
43 return A->getRangeMap();
46 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
55 #ifdef HAVE_MUELU_DEBUG
57 RCP<Matrix> A = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >(
"A");
65 "MueLu::XpetraOperator::apply: map of X is incompatible with range map of A");
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;
77 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
83 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
89 R.
update(STS::one(), B, STS::zero());
93 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
101 #endif // MUELU_XPETRAOPERATOR_DEF_HPP
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &map, size_t NumVectors, bool zeroOut=true)
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
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.
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
virtual ~XpetraOperator()
Destructor.
MueLu::DefaultScalar Scalar
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.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
virtual size_t getNumVectors() const =0
Provides methods to build a multigrid hierarchy and apply multigrid cycles.