Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Enumerations | Functions
Thyra/Epetra Operator/Vector Adapter Code
Collaboration diagram for Thyra/Epetra Operator/Vector Adapter Code:

Classes

class  Thyra::DiagonalEpetraLinearOpWithSolveFactory
 Create a DefaultDiagonalLinearOpWithSolve out of a diagonal Epetra_RowMatrix object. More...
 
class  Thyra::EpetraLinearOp
 Concrete LinearOpBase adapter subclass for Epetra_Operator object. More...
 
class  Thyra::EpetraLinearOpBase
 Abstract base class for all LinearOpBase objects that can return an Epetra_Operator view of themselves and details about how to apply the view. More...
 
class  Thyra::EpetraOperatorViewExtractorBase
 Strategy interface for extracting an Epetra_Operator view out of a Thyra::LinearOpBase<double> object. More...
 
class  Thyra::EpetraOperatorViewExtractorStd
 Standard strategy subclass for extracting an Epetra_Operator view out of a Thyra::LinearOpBase<double> object by dynamic casting to the EpetraLinearOpBase interface. More...
 
class  Thyra::EpetraOperatorWrapper
 Implements the Epetra_Operator interface with a Thyra LinearOperator. More...
 

Enumerations

enum  Thyra::EAdjointEpetraOp { Thyra::EPETRA_OP_ADJOINT_SUPPORTED, Thyra::EPETRA_OP_ADJOINT_UNSUPPORTED }
 Determine if adjoints are supported on Epetra_Opeator or not. More...
 
enum  Thyra::EApplyEpetraOpAs { Thyra::EPETRA_OP_APPLY_APPLY, Thyra::EPETRA_OP_APPLY_APPLY_INVERSE }
 Determine how the apply an Epetra_Operator as a linear operator. More...
 

Functions

RCP< const Teuchos::Comm
< Ordinal > > 
Thyra::create_Comm (const RCP< const Epetra_Comm > &epetraComm)
 Given an Epetra_Comm object, return an equivalent Teuchos::Comm object. More...
 
RCP< const VectorSpaceBase
< double > > 
Thyra::create_VectorSpace (const RCP< const Epetra_Map > &epetra_map)
 Create an VectorSpaceBase object given an Epetra_Map object. More...
 
RCP< const VectorSpaceBase
< double > > 
Thyra::create_LocallyReplicatedVectorSpace (const RCP< const VectorSpaceBase< double > > &parentSpace, const int dim)
 Create a VectorSpaceBase object that creates locally replicated vector objects. More...
 
RCP< VectorBase< double > > Thyra::create_Vector (const RCP< Epetra_Vector > &epetra_v, const RCP< const VectorSpaceBase< double > > &space=Teuchos::null)
 Create a non-const VectorBase object from a non-const Epetra_Vector object. More...
 
RCP< const VectorBase< double > > Thyra::create_Vector (const RCP< const Epetra_Vector > &epetra_v, const RCP< const VectorSpaceBase< double > > &space=Teuchos::null)
 Create an const VectorBase wrapper object for a const Epetra_Vector object. More...
 
RCP< MultiVectorBase< double > > Thyra::create_MultiVector (const RCP< Epetra_MultiVector > &epetra_mv, const RCP< const VectorSpaceBase< double > > &range=Teuchos::null, const RCP< const VectorSpaceBase< double > > &domain=Teuchos::null)
 Create a non-const MultiVectorBase object from a non-const Epetra_MultiVector object. More...
 
RCP< const MultiVectorBase
< double > > 
Thyra::create_MultiVector (const RCP< const Epetra_MultiVector > &epetra_mv, const RCP< const VectorSpaceBase< double > > &range=Teuchos::null, const RCP< const VectorSpaceBase< double > > &domain=Teuchos::null)
 Create an const MultiVectorBase wrapper object for a const Epetra_MultiVector object. More...
 
RCP< const Epetra_CommThyra::get_Epetra_Comm (const Teuchos::Comm< Ordinal > &comm)
 Get (or create) and Epetra_Comm given a Teuchos::Comm object. More...
 
RCP< const Epetra_MapThyra::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_Comm object. More...
 
RCP< Epetra_VectorThyra::get_Epetra_Vector (const Epetra_Map &map, const RCP< VectorBase< double > > &v)
 Get a non-const Epetra_Vector view from a non-const VectorBase object if possible. More...
 
