Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Related Functions | List of all members
Thyra::DefaultSpmdVectorSpace< Scalar > Class Template Reference

Concrete implementation of an SPMD vector space subclass which creates DefaultSpmdVector and DefaultSpmdMultiVector objects. More...

#include <Thyra_DefaultSpmdVectorSpace_decl.hpp>

Inheritance diagram for Thyra::DefaultSpmdVectorSpace< Scalar >:
Inheritance graph
[legend]

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP< DefaultSpmdVectorSpace
< Scalar > > 
defaultSpmdVectorSpace ()
 Nonmember consturctor that creats an uninitialized vector space. More...
 
template<class Scalar >
RCP< DefaultSpmdVectorSpace
< Scalar > > 
defaultSpmdVectorSpace (const Ordinal dim)
 Nonmember consturctor that creats a serial vector space. More...
 
template<class Scalar >
RCP< DefaultSpmdVectorSpace
< Scalar > > 
defaultSpmdVectorSpace (const RCP< const Teuchos::Comm< Ordinal > > &comm, const Ordinal localSubDim, const Ordinal globalDim, const bool isLocallyReplicated=false)
 Nonmember consturctor function that creates a distributed or locally-replicated parallel vector space. More...
 
template<class Scalar >
RCP< DefaultSpmdVectorSpace
< Scalar > > 
locallyReplicatedDefaultSpmdVectorSpace (const RCP< const Teuchos::Comm< Ordinal > > &comm, const Ordinal globalDim)
 Nonmember consturctor function that creates a locally-replicated parallel vector space. More...
 

Constructors and initializers

static RCP
< DefaultSpmdVectorSpace
< Scalar > > 
create ()
 Create with weak ownership to self. More...
 
void initialize (const Ordinal dim)
 Initialize a serial space. More...
 
void initialize (const RCP< const Teuchos::Comm< Ordinal > > &comm, const Ordinal localSubDim, const Ordinal globalDim, const bool isLocallyReplicated=false)
 Initialize an SPMD space. More...
 
void uninitialize ()
 Set to an uninitialized state. More...
 

Public overridden from VectorSpaceBase

bool hasInCoreView (const Range1D &rng, const EViewType viewType, const EStrideType strideType) const
 Returns true if all the elements in rng are in this process. More...
 
RCP< const VectorSpaceBase
< Scalar > > 
clone () const
 

Protected overridden from VectorSpaceBase

RCP< VectorBase< Scalar > > createMember () const
 
RCP< MultiVectorBase< Scalar > > createMembers (int numMembers) const
 
RCP< VectorBase< Scalar > > createMemberView (const RTOpPack::SubVectorView< Scalar > &raw_v) const
 
RCP< const VectorBase< Scalar > > createMemberView (const RTOpPack::ConstSubVectorView< Scalar > &raw_v) const
 
RCP< MultiVectorBase< Scalar > > createMembersView (const RTOpPack::SubMultiVectorView< Scalar > &raw_mv) const
 
RCP< const MultiVectorBase
< Scalar > > 
createMembersView (const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv) const
 

Public overridden from SpmdVectorSpaceDefaultBase

RCP< const Teuchos::Comm
< Ordinal > > 
getComm () const
 
Ordinal localSubDim () const
 

Additional Inherited Members

- Public Member Functions inherited from Thyra::SpmdVectorSpaceDefaultBase< Scalar >
 SpmdVectorSpaceDefaultBase ()
 
Ordinal localOffset () const
 
Ordinal mapCode () const
 
bool isLocallyReplicated () const
 Returns true if vector space is locally replicated space. More...
 
std::string description () const
 
Ordinal dim () const
 Returns the sum of the local number of elements on every process. More...
 
Teuchos::RCP< const
VectorSpaceFactoryBase< Scalar > > 
smallVecSpcFcty () const
 Returns a DefaultSpmdVectorSpaceFactory object that has been given getComm(). More...
 
bool isCompatible (const VectorSpaceBase< Scalar > &vecSpc) const
 Checks the general compatibility of parallel (or serial on one process) Spmd-based vector spaces. More...
 
- Public Member Functions inherited from Thyra::SpmdVectorSpaceBase< Scalar >
- Public Member Functions inherited from Thyra::VectorSpaceBase< Scalar >
void scalarProds (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds_out) const
 Return the scalar product of each column in two multi-vectors in the vector space. More...
 
