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::DefaultProductVectorSpace< Scalar > Class Template Reference

More...

#include <Thyra_DefaultProductMultiVector_decl.hpp>

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

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP< DefaultProductVectorSpace
< Scalar > > 
productVectorSpace ()
 Nonmember constructor that constructs to uninitialized. More...
 
template<class Scalar >
RCP< DefaultProductVectorSpace
< Scalar > > 
productVectorSpace (const ArrayView< RCP< const VectorSpaceBase< Scalar > > > &vecSpaces)
 Nonmember constructor that takes an array of vector spaces. More...
 
template<class Scalar >
RCP< DefaultProductVectorSpace
< Scalar > > 
productVectorSpace (const RCP< const VectorSpaceBase< Scalar > > &vecSpace, const int numBlocks)
 Nonmember constructor that duplicates a block vector space numBlock times to form a product space. More...
 

Constructors/initializers/accessors

 DefaultProductVectorSpace ()
 Default construct to uninitialized. More...
 
 DefaultProductVectorSpace (const ArrayView< const RCP< const VectorSpaceBase< Scalar > > > &vecSpaces)
 
virtual void initialize (const ArrayView< const RCP< const VectorSpaceBase< Scalar > > > &vecSpaces)
 Initialize with a list of constituent vector spaces. More...
 
bool hasBeenCloned () const
 Return if this vector space was cloned. More...
 
virtual void uninitialize (const ArrayView< RCP< const VectorSpaceBase< Scalar > > > &vecSpaces=Teuchos::null)
 Uninitialize. More...
 
virtual const RCP< const
VectorSpaceBase< Scalar > > * 
vecSpaces () const
 Returns a pointer to an array (of length this->numBlocks()) to the constituent vector spaces. More...
 
virtual const OrdinalvecSpacesOffsets () const
 Returns a pointer to an array (of length this->numBlocks()+1) of offset into each constituent vector space. More...
 
void getVecSpcPoss (Ordinal i, int *kth_vector_space, Ordinal *kth_global_offset) const
 Get the position of the vector space object and its offset into a composite vector that owns the ith index in the composite vector. More...
 

Overridden from DefaultProductVectorSpace

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

Overridden from VectorSpaceBase

Ordinal dim () const
 Returns the summation of the constituent vector spaces. More...
 
bool isCompatible (const VectorSpaceBase< Scalar > &vecSpc) const
 Returns true only if also a product vector space and all constituent vectors are compatible. More...
 
RCP< VectorBase< Scalar > > createMember () const
 Returns a DefaultProductVector object. More...
 
Scalar scalarProd (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y) const
 Returns the sum of the scalar products of the constituent vectors. More...
 
void scalarProdsImpl (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds) const
 Returns the sum of the scalar products of each of the columns of the constituent multi-vectors. More...
 
bool hasInCoreView (const Range1D &rng, const EViewType viewType, const EStrideType strideType) const
 Returns true if all of the constituent vector spaces return true. More...
 
RCP< const
VectorSpaceFactoryBase< Scalar > > 
smallVecSpcFcty () const
 Returns getBlock(0)->smallVecSpcFcty(). More...
 
RCP< MultiVectorBase< Scalar > > createMembers (int numMembers) const
 Returns a DefaultProductMultiVector object. More...
 
RCP< const VectorSpaceBase
< Scalar > > 
clone () const
 Clones the object as promised. More...
 

Overridden from Teuchos::Describable

std::string description () const
 Prints just the name DefaultProductVectorSpace along with the overall dimension and the number of blocks. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
 Prints the details about the constituent vector spaces. More...
 

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 bool isEuclidean () const
 Return if this vector space has a Euclidean (identity) basis in which case the scalar product is the same as the dot product. More...
 
- Protected Member Functions inherited from Thyra::VectorSpaceBase< Scalar >
- Protected Member Functions inherited from Thyra::VectorSpaceDefaultBase< Scalar >
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
 

Detailed Description

template<class Scalar>
class Thyra::DefaultProductVectorSpace< Scalar >

Standard concrete implementation of a product vector space.

This subclass allows VectorSpaceBase objects to be built out of one or more other vector space objects to form a product space. The specific type of vector created by this->createMember() is of type DefaultProductVector but the client need not ever know this or deal with this type directly. However, a client may want to rcp_dynamic_cast() to access the ProductVectorBase interface.

To demonstrate how to use this class, suppose one has p vector spaces V[k] for k = 0...p-1 and one wants to form a concatenated vector space Z containing all of these vector spaces stacked on top of each other to form:

     [ V[0]   ]
 Z = [ V[1]   ]
     [ .      ]
     [ V[p-1] ]

Such a vector space can be constructed out of an array of p constituent VectorSpaceBase objects as shown in the following function: x

void constructProductSpace(
const ArrayView<const RCP<const VectorSpaceBase<Scalar> > &Vs,
const Ptr<RCP<const VectorSpaceBase<Scalar> > > &Z
)
{
Z = productVectorSpace<Scalar>(Vs);
}

Or, a product space can be constructed out of p copies of the same VectorSpaceBase object as follows:

void constructProductSpace(
const RCP<const VectorSpaceBase<Scalar> > V, int p,
const Ptr<RCP<const VectorSpaceBase<Scalar> > > &Z
)
{
Array<RCP<const VectorSpaceBase<Scalar> > > vecSpaces;
for( int k = 0; k < p; ++k ) vecSpaces.push_back(V);
Z = productVectorSpace<Scalar>(vecSpaces());
}

Once a DefaultProductVectorSpace object is initialized, it can be used just like any other VectorSpaceBase object. The method createMember() will create DefaultProductVector objects containing vector members from the constituent vector spaces. The method createMembers() will create ProductMultiVector objects containing multi-vector members from the constituent vector spaces.

There are several methods that can be used by clients that need to work with the individual constituent vector spaces. The method numBlocks() give the number of constituent vector spaces while vecSpaces() returns a pointer to a copy of the array of the constituent vector spaces passed to initialize(). Some other useful utility methods are also defined. The method vecSpacesOffsets() returns a pointer to an array that gives the offset of each constituent vector in the overall composite product vector. For example, the (zero-based) kth vector space this->vecSpaces()[k] owns the element indexes this->vecSpacesOffsets()[k] to this->vecSpacesOffsets()[k+1]-1. Determining which constituent vector space owns a particular element index can be found out by calling getVecSpcPoss().

The default assignment operator is allowed since it has the correct semantics for shallow copy. The default copy constructor is also allowed but only performs a shallow copy of the constituent vector space objects. If you want to copy the constituent vector space objects also you need to use the clone() method.

Definition at line 58 of file Thyra_DefaultBlockedLinearOp_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::DefaultProductVectorSpace< Scalar >::DefaultProductVectorSpace ( )

Default construct to uninitialized.

Definition at line 61 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
Thyra::DefaultProductVectorSpace< Scalar >::DefaultProductVectorSpace ( const ArrayView< const RCP< const VectorSpaceBase< Scalar > > > &  vecSpaces)

Construct to an initialized state (calls initialize).

Definition at line 67 of file Thyra_DefaultProductVectorSpace_def.hpp.

Member Function Documentation

template<class Scalar >
void Thyra::DefaultProductVectorSpace< Scalar >::initialize ( const ArrayView< const RCP< const VectorSpaceBase< Scalar > > > &  vecSpaces)
virtual

Initialize with a list of constituent vector spaces.

Parameters
numBlocks[in] The number of constituent vector spaces.
vecSpaces[in] If vecSpaces!=NULL then vecSpaces must point to an array of length this->numBlocks and on output vecSpace[i] will be set to this->vecSpaces()[i] for i=0,..,this->numBlocks()-1.

Preconditions:

  • numBlocks > 0

  • vecSpaces != NULL

  • vecSpaces[i].get() != NULL, i=0,...,numBlocks

  • The vector space create by vecSpace[i]->smallVecSpcFcty()->createVecSpc(k) must be compatible. In other words,

     <tt>vecSpace[i]->smallVecSpcFcty()->createVecSpc(k)->isCompatible(
       *vecSpace[j]->smallVecSpcFcty()->createVecSpc(k)
       ) == true </tt> for all <tt>i=[0,numBlocks]</tt>, <tt>j=[0,numBlocks]</tt>
       and valid <tt>k > 0</tt>.
     This is required to insure that product multi-vectors can be created
     with constituent multi-vector blocks that have compatible <tt>domain()</tt>
     vector spaces.
    

Postconditions:

  • this->dim() == sum( vecSpaces[i]->dim(), i=0,...,numBlocks-1 )
  • this->numBlocks()==numBlocks
  • getBlock(i).get() == vecSpaces[i].get(), i=0,...,numBlocks-1
  • vecSpaces()[i].get() == vecSpaces[i].get(), i=0,...,numBlocks-1

Definition at line 77 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
bool Thyra::DefaultProductVectorSpace< Scalar >::hasBeenCloned ( ) const
inline

Return if this vector space was cloned.

If this function returns true then the client needs to be careful about how the constituent vector spaces returned from uninitialize() will be used.

Definition at line 455 of file Thyra_DefaultProductVectorSpace_decl.hpp.

template<class Scalar >
void Thyra::DefaultProductVectorSpace< Scalar >::uninitialize ( const ArrayView< RCP< const VectorSpaceBase< Scalar > > > &  vecSpaces = Teuchos::null)
virtual

Uninitialize.

Parameters
numBlocks[out] If numBlocks!=NULL then on output *numBlocks will be set to this->numBlocks().
vecSpaces[out] If vecSpaces!=NULL then vecSpaces must point to an array of length this->numBlocks and on output vecSpace[i] will be set to this->vecSpaces()[i] for i=0,..,this->numBlocks()-1.

Postconditions:

Warning! If this->hasBeenCloned()==true then the client had better not mess with the constituent vector spaces returned in vecSpaces[] since another DefaultProductVectorSpace object is still using them.

Definition at line 120 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
const RCP< const VectorSpaceBase< Scalar > > * Thyra::DefaultProductVectorSpace< Scalar >::vecSpaces ( ) const
inlinevirtual

Returns a pointer to an array (of length this->numBlocks()) to the constituent vector spaces.

Postconditions:

Definition at line 440 of file Thyra_DefaultProductVectorSpace_decl.hpp.

template<class Scalar >
const Ordinal * Thyra::DefaultProductVectorSpace< Scalar >::vecSpacesOffsets ( ) const
inlinevirtual

Returns a pointer to an array (of length this->numBlocks()+1) of offset into each constituent vector space.

Postconditions:

Definition at line 448 of file Thyra_DefaultProductVectorSpace_decl.hpp.

template<class Scalar >
void Thyra::DefaultProductVectorSpace< Scalar >::getVecSpcPoss ( Ordinal  i,
int *  kth_vector_space,
Ordinal kth_global_offset 
) const

Get the position of the vector space object and its offset into a composite vector that owns the ith index in the composite vector.

Parameters
i[in] The index offset of the element to find the vector space object for.
kth_vector_space[out] The index for this->vectorSpaces()[kth_vector_space] that owns the element i.
kth_global_offset[out] The global offset for this->vectorSpaces()[kth_vector_space] in the composite vector.

Preconditions:

Postconditions:

  • kth_global_offset == this->vecSpacesoffsets()[kth_vector-space]

  • kth_global_offset <= i <= kth_global_offset + this->vecSpaces()[kth_vector_space]->dim() - 1

Definition at line 134 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
int Thyra::DefaultProductVectorSpace< Scalar >::numBlocks ( ) const
virtual
template<class Scalar >
Teuchos::RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultProductVectorSpace< Scalar >::getBlock ( const int  k) const
virtual
template<class Scalar >
Ordinal Thyra::DefaultProductVectorSpace< Scalar >::dim ( ) const
virtual

Returns the summation of the constituent vector spaces.

Implements Thyra::VectorSpaceBase< Scalar >.

Definition at line 183 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
bool Thyra::DefaultProductVectorSpace< Scalar >::isCompatible ( const VectorSpaceBase< Scalar > &  vecSpc) const
virtual

Returns true only if also a product vector space and all constituent vectors are compatible.

Implements Thyra::VectorSpaceBase< Scalar >.

Definition at line 190 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
Teuchos::RCP< VectorBase< Scalar > > Thyra::DefaultProductVectorSpace< Scalar >::createMember ( ) const
virtual

Returns a DefaultProductVector object.

Implements Thyra::VectorSpaceBase< Scalar >.

Definition at line 228 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
Scalar Thyra::DefaultProductVectorSpace< Scalar >::scalarProd ( const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
) const
virtual

Returns the sum of the scalar products of the constituent vectors.

Implements Thyra::VectorSpaceBase< Scalar >.

Definition at line 235 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
void Thyra::DefaultProductVectorSpace< Scalar >::scalarProdsImpl ( const MultiVectorBase< Scalar > &  X,
const MultiVectorBase< Scalar > &  Y,
const ArrayView< Scalar > &  scalarProds 
) const
virtual

Returns the sum of the scalar products of each of the columns of the constituent multi-vectors.

Implements Thyra::VectorSpaceBase< Scalar >.

Definition at line 260 of file Thyra_DefaultProductVectorSpace_def.hpp.

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

Returns true if all of the constituent vector spaces return true.

Reimplemented from Thyra::VectorSpaceBase< Scalar >.

Definition at line 302 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
Teuchos::RCP< const VectorSpaceFactoryBase< Scalar > > Thyra::DefaultProductVectorSpace< Scalar >::smallVecSpcFcty ( ) const
virtual

Returns getBlock(0)->smallVecSpcFcty().

Implements Thyra::VectorSpaceBase< Scalar >.

Definition at line 339 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
Teuchos::RCP< MultiVectorBase< Scalar > > Thyra::DefaultProductVectorSpace< Scalar >::createMembers ( int  numMembers) const
virtual
template<class Scalar >
Teuchos::RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultProductVectorSpace< Scalar >::clone ( ) const
virtual

Clones the object as promised.

Reimplemented from Thyra::VectorSpaceBase< Scalar >.

Definition at line 358 of file Thyra_DefaultProductVectorSpace_def.hpp.

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

Prints just the name DefaultProductVectorSpace along with the overall dimension and the number of blocks.

Reimplemented from Teuchos::Describable.

Definition at line 378 of file Thyra_DefaultProductVectorSpace_def.hpp.

template<class Scalar >
void Thyra::DefaultProductVectorSpace< Scalar >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
virtual

Prints the details about the constituent vector spaces.

This function outputs different levels of detail based on the value passed in for verbLevel:

ToDo: Finish documentation!

Reimplemented from Teuchos::Describable.

Definition at line 391 of file Thyra_DefaultProductVectorSpace_def.hpp.

Friends And Related Function Documentation

template<class Scalar >
RCP< DefaultProductVectorSpace< Scalar > > productVectorSpace ( )
related

Nonmember constructor that constructs to uninitialized.

Definition at line 393 of file Thyra_DefaultProductVectorSpace_decl.hpp.

template<class Scalar >
RCP< DefaultProductVectorSpace< Scalar > > productVectorSpace ( const ArrayView< RCP< const VectorSpaceBase< Scalar > > > &  vecSpaces)
related

Nonmember constructor that takes an array of vector spaces.

Definition at line 406 of file Thyra_DefaultProductVectorSpace_decl.hpp.

template<class Scalar >
RCP< DefaultProductVectorSpace< Scalar > > productVectorSpace ( const RCP< const VectorSpaceBase< Scalar > > &  vecSpace,
const int  numBlocks 
)
related

Nonmember constructor that duplicates a block vector space numBlock times to form a product space.

Definition at line 422 of file Thyra_DefaultProductVectorSpace_decl.hpp.


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