RCP< const Epetra_VectorThyra::get_Epetra_Vector (const Epetra_Map &map, const RCP< const VectorBase< double > > &v)
 Get a const Epetra_Vector view from a const VectorBase object if possible. More...
 
RCP< Epetra_MultiVectorThyra::get_Epetra_MultiVector (const Epetra_Map &map, const RCP< MultiVectorBase< double > > &mv)
 Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible. More...
 
RCP< const Epetra_MultiVectorThyra::get_Epetra_MultiVector (const Epetra_Map &map, const RCP< const MultiVectorBase< double > > &mv)
 Get a const Epetra_MultiVector view from a const MultiVectorBase object if possible. More...
 
Teuchos::RCP< Epetra_MultiVectorThyra::get_Epetra_MultiVector (const Epetra_Map &map, MultiVectorBase< double > &mv)
 Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible where the client must maintain the memory of the input multivector. More...
 
Teuchos::RCP< const
Epetra_MultiVector
Thyra::get_Epetra_MultiVector (const Epetra_Map &map, const MultiVectorBase< double > &mv)
 Get a const Epetra_MultiVector view from a const MultiVectorBase object if possible where the client must maintain the memory of the input multivector. More...
 
const std::string Thyra::toString (const EAdjointEpetraOp adjointEpetraOp)
 
const std::string Thyra::toString (const EApplyEpetraOpAs applyEpetraOpAs)
 
template<class Scalar >
Teuchos::RCP< Epetra_OperatorThyra::get_Epetra_Operator (LinearOpBase< Scalar > &)
 Get smart pointer to non-const Epetra_Operator object from reference to a non-const EpetraLinearOp accessed through its LinearOpBase interface. More...
 
template<class Scalar >
Teuchos::RCP< const
Epetra_Operator
Thyra::get_Epetra_Operator (const LinearOpBase< Scalar > &)
 Get smart pointer to const Epetra_Operator object from reference to a const EpetraLinearOp accessed through its LinearOpBase interface. More...
 

Detailed Description

The following functions and classes are used to create Thyra objects which wrap (or adapt) Epetra objects:

The above adapter code is based directly from the general Thyra Operator/Vector Subclasses for SPMD Distributed-Memory Platforms. Therefore these Epetra adapted objects are automatically compatible with any other such MPI-based SPMD adapter subclasses.

There is, however, one issue that requires a little care and that is using arbitrary Thyra::VectorBase and Thyra::MultiVectorBase objects with the vector and multi-vector versions of Thyra::EpetraLinearOp::apply(). The issue is that the underlying Epetra_Operator::Apply() function can only accept Epetra_MultiVector objects. The utility functions Thyra::get_Epetra_MultiVector() return an Epetra_MultiVector view of any Thyra::MultiVectorBase object with a compatible range space. Studying the implementations of these utility functions will show you how simple it is to provide for this type of interoperabiity. This type of interoperabiity machinary should also be used for other types of concrete adapter subclasses.

Other types of code will need to extract an Epetra_Vector view of a Thyra::VectorBase object. For this purpose the Thyra::get_Epetra_Vector() functions are provided.

The utility functions Thyra::get_Epetra_Operator() are also provided that encapsulate the extraction of an "adapted" Epetra_Operator object out of a Thyra::EpetraLinearOp object through its Thyra::LinearOpBase base interface. This is a common type of activity in object-oriented programming.

Examples

Enumeration Type Documentation

Determine if adjoints are supported on Epetra_Opeator or not.

Enumerator
EPETRA_OP_ADJOINT_SUPPORTED 

Adjoint supported.

EPETRA_OP_ADJOINT_UNSUPPORTED 

Adjoint not supported.

Definition at line 64 of file Thyra_EpetraTypes.hpp.

Determine how the apply an Epetra_Operator as a linear operator.

Enumerator
EPETRA_OP_APPLY_APPLY 

Apply using Epetra_Operator::Apply(...)

EPETRA_OP_APPLY_APPLY_INVERSE 

Apply using Epetra_Operator::ApplyInverse(...)

Definition at line 93 of file Thyra_EpetraTypes.hpp.

Function Documentation

Teuchos::RCP< const Teuchos::Comm< Thyra::Ordinal > > Thyra::create_Comm ( const RCP< const Epetra_Comm > &  epetraComm)

Given an Epetra_Comm object, return an equivalent Teuchos::Comm object.