- Public Member Functions inherited from Thyra::ScalarProdVectorSpaceBase< Scalar >
 ScalarProdVectorSpaceBase ()
 Construct to use dot product as the default. More...
 
 ScalarProdVectorSpaceBase (const RCP< const ScalarProdBase< Scalar > > &scalarProd)
 Construct with a different scalar product. More...
 
virtual void setScalarProd (const RCP< const ScalarProdBase< Scalar > > &scalarProd)
 Set a different scalar product. More...
 
RCP< const ScalarProdBase
< Scalar > > 
getScalarProd () const
 Return the current scalar product. More...
 
bool isEuclidean () const
 Returns getScalarProd()->isEuclidean() More...
 
Scalar scalarProd (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y) const
 Returns getScalarProd()->scalarProd(x,y) More...
 
void scalarProdsImpl (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds_out) const
 Calls getScalarProd()->scalarProds(X,Y,scalar_prods) More...
 
- Protected Member Functions inherited from Thyra::SpmdVectorSpaceDefaultBase< Scalar >
virtual void updateState (const Ordinal globalDim, const bool isLocallyReplicated=false)
 This function must be called whenever the state of this changes and some internal state must be updated. More...
 
- Protected Member Functions inherited from Thyra::VectorSpaceBase< Scalar >
- Protected Member Functions inherited from Thyra::VectorSpaceDefaultBase< Scalar >

Detailed Description

template<class Scalar>
class Thyra::DefaultSpmdVectorSpace< Scalar >

Concrete implementation of an SPMD vector space subclass which creates DefaultSpmdVector and DefaultSpmdMultiVector objects.

This is a simple but yet fully general and functional concrete subclass of SpmdVectorSpaceBase that returns DefaultSpmdMultiVector objects from createMembers() and DefaultSpmdVector objects from createMember().

See the function initialize() that describes the different kinds of distributions this class can handle.

Definition at line 69 of file Thyra_DefaultSpmdVectorSpace_decl.hpp.

Member Function Documentation

template<class Scalar >
RCP< DefaultSpmdVectorSpace< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::create ( )
static

Create with weak ownership to self.

Definition at line 58 of file Thyra_DefaultSpmdVectorSpace_def.hpp.

template<class Scalar >
void Thyra::DefaultSpmdVectorSpace< Scalar >::initialize ( const Ordinal  dim)

Initialize a serial space.

Parameters
dim[in] Gives the dimension of the vector space.

Equivalent to calling this->initialize(Teuchos::null,dim,dim)

Definition at line 67 of file Thyra_DefaultSpmdVectorSpace_def.hpp.

template<class Scalar >
void Thyra::DefaultSpmdVectorSpace< Scalar >::initialize ( const RCP< const Teuchos::Comm< Ordinal > > &  comm,
const Ordinal  localSubDim,
const Ordinal  globalDim,
const bool  isLocallyReplicated = false 
)

Initialize an SPMD space.

Parameters
comm[in] The communicator. This object must be maintained by the client the entire time that this is in use.
localSubDim[in] The number of elements in the local process. This number can be different in every process.
globalDim[in] Gives the number of global elements in the vector if globalDim > 0. If globalDim < 0 then the global dimension is determined by the above argument localSubDim but requires a global communication to do so (i.e. Spmd_Allreduce()).

Preconditions:

  • localSubDim > 0
  • globalDim != 0
  • [comm.get() != NULL && globalDim > 0] globalDim >= localSubDim

Postconditions:

  • this->getComm().get() == comm.get()
  • this->localSubDim() == localSubDim
  • [comm.get() == NULL] this->dim() == localSubDim
  • [comm.get() != NULL && globalDim > 0] this->dim() == globalDim
  • [comm.get() != NULL && globalDim < 0] this->dim() == sum(localSubDim[i],i=0...numProc-1)

