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::VectorSpaceBase< Scalar > Class Template Referenceabstract

Abstract interface for objects that represent a space for vectors. More...

#include <Thyra_VectorSpaceBase_decl.hpp>

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

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP< const VectorSpaceBase
< Scalar > > 
makeHaveOwnership (const RCP< const VectorSpaceBase< Scalar > > &vs)
 Helper function that clones a VectorSpaceBase object if the RCP does not have ownership. More...
 
template<class Scalar >
RCP< VectorBase< Scalar > > createMember (const RCP< const VectorSpaceBase< Scalar > > &vs, const std::string &label="")
 Create a vector member from the vector space. More...
 
template<class Scalar >
RCP< VectorBase< Scalar > > createMember (const VectorSpaceBase< Scalar > &vs, const std::string &label="")
 Calls createMember(Teuchos::rcp(&vs,false)). More...
 
template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembers (const RCP< const VectorSpaceBase< Scalar > > &vs, int numMembers, const std::string &label="")
 Create a set of vector members (a MultiVectorBase) from the vector space. More...
 
template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembers (const RCP< const VectorSpaceBase< Scalar > > &vs, const RCP< const VectorSpaceBase< Scalar > > &domain, const std::string &label="")
 Create a set of vector members (a MultiVectorBase) from the vector space. More...
 
template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembers (const VectorSpaceBase< Scalar > &vs, int numMembers, const std::string &label="")
 Calls createMembers(Teuchos::rcp(&vs,false),numMembers). More...
 
template<class Scalar >
RCP< VectorBase< Scalar > > createMemberView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::SubVectorView< Scalar > &raw_v, const std::string &label="")
 Create a vector member that is a non-const view of raw data. More...
 
template<class Scalar >
RCP< VectorBase< Scalar > > createMemberView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::SubVectorView< Scalar > &raw_v, const std::string &label="")
 Calls createMemberView(Teuchos::rcp(&vs,false),raw_v). More...
 
template<class Scalar >
RCP< const VectorBase< Scalar > > createMemberView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::ConstSubVectorView< Scalar > &raw_v, const std::string &label="")
 Create a vector member that is a const view of raw data. More...
 
template<class Scalar >
RCP< const VectorBase< Scalar > > createMemberView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::ConstSubVectorView< Scalar > &raw_v, const std::string &label="")
 Calls createMemberView(Teuchos::rcp(&vs,false),raw_v). More...
 
template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembersView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::SubMultiVectorView< Scalar > &raw_mv, const std::string &label="")
 Create a multi-vector member that is a non-const view of raw data. More...
 
template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembersView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::SubMultiVectorView< Scalar > &raw_mv, const std::string &label="")
 Calls createMembersView(Teuchos::rcp(&vs,false),raw_mv). More...
 
template<class Scalar >
RCP< const MultiVectorBase
< Scalar > > 
createMembersView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv, const std::string &label="")
 Create a multi-vector member that is a const view of raw data. More...
 
template<class Scalar >
RCP< const MultiVectorBase
< Scalar > > 
createMembersView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv, const std::string &label="")
 Calls createMembersView(Teuchos::rcp(&vs,false),raw_mv). More...
 

Public pure virtual functions that must be overridden

virtual Ordinal dim () const =0
 Return the dimension of the vector space. More...
 
virtual bool isCompatible (const VectorSpaceBase< Scalar > &vecSpc) const =0
 Compare the compatibility of two vector spaces. More...
 
virtual RCP< const
VectorSpaceFactoryBase< Scalar > > 
smallVecSpcFcty () const =0
 Return a VectorSpaceFactoryBase object for the creation of (usually serial) vector spaces with a small dimension. More...
 
virtual Scalar scalarProd (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y) const =0
 Return the scalar product of two vectors in the vector space. More...
 
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 virtual functions with default implementations

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...
 
virtual bool hasInCoreView (const Range1D &rng=Range1D(), const EViewType viewType=VIEW_TYPE_DETACHED, const EStrideType strideType=STRIDE_TYPE_NONUNIT) const
 Returns true if this->acquireDetachedView(rng,...) returns a direct view of the range of data requested. More...
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
clone () const
 Clone this object (if supported). More...
 

Protected pure virtual functions that must be overridden

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...
 
virtual RCP< VectorBase< Scalar > > createMember () const =0
 Create a vector member from the vector space. More...
 
virtual RCP< MultiVectorBase
< Scalar > > 
createMembers (int numMembers) const =0
 Create a set of vector members (a MultiVectorBase) from the vector space. More...
 
virtual RCP< VectorBase< Scalar > > createMemberView (const RTOpPack::SubVectorView< Scalar > &raw_v) const =0
 Create a vector member that is a non-const view of raw vector data. More...
 
virtual RCP< const VectorBase
< Scalar > > 
createMemberView (const RTOpPack::ConstSubVectorView< Scalar > &raw_v) const =0
 Create a vector member that is a const view of raw vector data. More...
 
virtual RCP< MultiVectorBase
< Scalar > > 
createMembersView (const RTOpPack::SubMultiVectorView< Scalar > &raw_mv) const =0
 Create a multi-vector member that is a non-const view of raw multi-vector data. More...
 
virtual RCP< const
MultiVectorBase< Scalar > > 
createMembersView (const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv) const =0
 Create a multi-vector member that is a const view of raw multi-vector data. More...
 

Protected virtual funtions.

virtual void scalarProdsImpl (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds) const =0
 

Detailed Description

template<class Scalar>
class Thyra::VectorSpaceBase< Scalar >

Abstract interface for objects that represent a space for vectors.

This interface acts primarily as an "Abstract Factory" interface for creating VectorBase objects using the non-member function Thyra::createMember(). A VectorSpaceBase can also create MultiVectorBase objects which represent a compact collection of vectors using the non-member function Thyra::createMembers(). A secondary role for VectorSpaceBase objects is to test for compatibility of vector space objects using the isCompatible() method and to apply the space's scalar (inner) product.

Clients can not directly create VectorBase and MultiVectorBase objects using the member functions createMember(), createMembers(), createMemberView(), and createMembersView() but instead must use the non-member Thyra_Op_Vec_createMember_grp.

Note that a VectorSpaceBase object must also be able to create MultiVectorBase objects with any number of column vectors using the Thyra::createMembers() function. Once created, the LinearOpBase::domain() function supported by a created MultiVectorBase object returns a vector space of dimension equal to the number of columns in the multi-vector. An interesting implication of this design is that the creation of a multi-vector provides a way for clients to create vector spaces of any arbitrary (although small usually) dimension. In order to give the client the same ability to create smaller vector spaces without having to create a dummy multi-vector object first, the method smallVecSpcFcty() is included. The method smallVecSpcFcty() returns a VectorSpaceFactoryBase object that can create (typically serial) vector spaces of any small dimension that are compatible with the domain spaces of MultiVectorBase objects created by the vector space object.

A vector space is also where the scalar product for the space is defined which is computed by the scalarProd() method. A scalar product allows the vector space to introduce scaling into many different types of numerical algorithms.

If the underlying object is not initialized, then dim()==0 will be true and none of the other methods should be called or exceptions will be thrown.

Notes for subclass developers

This is a fairly bare-bones interface class without much in the way of default function implementations. The subclass VectorSpaceDefaultBase provides a default multi-vector implementation and should the the first choice for subclass implementations.

Definition at line 367 of file Thyra_OperatorVectorTypes.hpp.

Member Function Documentation

template<class Scalar>
virtual Ordinal Thyra::VectorSpaceBase< Scalar >::dim ( ) const
pure virtual

Return the dimension of the vector space.

If the underlying object is not initialized, then dim()==0 will be true and none of the other methods should be called.

Implemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, Thyra::SpmdVectorSpaceDefaultBase< Scalar >, and Thyra::DefaultMultiVectorProductVectorSpace< Scalar >.

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

Compare the compatibility of two vector spaces.

If this function returns true, then vectors created from either of the vector spaces will be compatible and can be combined in vector operations.

Preconditions:

Postconditions:

  • [this->dim() != vecSpc.dim()] returnVal == false

Invariants:

  • [this->isCompatible(vecSpc) == true] vecSpc.isCompatible(*this) == true

Implemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::SpmdVectorSpaceDefaultBase< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, and Thyra::DefaultMultiVectorProductVectorSpace< Scalar >.

template<class Scalar>
virtual RCP< const VectorSpaceFactoryBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::smallVecSpcFcty ( ) const
pure virtual
template<class Scalar>
virtual Scalar Thyra::VectorSpaceBase< Scalar >::scalarProd ( const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
) const
pure virtual

Return the scalar product of two vectors in the vector space.

Preconditions:

Implemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, Thyra::ScalarProdVectorSpaceBase< Scalar >, and Thyra::DefaultMultiVectorProductVectorSpace< Scalar >.

template<class Scalar>
void Thyra::VectorSpaceBase< Scalar >::scalarProds ( const MultiVectorBase< Scalar > &  X,
const MultiVectorBase< Scalar > &  Y,
const ArrayView< Scalar > &  scalarProds_out 
) const
inline

Return the scalar product of each column in two multi-vectors in the vector space.

Parameters
X[in] Multi-vector.
Y[in] Multi-vector.
scalarProds_out[out] Array (length X.domain()->dim()) containing the scalar products scalarProds_out[j] = this->scalarProds_out(*X.col(j),*Y.col(j)), for j = 0 ... X.domain()->dim()-1.

Preconditions:

Definition at line 387 of file Thyra_VectorSpaceBase_decl.hpp.

template<class Scalar >
bool Thyra::VectorSpaceBase< Scalar >::isEuclidean ( ) const
virtual

Return if this vector space has a Euclidean (identity) basis in which case the scalar product is the same as the dot product.

The default implementation returns false (even though on average the Euclidean scalar product is used).

Reimplemented in Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, and Thyra::ScalarProdVectorSpaceBase< Scalar >.

Definition at line 69 of file Thyra_VectorSpaceBase_def.hpp.

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

Returns true if this->acquireDetachedView(rng,...) returns a direct view of the range of data requested.

Parameters
rng[in] The range of elements for the view (see acquireDetachedView()). The default value is Range1D() (i.e. All of the elements in the vector).
viewType[in] The type of view allowed.
strideType[in] The type of stride the view is allowed to be.

Preconditions:

There are three different questions about the behavior of the acquireDetachedView() that this query function can answer:

  • The elements in rng are fairly cheaply accessble in local (i.e. in-core) memory if this->hasInCoreView(rng)==true. Note that this also allows for detached temporary copies of data.

  • A direct view of the elements in rng is available in local (i.e. in-core) memory if this->hasInCoreView(rng,VIEW_TYPE_DIRECT)==true. No copy of data is allowed here.

  • A direct view of the elements in rng with unit stride is available in local (i.e. in-core) memory if this->hasInCoreView(rng,VIEW_TYPE_DIRECT,STRIDE_TYPE_UNIT)==true No copy of data is allowed here.

The default implementation returns false (i.e. by default we do not assume that any direct and/or contiguous views of any range of data are provided).

Reimplemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, Thyra::DefaultSpmdVectorSpace< Scalar >, Thyra::DefaultMultiVectorProductVectorSpace< Scalar >, and Thyra::TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 76 of file Thyra_VectorSpaceBase_def.hpp.

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

Clone this object (if supported).

It is allowed for returnVal.get()==NULL which means that this capability is optional.

Preconditions:

The default implementation returns returnVal.get()==NULL.

Reimplemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, Thyra::DefaultMultiVectorProductVectorSpace< Scalar >, Thyra::DefaultSpmdVectorSpace< Scalar >, and Thyra::TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 85 of file Thyra_VectorSpaceBase_def.hpp.

template<class Scalar>
virtual RCP< VectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createMember ( ) const
protectedpure virtual

Create a vector member from the vector space.

Preconditions:

Postconditions:

  • returnVal.get() != NULL

  • returnVal->space()->isCompatible(*this) == true

Note: This function is not to be called directly since it is protected! See the Thyra_Op_Vec_createMember_grp.

Returns
A smart reference counted pointer to a dynamically allocated vector object. After construction the values in the vector *returnVal are unspecified (uninitialized). This allows for faster execution times. Note that returnVal->space().get() == this need not be true.

Implemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, Thyra::DefaultSpmdVectorSpace< Scalar >, Thyra::DefaultMultiVectorProductVectorSpace< Scalar >, and Thyra::TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

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

Create a set of vector members (a MultiVectorBase) from the vector space.

Preconditions:

  • this->dim() > 0

  • num_vecs >= 1

Postconditions:

  • returnVal->range()->isCompatible(*this) == true

  • returnVal->domain()->dim() == numMembers

Returns
A smart reference-counted pointer to a dynamically allocated multi-vector object. After construction, the values in *returnVal are unspecified (uninitialized). This allows for faster execution times. Note that returnVal->range().get()==this does not have to be true but will be in may cases.

Implemented in Thyra::DefaultProductVectorSpace< Scalar >, Thyra::DefaultClusteredSpmdProductVectorSpace< Scalar >, Thyra::DefaultSpmdVectorSpace< Scalar >, Thyra::DefaultMultiVectorProductVectorSpace< Scalar >, Thyra::TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Thyra::VectorSpaceDefaultBase< Scalar >.

template<class Scalar>
virtual RCP<VectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createMemberView ( const RTOpPack::SubVectorView< Scalar > &  raw_v) const
protectedpure virtual

Create a vector member that is a non-const view of raw vector data.

Parameters
raw_v[in] On input contains pointer (i.e. raw_v.values()) to array that the returned VectorBase will be a view of. The data pointed to by raw_v.values() must remain valid until the returned VectorBase object is destroyed.

Preconditions:

  • raw_v has been initialized to memory (i.e. raw_v.subDim()!=0 && raw_v.values()!=NULL).

  • raw_v is consistent with the local storage of this vector spaces vector data. This precondition is purposefully vague since this function can be used an variety of specialized use-cases.

Postconditions:

It is stated here that the client can not expect that the values pointed to by raw_v.values() to be changed until the smart pointer returnVal goes out of scope. This is to allow an implementation that temporarily copies data into and out of a VectorBase object using explicit vector access.

Implemented in Thyra::DefaultSpmdVectorSpace< Scalar >, and Thyra::VectorSpaceDefaultBase< Scalar >.

template<class Scalar>
virtual RCP<const VectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createMemberView ( const RTOpPack::ConstSubVectorView< Scalar > &  raw_v) const
protectedpure virtual

Create a vector member that is a const view of raw vector data.

Parameters
raw_v[in] On input contains pointer (i.e. raw_v.values()) to array that the returned VectorBase will be a view of. The data pointed to by raw_v.values() must remain valid until the returned VectorBase object is destroyed.

This function works exactly the same as the previous version that takes a RTOpPack::SubVectorView object except that this version takes a RTOpPack::ConstSubVectorView and returns a smart pointer to a const VectorBase object.

Preconditions:

Postconditions:

Implemented in Thyra::DefaultSpmdVectorSpace< Scalar >, and Thyra::VectorSpaceDefaultBase< Scalar >.

template<class Scalar>
virtual RCP<MultiVectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createMembersView ( const RTOpPack::SubMultiVectorView< Scalar > &  raw_mv) const
protectedpure virtual