If a successful conversion could not be performed then return.get()==NULL.

Definition at line 102 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Thyra::VectorSpaceBase< double > > Thyra::create_VectorSpace ( const RCP< const Epetra_Map > &  epetra_map)

Create an VectorSpaceBase object given an Epetra_Map object.

Parameters
epetra_map[in] The Epetra map defining the partitioning of elements to processors.

Preconditions:

  • epetra_map.get() != NULL

Postconditions:

  • return.get() != NULL
  • The RCP object epetra_map is copied into the return object and therefore a memory of epetra_map is kept.

This uses an Epetra_Map object to initialize a compatible VectorSpaceBase object.

The fact that this function only accepts an Epetra_Map object means that only maps that have elements of size one can be used to define a vector space. General Epetra_BlockMaps can not be used. This is not a serious limitation since Epetra_Operator's domain and range maps are of type Epetra_Map.

This function works properly even if Epetra is not compiled with support for SPMD (i.e. HAVE_MPI is not defined when compiling and linking). If SPMD support is not compiled into Epetra, then a serial implementation of the communication is used instead.

Definition at line 139 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Thyra::VectorSpaceBase< double > > Thyra::create_LocallyReplicatedVectorSpace ( const RCP< const VectorSpaceBase< double > > &  parentSpace,
const int  dim 
)

Create a VectorSpaceBase object that creates locally replicated vector objects.

Parameters
parentSpace[in] The vector space that will be used to create the smaller locally-replicated vector space.
dim[in] The dimension of the locally replicated vector space.

Note: This vector space will be compatible with the domain space of a multivector. which has the range space parentSpace.

Definition at line 177 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< Thyra::VectorBase< double > > Thyra::create_Vector ( const RCP< Epetra_Vector > &  epetra_v,
const RCP< const VectorSpaceBase< double > > &  space = Teuchos::null 
)

Create a non-const VectorBase object from a non-const Epetra_Vector object.

Parameters
epetra_v[in] Smart pointer to the Epetra_Vector object to wrap.
space[in] The vector space that is compatible with epetra_v->Map().

Precondiitions:

  • [epetra_v.get()!=NULL] space.get()!=NULL

Postconditions:

  • [epetra_v.get()==NULL] return.get()==NULL
  • [epetra_v.get()!=NULL] return.get()!=NULL

Returns
The returned RCP object contains a copy of the input RCP<Epetra_Vector> wrapped Epetra_Vector object. It is also stated that *epetra_v will only be guaranteed to be modifed after the last RCP to the returned VectorBase is destroyed. In addition, *return is only valid as long as one RefCoutPtr wrapper object still exits.

Definition at line 193 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Thyra::VectorBase< double > > Thyra::create_Vector ( const RCP< const Epetra_Vector > &  epetra_v,
const RCP< const VectorSpaceBase< double > > &  space = Teuchos::null 
)

Create an const VectorBase wrapper object for a const Epetra_Vector object.

Parameters
epetra_v[in] Smart pointer to the Epetra_Vector object to wrap.
space[in] The vector space that is compatible with epetra_v->Map().

Precondiitions:

  • [epetra_v.get()!=NULL] space.get()!=NULL

Postconditions:

  • [epetra_v.get()==NULL] return.get()==NULL
  • [epetra_v.get()!=NULL] return.get()!=NULL

Returns
The returned RCP object contains a copy of the input RCP<Epetra_Vector> wrapped Epetra_Vector object. In addition, *return is only valid as long as one RefCoutPtr wrapper object still exits.

Definition at line 229 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< Thyra::MultiVectorBase< double > > Thyra::create_MultiVector ( const RCP< Epetra_MultiVector > &  epetra_mv,
const RCP< const VectorSpaceBase< double > > &  range = Teuchos::null,
const RCP< const VectorSpaceBase< double > > &  domain = Teuchos::null 
)

Create a non-const MultiVectorBase object from a non-const Epetra_MultiVector object.

Parameters
epetra_mv[in] Smart pointer to the Epetra_MultiVector object to wrap.
range[in] The vector space that is compatible with epetra_mv->Map().
domain[in] The vector space that is compatible with epetra_mv.NumVectors. If domain.get()==NULL, then a space will be created internally.

Precondiitions:

  • epetra_mv.get()!=NULL
  • range.get()!=NULL
