47 #ifndef THYRA_XPETRA_LINEAR_OP_DECL_HPP
48 #define THYRA_XPETRA_LINEAR_OP_DECL_HPP
50 #include "Thyra_LinearOpDefaultBase.hpp"
51 #include "Xpetra_Operator.hpp"
52 #include "Teuchos_ConstNonconstObjectContainer.hpp"
65 class Node=KokkosClassic::DefaultNode::DefaultNodeType>
67 :
virtual public Thyra::LinearOpDefaultBase<Scalar>
79 const RCP<
const VectorSpaceBase<Scalar> > &rangeSpace,
80 const RCP<
const VectorSpaceBase<Scalar> > &domainSpace,
81 const RCP<Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &xpetraOperator
86 const RCP<
const VectorSpaceBase<Scalar> > &rangeSpace,
87 const RCP<
const VectorSpaceBase<Scalar> > &domainSpace,
88 const RCP<
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &xpetraOperator
122 const Thyra::EOpTransp M_trans,
123 const Thyra::MultiVectorBase<Scalar> &X_in,
124 const Teuchos::Ptr<Thyra::MultiVectorBase<Scalar> > &Y_inout,
142 template<
class XpetraOperator_t>
144 const RCP<
const VectorSpaceBase<Scalar> > &rangeSpace,
145 const RCP<
const VectorSpaceBase<Scalar> > &domainSpace,
156 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
159 const RCP<
const VectorSpaceBase<Scalar> > &rangeSpace,
160 const RCP<
const VectorSpaceBase<Scalar> > &domainSpace,
161 const RCP<Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &xpetraOperator
166 op->initialize(rangeSpace, domainSpace, xpetraOperator);
175 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
178 const RCP<
const VectorSpaceBase<Scalar> > &rangeSpace,
179 const RCP<
const VectorSpaceBase<Scalar> > &domainSpace,
180 const RCP<
const Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &xpetraOperator
185 op->constInitialize(rangeSpace, domainSpace, xpetraOperator);
191 #endif // THYRA_XPETRA_LINEAR_OP_DECL_HPP
RCP< const VectorSpaceBase< Scalar > > domainSpace_
MueLu::DefaultLocalOrdinal LocalOrdinal
RCP< const Thyra::VectorSpaceBase< Scalar > > range() const
RCP< const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstXpetraOperator() const
Get embedded const Xpetra::Operator.
RCP< const XpetraLinearOp< Scalar, LocalOrdinal, GlobalOrdinal, Node > > constXpetraLinearOp(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &xpetraOperator)
Nonmmeber constructor for XpetraLinearOp.
bool opSupportedImpl(Thyra::EOpTransp M_trans) const
RCP< const Thyra::VectorSpaceBase< Scalar > > domain() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void initializeImpl(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< XpetraOperator_t > &xpetraOperator)
XpetraLinearOp()
Construct to uninitialized.
void applyImpl(const Thyra::EOpTransp M_trans, const Thyra::MultiVectorBase< Scalar > &X_in, const Teuchos::Ptr< Thyra::MultiVectorBase< Scalar > > &Y_inout, const Scalar alpha, const Scalar beta) const
RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getXpetraOperator()
Get embedded non-const Xpetra::Operator.
Concrete Thyra::LinearOpBase subclass for Xpetra::Operator.
void initialize(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &xpetraOperator)
Initialize.
RCP< const VectorSpaceBase< Scalar > > rangeSpace_
RCP< XpetraLinearOp< Scalar, LocalOrdinal, GlobalOrdinal, Node > > xpetraLinearOp(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &xpetraOperator)
Nonmmeber constructor for XpetraLinearOp.
Teuchos::ConstNonconstObjectContainer< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > xpetraOperator_
void constInitialize(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &xpetraOperator)
Initialize.