Create a multi-vector member that is a non-const view of raw multi-vector data.

Parameters
raw_mv[in] On input contains pointer (i.e. raw_mv.values()) to array that the returned MultiVectorBase will be a view of.

Preconditions:

  • raw_mv has been initialized to memory (i.e. raw_mv.subDim()!=0 && raw_mv.values()!=NULL).

  • raw_mv is consistent with the local storage of this spaces vector data. This precondition is purposefully vague since this function can be used an variety of specialized use-cases.

Postconditions:

It is stated here that the client can not expect that the values pointed to by raw_mv.values() to be changed until the smart pointer returnVal goes out of scope. This is to allow for an implementation that temporarily copies data into and out of a MultiVectorBase object using explicit vector access.

Implemented in Thyra::DefaultSpmdVectorSpace< Scalar >, and Thyra::VectorSpaceDefaultBase< Scalar >.

template<class Scalar>
virtual RCP<const MultiVectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createMembersView ( const RTOpPack::ConstSubMultiVectorView< Scalar > &  raw_mv) const
protectedpure virtual

Create a multi-vector member that is a const view of raw multi-vector data.

Parameters
raw_mv[in] On input contains pointer (i.e. raw_mv.values()) to array that the returned MultiVectorBase will be a view of. The data pointed to by raw_mv.values() must remain valid until the returned MultiVectorBase object is destroyed.

This function works exactly the same as the previous version that takes a RTOpPack::SubMultiVectorView object except that this version takes a RTOpPack::ConstSubMultiVectorView object and returns a smart pointer to a const MultiVectorBase object.

Preconditions:

Postconditions:

Implemented in Thyra::DefaultSpmdVectorSpace< Scalar >, and Thyra::VectorSpaceDefaultBase< Scalar >.

template<class Scalar>
virtual RCP<MultiVectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createCachedMembersView ( const RTOpPack::SubMultiVectorView< Scalar > &  raw_mv) const
inlinevirtual

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.

Parameters
raw_mv[in] On input contains pointer (i.e. raw_mv.values()) to array that the returned MultiVectorBase will be a view of.

Preconditions:

  • raw_mv has been initialized to memory (i.e. raw_mv.subDim()!=0 && raw_mv.values()!=NULL).

  • raw_mv is consistent with the local storage of this spaces vector data. This precondition is purposefully vague since this function can be used an variety of specialized use-cases.

Postconditions:

It is stated here that the client can not expect that the values pointed to by raw_mv.values() to be changed until the smart pointer returnVal goes out of scope. This is to allow for an implementation that temporarily copies data into and out of a MultiVectorBase object using explicit vector access.

Reimplemented in Thyra::TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 748 of file Thyra_VectorSpaceBase_decl.hpp.

template<class Scalar>
virtual RCP<const MultiVectorBase<Scalar> > Thyra::VectorSpaceBase< Scalar >::createCachedMembersView ( const RTOpPack::ConstSubMultiVectorView< Scalar > &  raw_mv) const
inlinevirtual

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.

Parameters
raw_mv[in] On input contains pointer (i.e. raw_mv.values()) to array that the returned MultiVectorBase will be a view of. The data pointed to by raw_mv.values() must remain valid until the returned MultiVectorBase object is destroyed.

This function works exactly the same as the previous version that takes a RTOpPack::SubMultiVectorView object except that this version takes a RTOpPack::ConstSubMultiVectorView object and returns a smart pointer to a const MultiVectorBase object.

Preconditions:

Postconditions:

Reimplemented in Thyra::TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 779 of file Thyra_VectorSpaceBase_decl.hpp.

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

Friends And Related Function Documentation

template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > makeHaveOwnership ( const RCP< const VectorSpaceBase< Scalar > > &  vs)
related

Helper function that clones a VectorSpaceBase object if the RCP does not have ownership.

template<class Scalar >
RCP< VectorBase< Scalar > > createMember ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
const std::string &  label = "" 
)
related