This function supports three different types of use-cases:

  • comm.get()==NULL : Serial (i.e. single process) vectors where this->dim() == localSubDim.

  • comm.get()!=NULL && globalDim < 0 : Distributed-memory vectors where this->dim() is equal to the sum of the localSubDim arguments in each process. This will result in a call to Spmd_Allreduce() inside of this function.

  • comm.get()!=NULL && globalDim > 0 : Distributed-memory vectors where this->dim() returns globalDim. This will not result in a call to Teuchos::reduceAll() inside this function and therefore the client had better be sure that globalDim is consistent with localSubDim in each process.

  • comm.get()!=NULL && globalDim == localSubDim : Locally-replicated distributed-memory vectors where this->dim() == globalDim == localSubDim.

Definition at line 76 of file Thyra_DefaultSpmdVectorSpace_def.hpp.

template<class Scalar >
void Thyra::DefaultSpmdVectorSpace< Scalar >::uninitialize ( )

Set to an uninitialized state.

Postconditions:

Definition at line 100 of file Thyra_DefaultSpmdVectorSpace_def.hpp.

template<class Scalar >
bool Thyra::DefaultSpmdVectorSpace< Scalar >::hasInCoreView ( const Range1D rng,
const EViewType  viewType,
const EStrideType  strideType 
) const
virtual

Returns true if all the elements in rng are in this process.

Reimplemented from Thyra::VectorSpaceBase< Scalar >.

Definition at line 225 of file Thyra_DefaultSpmdVectorSpace_def.hpp.

template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::clone ( ) const
virtual

Reimplemented from Thyra::VectorSpaceBase< Scalar >.

Definition at line 237 of file Thyra_DefaultSpmdVectorSpace_def.hpp.

template<class Scalar >
RCP< VectorBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::createMember ( ) const
protectedvirtual
template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::createMembers ( int  numMembers) const
protectedvirtual
template<class Scalar >
RCP< VectorBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::createMemberView ( const RTOpPack::SubVectorView< Scalar > &  raw_v) const
protectedvirtual
template<class Scalar >
RCP< const VectorBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::createMemberView ( const RTOpPack::ConstSubVectorView< Scalar > &  raw_v) const
protectedvirtual
template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::createMembersView ( const RTOpPack::SubMultiVectorView< Scalar > &  raw_mv) const
protectedvirtual
template<class Scalar >
RCP< const MultiVectorBase< Scalar > > Thyra::DefaultSpmdVectorSpace< Scalar >::createMembersView ( const RTOpPack::ConstSubMultiVectorView< Scalar > &  raw_mv) const
protectedvirtual
template<class Scalar >
RCP< const Teuchos::Comm< Ordinal > > Thyra::DefaultSpmdVectorSpace< Scalar >::getComm ( ) const
virtual
template<class Scalar >
Ordinal Thyra::DefaultSpmdVectorSpace< Scalar >::localSubDim ( ) const
virtual

Friends And Related Function Documentation

template<class Scalar >
RCP< DefaultSpmdVectorSpace< Scalar > > defaultSpmdVectorSpace ( )
related

Nonmember consturctor that creats an uninitialized vector space.

Definition at line 238 of file Thyra_DefaultSpmdVectorSpace_decl.hpp.

template<class Scalar >
RCP< DefaultSpmdVectorSpace< Scalar > > defaultSpmdVectorSpace ( const Ordinal  dim)
related

Nonmember consturctor that creats a serial vector space.

Definition at line 250 of file Thyra_DefaultSpmdVectorSpace_decl.hpp.

template<class Scalar >
RCP< DefaultSpmdVectorSpace< Scalar > > defaultSpmdVectorSpace ( const RCP< const Teuchos::Comm< Ordinal > > &  comm,
const Ordinal  localSubDim,
const Ordinal  globalDim,
const bool  isLocallyReplicated = false 
)
related

Nonmember consturctor function that creates a distributed or locally-replicated parallel vector space.

Definition at line 266 of file Thyra_DefaultSpmdVectorSpace_decl.hpp.

template<class Scalar >
RCP< DefaultSpmdVectorSpace< Scalar > > locallyReplicatedDefaultSpmdVectorSpace ( const RCP< const Teuchos::Comm< Ordinal > > &  comm,
const Ordinal  globalDim 
)
related

Nonmember consturctor function that creates a locally-replicated parallel vector space.

Definition at line 286 of file Thyra_DefaultSpmdVectorSpace_decl.hpp.


The documentation for this class was generated from the following files: