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

More...

#include <Thyra_DefaultClusteredSpmdProductVector_decl.hpp>

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

Constructors/Intializers/Accessors

 DefaultClusteredSpmdProductVectorSpace ()
 
 DefaultClusteredSpmdProductVectorSpace (const RCP< const Teuchos::Comm< Ordinal > > &intraClusterComm, const int clusterRootRank, const RCP< const Teuchos::Comm< Ordinal > > &interClusterComm, const int numBlocks, const RCP< const VectorSpaceBase< Scalar > > vecSpaces[])
 
void initialize (const RCP< const Teuchos::Comm< Ordinal > > &intraClusterComm, const int clusterRootRank, const RCP< const Teuchos::Comm< Ordinal > > &interClusterComm, const int numBlocks, const RCP< const VectorSpaceBase< Scalar > > vecSpaces[])
 Initalize. More...
 
RCP< const Teuchos::Comm
< Ordinal > > 
intraClusterComm () const
 
int clusterRootRank () const
 
RCP< const Teuchos::Comm
< Ordinal > > 
interClusterComm () const
 
int clusterSubDim () const
 
int clusterOffset () const
 

Overridden form Teuchos::Describable

std::string description () const
 

Public overridden from VectorSpaceBase

Ordinal dim () const
 
bool isCompatible (const VectorSpaceBase< Scalar > &vecSpc) const
 
RCP< const
VectorSpaceFactoryBase< Scalar > > 
smallVecSpcFcty () const
 
Scalar scalarProd (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y) const
 
void scalarProdsImpl (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds) const
 
bool isEuclidean () const
 
bool hasInCoreView (const Range1D &rng, const EViewType viewType, const EStrideType strideType) const
 
RCP< const VectorSpaceBase
< Scalar > > 
clone () const
 

Protected overridden from ProductVectorSpaceBase

int numBlocks () const
 
RCP< const VectorSpaceBase
< Scalar > > 
getBlock (const int k) const
 

Protected overridden from VectorSpaceBase

RCP< VectorBase< Scalar > > createMember () const
 
RCP< MultiVectorBase< Scalar > > createMembers (int numMembers) const
 

Additional Inherited Members

- 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...
 
virtual RCP< MultiVectorBase
< Scalar > > 
createCachedMembersView (const RTOpPack::SubMultiVectorView< Scalar > &raw_mv) const
 Create a (possibly) cached multi-vector member that is a non-const view of raw multi-vector data. The caching mechanism must be implemented by child classes, by default this just calls the regular createMembersView. More...
 
virtual RCP< const
MultiVectorBase< Scalar > > 
createCachedMembersView (const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv) const
 Create a (possibly) cached multi-vector member that is a const view of raw multi-vector data. The caching mechanism must be implemented by child classes, by default this just calls the regular createMembersView. More...
 
- Protected Member Functions inherited from Thyra::VectorSpaceBase< Scalar >
- Protected Member Functions inherited from Thyra::VectorSpaceDefaultBase< Scalar >
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
 

Detailed Description

template<class Scalar>
class Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >

Concrete subclass of VectorSpaceBase that takes a collection of individual VectorSpaceBase objects distributed over many different processes and creates a single vector space.

What this class does is to take different vector space objects distributed over a set of different clusters of processes and then creates a single vector space object.

Let numClusters be the number of clusters that the processes in the global communicator are partitioned into. Each cluster of processes is represented by another communicator known in this process as intraClusterComm. The communicator that links the root of each cluster is called interClusterComm. All communication is performed with just these two communicators. There is no overall global communicator that encompasses all of the processes.

Within a cluster of processes, the number of constituent vector space objects must be the same. Let v[0]..v[numBlocks-1] be the numBlocks VectorBase objects for constituent vectors in this cluster of processes. There is no assumption make whatsoever about the natrue of the VectorSpaceBase objects or the [Multi]VectorBase objects that they create.

ToDo: Finish documentation!

The default copy constructor and assign operator are allowed and they work correctly and perform shallow copies of the constituent vector spaces!

Definition at line 51 of file Thyra_DefaultClusteredSpmdProductVector_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::DefaultClusteredSpmdProductVectorSpace ( const RCP< const Teuchos::Comm< Ordinal > > &  intraClusterComm,
const int  clusterRootRank,
const RCP< const Teuchos::Comm< Ordinal > > &  interClusterComm,
const int  numBlocks,
const RCP< const VectorSpaceBase< Scalar > >  vecSpaces[] 
)

Member Function Documentation

template<class Scalar >
void Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::initialize ( const RCP< const Teuchos::Comm< Ordinal > > &  intraClusterComm,
const int  clusterRootRank,
const RCP< const Teuchos::Comm< Ordinal > > &  interClusterComm,
const int  numBlocks,
const RCP< const VectorSpaceBase< Scalar > >  vecSpaces[] 
)

Initalize.

Parameters
intraClusterComm[in] The communicator over just this cluster.
clusterRootRank[in] The rank of the root process in *interClusterComm for this cluster of processes. This is also the process that in included in *interClusterComm (which has a different rank ther obviously).
interClusterComm[in] Defines the communicator between the root processes of all of the clusters. For the root process of each cluster *interClusterComm!=Spmd_COMM_NULL, otherwise interClusterComm==Teuchos::null or *interClusterComm==Spmd_COMM_NULL.
numBlocks[in] Number of vector space blocks for this cluster of processes.
vecSpaces[in] Array (length numBlocks) of the vector spaces for this cluster of processes.

Definition at line 80 of file Thyra_DefaultClusteredSpmdProductVectorSpace_def.hpp.

template<class Scalar >
RCP< const Teuchos::Comm< Ordinal > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::intraClusterComm ( ) const
template<class Scalar >
int Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::clusterRootRank ( ) const
template<class Scalar >
RCP< const Teuchos::Comm< Ordinal > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::interClusterComm ( ) const
template<class Scalar >
int Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::clusterSubDim ( ) const

The sum of the dimensions of the block vector spaces in this cluster.

Definition at line 261 of file Thyra_DefaultClusteredSpmdProductVectorSpace_decl.hpp.

template<class Scalar >
int Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::clusterOffset ( ) const

The offset of the first element in the first constituent vector in this cluster in the w.r.t. the global vector.

Definition at line 268 of file Thyra_DefaultClusteredSpmdProductVectorSpace_decl.hpp.

template<class Scalar >
std::string Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::description ( ) const
virtual

Reimplemented from Teuchos::Describable.

Definition at line 131 of file Thyra_DefaultClusteredSpmdProductVectorSpace_def.hpp.

template<class Scalar >
Ordinal Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::dim ( ) const
virtual
template<class Scalar >
bool Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::isCompatible ( const VectorSpaceBase< Scalar > &  vecSpc) const
virtual
template<class Scalar >
Teuchos::RCP< const VectorSpaceFactoryBase< Scalar > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::smallVecSpcFcty ( ) const
virtual
template<class Scalar >
Scalar Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::scalarProd ( const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
) const
virtual
template<class Scalar >
void Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::scalarProdsImpl ( const MultiVectorBase< Scalar > &  X,
const MultiVectorBase< Scalar > &  Y,
const ArrayView< Scalar > &  scalarProds 
) const
virtual
template<class Scalar >
bool Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::isEuclidean ( ) const
virtual
template<class Scalar >
bool Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::hasInCoreView ( const Range1D rng,
const EViewType  viewType,
const EStrideType  strideType 
) const
virtual
template<class Scalar >
Teuchos::RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::clone ( ) const
virtual
template<class Scalar >
int Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::numBlocks ( ) const
virtual
template<class Scalar >
Teuchos::RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::getBlock ( const int  k) const
virtual
template<class Scalar >
Teuchos::RCP< VectorBase< Scalar > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::createMember ( ) const
protectedvirtual
template<class Scalar >
Teuchos::RCP< MultiVectorBase< Scalar > > Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >::createMembers ( int  numMembers) const
protectedvirtual

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