Returns
The returned RCP object contains a copy of the input RCP<Epetra_MultiVector> wrapped Epetra_MultiVector object. It is also stated that *epetra_mv will only be guaranteed to be modifed after the last RCP to the returned MultiVectorBase is destroyed. In addition, *return is only valid as long as one RefCoutPtr wrapper object still exits.

Definition at line 265 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Thyra::MultiVectorBase< double > > Thyra::create_MultiVector ( const RCP< const Epetra_MultiVector > &  epetra_mv,
const RCP< const VectorSpaceBase< double > > &  range = Teuchos::null,
const RCP< const VectorSpaceBase< double > > &  domain = Teuchos::null 
)

Create an const MultiVectorBase wrapper object for a const Epetra_MultiVector object.

Parameters
epetra_mv[in] Smart pointer to the Epetra_MultiVector object to wrap.
range[in] The vector space that is compatible with epetra_mv->Map().
domain[in] The vector space that is compatible with epetra_mv.NumVectors. If domain.get()==NULL, then a space will be created internally.

Preconditions:

  • epetra_mv.get()!=NULL
  • range.get()!=NULL
Returns
The returned RCP object contains a copy of the input RCP<Epetra_MultiVector> wrapped Epetra_MultiVector object. In addition, *return is only valid as long as one RefCoutPtr wrapper object still exits.

Definition at line 321 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Epetra_Comm > Thyra::get_Epetra_Comm ( const Teuchos::Comm< Ordinal > &  comm)

Get (or create) and Epetra_Comm given a Teuchos::Comm object.

This function returns a new Epetra_Comm object from the input Teuchos::Comm object. If an Epetra_Comm was used to create this Tpetra::Comm object using the function create_Comm() then the Epetra_Comm object returned from this function will not be the same!

NOTE: The behavior of the implementation of this function is likely not the desired behavior and would likely surpise clients and is not consistent with the other functions that actually return the same Epetra object passed in.

Definition at line 377 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Epetra_Map > Thyra::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_Comm object.

ToDo: Finish documentation!

Definition at line 431 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< Epetra_Vector > Thyra::get_Epetra_Vector ( const Epetra_Map map,
const RCP< VectorBase< double > > &  v 
)

Get a non-const Epetra_Vector view from a non-const VectorBase object if possible.

Preconditions:

  • v.get()!=NULL
  • map must be compatible with *v.space()

If a RCP<Epetra_Vector> object is already attached to the node of the smart pointer for mv then this is returned directly. If not, then a view of the data in *v is created and returned. In the latter case the smart pointer v is copied and attached to the returned RCP object. Therefore, a temporary VectorBase object can be created in the call to this function and the view in return will persist until all of the RCP objects to the returned Epetra_Vector object go away.

Note: the v object is not guaranteed to be modified until the last smart pointer to the returned Epetra_Vector object is destroyed.

Definition at line 538 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Epetra_Vector > Thyra::get_Epetra_Vector ( const Epetra_Map map,
const RCP< const VectorBase< double > > &  v 
)

Get a const Epetra_Vector view from a const VectorBase object if possible.

Preconditions:

  • v.get()!=NULL
  • map must be compatible with *v->space()

If a RCP<Epetra_Vector> object is already attached to the node of the smart pointer for mv then this is returned directly. If not, then a view of the data in *v is created and returned. In the latter case the smart pointer v is copied and attached to the returned RCP object. Therefore, a temporary VectorBase object can be created in the call to this function and the view in return will persist until all of the RCP objects to the returned Epetra_Vector object go away.

Definition at line 645 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< Epetra_MultiVector > Thyra::get_Epetra_MultiVector ( const Epetra_Map map,
const RCP< MultiVectorBase< double > > &  mv 
)

Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible.

Preconditions:

  • mv.get()!=NULL
  • map must be compatible with *mv->range()

If a RCP<Epetra_MultiVector> object is already attached to the node of the smart pointer for mv then this is returned directly. If not, then a view of the data in *mv is created and returned. In the latter case the smart pointer mv is copied and attached to the returned RCP object. Therefore, a temporary MultiVectorBase object can be created in the call to this function and the view in return will persist until all of the RCP objects to the returned Epetra_MultiVector object go away.

Note: the mv object is not guaranteed to be modified until the last smart pointer to the returned Epetra_MultiVector object is destroyed.

Definition at line 744 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Epetra_MultiVector > Thyra::get_Epetra_MultiVector ( const Epetra_Map map,
const RCP< const MultiVectorBase< double > > &  mv 
)

Get a const Epetra_MultiVector view from a const MultiVectorBase object if possible.

Preconditions:

  • mv.get()!=NULL
  • map must be compatible with *mv.range()

If a RCP<const Epetra_MultiVector> object is already attached to the node of the smart pointer for mv then this is returned directly. If not, then a view of the data in *mv is created and returned. In the latter case the smart pointer mv is copied and attached to the returned RCP object. Therefore, a temporary MultiVectorBase object can be created in the call to this function and the view in return will persist until all of the RCP objects to the returned Epetra_MultiVector object go away.

Definition at line 852 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< Epetra_MultiVector > Thyra::get_Epetra_MultiVector ( const Epetra_Map map,
MultiVectorBase< double > &  mv 
)

Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible where the client must maintain the memory of the input multivector.

Preconditions:

  • map must be compatible with *mv.range()

This function trys to dynamic cast some some known interfaces classes where data can be directly accessed and no RCP magic needs to be used. This results in improved performance in time-critical use cases (like when called for EpetraLinearOp in the inner loop of a Krylov solver).

If this function can not dynamic cast to the direct data access interfaces it punts and calls the more general (but more expensive) get_Epetra_MultiVector() function.

Note: the mv object is not guaranteed to be modified until the last smart pointer to the returned Epetra_MultiVector object is destroyed.

Definition at line 957 of file Thyra_EpetraThyraWrappers.cpp.

Teuchos::RCP< const Epetra_MultiVector > Thyra::get_Epetra_MultiVector ( const Epetra_Map map,
const MultiVectorBase< double > &  mv 
)

Get a const Epetra_MultiVector view from a const MultiVectorBase object if possible where the client must maintain the memory of the input multivector.

Preconditions:

  • map must be compatible with *mv.range()

This function trys to dynamic cast some some known interfaces classes where data can be directly accessed and no RCP magic needs to be used. This results in improved performance in time-critical use cases (like when called for EpetraLinearOp in the inner loop of a Krylov solver).

If this function can not dynamic cast to the direct data access interfaces it punts and calls the more general (but more expensive) get_Epetra_MultiVector() function.

Definition at line 986 of file Thyra_EpetraThyraWrappers.cpp.

const std::string Thyra::toString ( const EAdjointEpetraOp  adjointEpetraOp)
inline

Definition at line 75 of file Thyra_EpetraTypes.hpp.

const std::string Thyra::toString ( const EApplyEpetraOpAs  applyEpetraOpAs)
inline

Definition at line 104 of file Thyra_EpetraTypes.hpp.

template<class Scalar >
Teuchos::RCP<Epetra_Operator> Thyra::get_Epetra_Operator ( LinearOpBase< Scalar > &  )

Get smart pointer to non-const Epetra_Operator object from reference to a non-const EpetraLinearOp accessed through its LinearOpBase interface.

Parameters
op[in] Reference to operator to extract Epetra_Operator out of.

Preconditions:

  • dynamic_cast<EpetraLinearOp*>(&op) != NULL

This function is designed to provide an easy way for non-C++ experts to get at the Epetra_Operator object that was stuck into an EpetraLinearOp object.

If the dynamic cast fails then a std::bad_cast exception is thrown containing a very detailed error message as to why the cast failed.

This function is simple enough and developers can see what needs to be done to accomplish this type of access by looking at the source code by clicking on:

Definition at line 80 of file Thyra_get_Epetra_Operator.hpp.

template<class Scalar >
Teuchos::RCP<const Epetra_Operator> Thyra::get_Epetra_Operator ( const LinearOpBase< Scalar > &  )

Get smart pointer to const Epetra_Operator object from reference to a const EpetraLinearOp accessed through its LinearOpBase interface.

Parameters
op[in] Reference to operator to extract Epetra_Operator out of.

Preconditions:

  • dynamic_cast<const EpetraLinearOp*>(&op) != NULL

This function is designed to provide an easy way for non-C++ experts to get at the Epetra_Operator object that was stuck into an EpetraLinearOp object.

If the dynamic cast fails then a std::bad_cast exception is thrown containing a very detailed error message as to why the cast failed.

This function is simple enough and developers can see what needs to be done to accomplish this type of access by looking at the source code by clicking on:

Definition at line 122 of file Thyra_get_Epetra_Operator.hpp.