Create a vector member from the vector space.

Calls VectorSpaceBase::createMember() on vs but the returned VectorBase object can live past vs.

template<class Scalar >
RCP< VectorBase< Scalar > > createMember ( const VectorSpaceBase< Scalar > &  vs,
const std::string &  label = "" 
)
related

Calls createMember(Teuchos::rcp(&vs,false)).

template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembers ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
int  numMembers,
const std::string &  label = "" 
)
related

Create a set of vector members (a MultiVectorBase) from the vector space.

Calls VectorSpaceBase::createMembers() on vs but the returned MultiVectorBase object can live past vs.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembers ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
const RCP< const VectorSpaceBase< Scalar > > &  domain,
const std::string &  label = "" 
)
related

Create a set of vector members (a MultiVectorBase) from the vector space.

Calls VectorSpaceBase::createMembers() on vs but the returned MultiVectorBase object can live past vs.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembers ( const VectorSpaceBase< Scalar > &  vs,
int  numMembers,
const std::string &  label = "" 
)
related

Calls createMembers(Teuchos::rcp(&vs,false),numMembers).

template<class Scalar >
RCP< VectorBase< Scalar > > createMemberView ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
const RTOpPack::SubVectorView< Scalar > &  raw_v,
const std::string &  label = "" 
)
related

Create a vector member that is a non-const view of raw data.

Calls VectorSpaceBase::createMemberView() on vs but the returned VectorBase object can live past vs.

template<class Scalar >
RCP< VectorBase< Scalar > > createMemberView ( const VectorSpaceBase< Scalar > &  vs,
const RTOpPack::SubVectorView< Scalar > &  raw_v,
const std::string &  label = "" 
)
related

Calls createMemberView(Teuchos::rcp(&vs,false),raw_v).

template<class Scalar >
RCP< const VectorBase< Scalar > > createMemberView ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
const RTOpPack::ConstSubVectorView< Scalar > &  raw_v,
const std::string &  label = "" 
)
related

Create a vector member that is a const view of raw data.

Calls VectorSpaceBase::createMemberView() on vs but the returned VectorBase object can live past vs.

template<class Scalar >
RCP< const VectorBase< Scalar > > createMemberView ( const VectorSpaceBase< Scalar > &  vs,
const RTOpPack::ConstSubVectorView< Scalar > &  raw_v,
const std::string &  label = "" 
)
related

Calls createMemberView(Teuchos::rcp(&vs,false),raw_v).

template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembersView ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
const RTOpPack::SubMultiVectorView< Scalar > &  raw_mv,
const std::string &  label = "" 
)
related

Create a multi-vector member that is a non-const view of raw data.

Calls VectorSpaceBase::createMembersView() on vs but the returned MultiVectorBase object can live past vs.

template<class Scalar >
RCP< MultiVectorBase< Scalar > > createMembersView ( const VectorSpaceBase< Scalar > &  vs,
const RTOpPack::SubMultiVectorView< Scalar > &  raw_mv,
const std::string &  label = "" 
)
related

Calls createMembersView(Teuchos::rcp(&vs,false),raw_mv).

template<class Scalar >
RCP< const MultiVectorBase< Scalar > > createMembersView ( const RCP< const VectorSpaceBase< Scalar > > &  vs,
const RTOpPack::ConstSubMultiVectorView< Scalar > &  raw_mv,
const std::string &  label = "" 
)
related

Create a multi-vector member that is a const view of raw data.

Calls VectorSpaceBase::createMembersView() on vs but the returned MultiVectorBase object can live past vs.

template<class Scalar >
RCP< const MultiVectorBase< Scalar > > createMembersView ( const VectorSpaceBase< Scalar > &  vs,
const RTOpPack::ConstSubMultiVectorView< Scalar > &  raw_mv,
const std::string &  label = "" 
)
related

Calls createMembersView(Teuchos::rcp(&vs,false),raw_mv).


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