42 #include "Thyra_EpetraOperatorWrapper.hpp" 
   44 #include "Thyra_DetachedSpmdVectorView.hpp" 
   45 #include "Thyra_DefaultProductVector.hpp" 
   46 #include "Thyra_ProductVectorSpaceBase.hpp" 
   47 #include "Thyra_SpmdVectorBase.hpp" 
   48 #include "Thyra_EpetraLinearOp.hpp" 
   51 #  include "Epetra_MpiComm.h" 
   53 #include "Epetra_SerialComm.h" 
   54 #include "Epetra_Vector.h" 
   57 #  include "Teuchos_DefaultMpiComm.hpp" 
   59 #include "Teuchos_DefaultSerialComm.hpp" 
   71   : useTranspose_(false),
 
   73     range_(thyraOp->range()),
 
   74     domain_(thyraOp->domain()),
 
   75     comm_(getEpetraComm(*thyraOp)),
 
   78     label_(thyraOp->description())
 
   86   using Teuchos::rcpFromPtr;
 
   87   using Teuchos::rcp_dynamic_cast;
 
   89   const int numVecs = x.NumVectors();
 
   92     "epetraToThyra does not work with MV dimension != 1");
 
   95     castOrCreateNonconstProductVectorBase(rcpFromPtr(thyraVec));
 
  103   const int numBlocks = prodThyraVec->productSpace()->numBlocks();
 
  104   for (
int b = 0; b < numBlocks; ++b) {
 
  110     const int localNumElems = spmd_vs_b->localSubDim();
 
  111     for (
int i=0; i < localNumElems; ++i) {
 
  112       thyraData[i] = epetraData[i+offset];
 
  114     offset += localNumElems;
 
  124   using Teuchos::rcpFromRef;
 
  125   using Teuchos::rcp_dynamic_cast;
 
  127   const int numVecs = x.NumVectors();
 
  130     "epetraToThyra does not work with MV dimension != 1");
 
  133     castOrCreateProductVectorBase(rcpFromRef(thyraVec));
 
  139   const int numBlocks = prodThyraVec->productSpace()->numBlocks();
 
  140   for (
int b = 0; b < numBlocks; ++b) {
 
  146     const int localNumElems = spmd_vs_b->localSubDim();
 
  147     for (
int i=0; i < localNumElems; ++i) {
 
  148       epetraData[i+offset] = thyraData[i];
 
  150     offset += localNumElems;
 
  164     opRange = ( !useTranspose_ ? range_ : domain_ ),
 
  165     opDomain = ( !useTranspose_ ? domain_ : range_ );
 
  185     "EpetraOperatorWrapper::ApplyInverse not implemented");
 
  193     "EpetraOperatorWrapper::NormInf not implemated");
 
  206   using Teuchos::rcp_dynamic_cast;
 
  209   using Teuchos::MpiComm;
 
  218     vs = prod_vs->getBlock(0);
 
  232 Thyra::makeEpetraWrapper(
const RCP<
const LinearOpBase<double> > &thyraOp)
 
Create an explicit detached mutable (non-const) view of all of the local elements on this process of ...
virtual RCP< const VectorSpaceBase< Scalar > > range() const =0
Return a smart pointer for the range space for this operator. 
void copyEpetraIntoThyra(const Epetra_MultiVector &x, const Ptr< VectorBase< double > > &thyraVec) const 
EpetraOperatorWrapper(const RCP< const LinearOpBase< double > > &thyraOp)
RCP< const Epetra_Map > get_Epetra_Map(const VectorSpaceBase< double > &vs, const RCP< const Epetra_Comm > &comm)
Get (or create) an Epetra_Map object given an VectorSpaceBase object an optionally an extra Epetra_Co...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Use the non-transposed operator. 
Use the transposed operator with complex-conjugate clements (same as TRANS for real scalar types)...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Abstract interface for finite-dimensional dense vectors. 
RCP< const EpetraLinearOp > epetraLinearOp(const RCP< const Epetra_Operator > &op, EOpTransp opTrans=NOTRANS, EApplyEpetraOpAs applyAs=EPETRA_OP_APPLY_APPLY, EAdjointEpetraOp adjointSupport=EPETRA_OP_ADJOINT_SUPPORTED, const RCP< const VectorSpaceBase< double > > &range=Teuchos::null, const RCP< const VectorSpaceBase< double > > &domain=Teuchos::null)
Dynamically allocate a nonconst EpetraLinearOp to wrap a const Epetra_Operator object. 
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const 
RCP< const Epetra_Comm > get_Epetra_Comm(const Teuchos::Comm< Ordinal > &comm)
Get (or create) and Epetra_Comm given a Teuchos::Comm object. 
void copyThyraIntoEpetra(const VectorBase< double > &thyraVec, Epetra_MultiVector &x) const 
bool nonnull(const boost::shared_ptr< T > &p)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const 
Create an explicit detached non-mutable (const) view of all of the local elements on this process of ...
Base abstract VectorSpaceBase class for all SPMD-based vector spaces.