46 #ifndef XPETRA_EPETRAOPERATOR_HPP
47 #define XPETRA_EPETRAOPERATOR_HPP
51 #include <Epetra_Operator.h>
52 #include <Epetra_Map.h>
54 #include "Xpetra_Map.hpp"
56 #include "Xpetra_MultiVector.hpp"
64 template<
class EpetraGlobalOrdinal,
class Node>
76 virtual Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
getDomainMap()
const
79 return toXpetra<GlobalOrdinal,Node>(
op_->OperatorDomainMap());
84 virtual Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
getRangeMap()
const
87 return toXpetra<GlobalOrdinal,Node>(
op_->OperatorRangeMap());
100 Teuchos::ETransp mode = Teuchos::NO_TRANS,
101 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
102 Scalar beta = Teuchos::ScalarTraits<Scalar>::zero())
const
110 "Xpetra::EpetraOperator->apply(): can only accept mode == NO_TRANS or mode == TRANS");
112 "Xpetra::EpetraOperator->apply(): cannot apply transpose as underlying Epetra operator does not support it");
115 RCP<Epetra_MultiVector> epY = eY.getEpetra_MultiVector();
116 RCP<Epetra_MultiVector> tmp = Teuchos::rcp(
new Epetra_MultiVector(*epY));
119 op_->SetUseTranspose(mode == Teuchos::TRANS);
131 int err =
op_->SetUseTranspose(
true);
132 op_->SetUseTranspose(
false);
145 void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default)
const
148 out <<
"Epetra_Operator" << std::endl;
163 using STS = Teuchos::ScalarTraits<Scalar>;
164 R.
update(STS::one(),B,STS::zero());
165 this->
apply (X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
172 RCP< Epetra_Operator>
op_;
179 #endif // XPETRA_EPETRAOPERATOR_HPP
std::string description() const
A simple one-line description of this object.
virtual Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
The Map associated with the range of this operator, which must be compatible with Y...
Exception throws to report errors in the internal logical of the program.
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
#define XPETRA_ERR_CHECK(arg)
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
EpetraOperator(const Teuchos::RCP< Epetra_Operator > &op)
EpetraOperator constructor to wrap a Epetra_Operator object.
Exception throws when you call an unimplemented method of Xpetra.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
virtual void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Computes the operator-multivector application.
virtual bool hasTransposeApply() const
Whether this operator supports applying the transpose or conjugate transpose.
EpetraGlobalOrdinal GlobalOrdinal
RCP< Epetra_Operator > op_
The Tpetra::Operator which this class wraps.
#define XPETRA_MONITOR(funcName)
virtual Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
The Map associated with the domain of this operator, which must be compatible with X...