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

Base interface for product vectors. More...

#include <Thyra_ProductVectorBase.hpp>

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

Public Member Functions

virtual RCP< VectorBase< Scalar > > getNonconstVectorBlock (const int k)=0
 Returns a non-persisting non-const view of the (zero-based) kth block vector. More...
 
virtual RCP< const VectorBase
< Scalar > > 
getVectorBlock (const int k) const =0
 Returns a non-persisting const view of the (zero-based) kth block vector. More...
 
- Public Member Functions inherited from Thyra::VectorBase< Scalar >
void assign (const VectorBase< Scalar > &x)
 Vector assignment: More...
 
void randomize (Scalar l, Scalar u)
 Random vector generation: More...
 
void update (Scalar alpha, const VectorBase< Scalar > &x)
 AXPY: More...
 
void linear_combination (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const VectorBase< Scalar > > > &x, const Scalar &beta)
 Linear combination: More...
 
Scalar dot (const VectorBase< Scalar > &x) const
 Euclidean dot product: result = x^H * this. More...
 
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_1 () const
 One (1) norm: result = ||v||1. More...
 
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_2 () const
 Euclidean (2) norm: result = ||v||2. More...
 
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_2 (const VectorBase< Scalar > &x) const
 Weighted Euclidean (2) norm: result = ||v||2. More...
 
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_inf () const
 Infinity norm: result = ||v||inf. More...
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
space () const =0
 Return a smart pointer to the vector space that this vector belongs to. More...
 
void applyOp (const RTOpPack::RTOpT< Scalar > &op, const ArrayView< const Ptr< const VectorBase< Scalar > > > &vecs, const ArrayView< const Ptr< VectorBase< Scalar > > > &targ_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal global_offset) const
 Calls applyOpImpl(). More...
 
virtual RCP< VectorBase< Scalar > > clone_v () const =0
 Returns a seprate cloned copy of *this vector with the same values but different storage. More...
 
void acquireDetachedView (const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
 Calls acquireDetachedVectorViewImpl(). More...
 
void releaseDetachedView (RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
 Calls releaseDetachedVectorViewImpl(). More...
 
void acquireDetachedView (const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
 Calls acquireNonconstDetachedVectorViewImpl(). More...
 
void commitDetachedView (RTOpPack::SubVectorView< Scalar > *sub_vec)
 Calls commitDetachedView(). More...
 
void setSubVector (const RTOpPack::SparseSubVectorT< Scalar > &sub_vec)
 Calls setSubVectorImpl(). More...
 
- Public Member Functions inherited from Thyra::MultiVectorBase< Scalar >
void assign (Scalar alpha)
 V = alpha. More...
 
void assign (const MultiVectorBase< Scalar > &mv)
 V = mv. More...
 
void scale (Scalar alpha)
 
void update (Scalar alpha, const MultiVectorBase< Scalar > &mv)
 
void linear_combination (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
 Y.col(j)(i) = beta*Y.col(j)(i) + sum( alpha[k]*X[k].col(j)(i), More...
 
void dots (const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
 Column-wise Euclidean dot product. More...
 
void norms_1 (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
 Column-wise 1-norms. More...
 
void norms_2 (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
 Column-wise 2-norms. More...
 
void norms_inf (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
 Column-wise infinity-norms. More...
 
RCP< const VectorBase< Scalar > > col (Ordinal j) const
 Calls colImpl(). More...
 
RCP< VectorBase< Scalar > > col (Ordinal j)
 Calls nonconstColImpl(). More...
 
RCP< const MultiVectorBase
< Scalar > > 
subView (const Range1D &colRng) const
 Calls contigSubViewImpl(). More...
 
RCP< MultiVectorBase< Scalar > > subView (const Range1D &colRng)
 Calls nonconstContigSubViewImpl(). More...
 
RCP< const MultiVectorBase
< Scalar > > 
subView (const ArrayView< const int > &cols) const
 nonContigSubViewImpl(). More...
 
RCP< MultiVectorBase< Scalar > > subView (const ArrayView< const int > &cols)
 nonconstNonContigSubViewImpl(). More...
 
void applyOp (const RTOpPack::RTOpT< Scalar > &primary_op, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &multi_vecs, const ArrayView< const Ptr< MultiVectorBase< Scalar > > > &targ_multi_vecs, const ArrayView< const Ptr< RTOpPack::ReductTarget > > &reduct_objs, const Ordinal primary_global_offset) const
 Calls mvMultiReductApplyOpImpl(). More...
 
void applyOp (const RTOpPack::RTOpT< Scalar > &primary_op, const RTOpPack::RTOpT< Scalar > &secondary_op, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &multi_vecs, const ArrayView< const Ptr< MultiVectorBase< Scalar > > > &targ_multi_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal primary_global_offset) const
 mvSingleReductApplyOpImpl(). More...
 
void acquireDetachedView (const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
 Calls acquireDetachedMultiVectorViewImpl(). More...
 
void releaseDetachedView (RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
 Calls releaseDetachedMultiVectorViewImpl(). More...
 
void acquireDetachedView (const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
 Calls acquireNonconstDetachedMultiVectorViewImpl(). More...
 
void commitDetachedView (RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
 Calls commitNonconstDetachedMultiVectorViewImpl(). More...
 
virtual RCP< MultiVectorBase
< Scalar > > 
clone_mv () const =0
 Clone the multi-vector object (if supported). More...
 
RCP< const LinearOpBase< Scalar > > clone () const
 This function is simply overridden to return this->clone_mv(). More...
 
- Public Member Functions inherited from Thyra::LinearOpBase< Scalar >
virtual RCP< const
VectorSpaceBase< Scalar > > 
range () const =0
 Return a smart pointer for the range space for this operator. More...
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
domain () const =0
 Return a smart pointer for the domain space for this operator. More...
 
bool opSupported (EOpTransp M_trans) const
 Return if the M_trans operation of apply() is supported or not. More...
 
void apply (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
 Apply the linear operator to a multi-vector : Y = alpha*op(M)*X + beta*Y. More...
 
- Public Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar >
bool rowStatIsSupported (const RowStatLinearOpBaseUtils::ERowStat rowStat) const
 Determine if a given row stat is supported. More...
 
void getRowStat (const RowStatLinearOpBaseUtils::ERowStat rowStat, const Ptr< VectorBase< Scalar > > &rowStatVec) const
 Get some statistics about a supported row. More...
 
- Public Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar >
bool supportsScaleLeft () const
 Determines if this objects supports left scaling. More...
 
bool supportsScaleRight () const
 Determines if this objects supports right scaling. More...
 
void scaleLeft (const VectorBase< Scalar > &row_scaling)
 Left scales operator with diagonal scaling operator. More...
 
void scaleRight (const VectorBase< Scalar > &col_scaling)
 Right scales operator with diagonal scaling operator. More...
 
- Public Member Functions inherited from Thyra::ProductMultiVectorBase< Scalar >
virtual Teuchos::RCP< const
ProductVectorSpaceBase< Scalar > > 
productSpace () const =0
 Returns the associated product vector space that represents the range. More...
 
virtual bool blockIsConst (const int k) const =0
 Return if the kth multi-vector block is const-only. More...
 
virtual Teuchos::RCP
< MultiVectorBase< Scalar > > 
getNonconstMultiVectorBlock (const int k)=0
 Returns a non-persisting non-const view of the zero-based kth block multi-vector. More...
 
virtual Teuchos::RCP< const
MultiVectorBase< Scalar > > 
getMultiVectorBlock (const int k) const =0
 Returns a non-persisting const view of the (zero-based) kth block multi-vector. More...
 

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP< Thyra::ProductVectorBase
< Scalar > > 
nonconstProductVectorBase (const RCP< Thyra::VectorBase< Scalar > > &v)
 Dynamic cast from a VectorBase to a ProductVectorBase object and thow exception if this fails. More...
 
template<class Scalar >
RCP< const
Thyra::ProductVectorBase
< Scalar > > 
productVectorBase (const RCP< const Thyra::VectorBase< Scalar > > &v)
 Dynamic cast from a const VectorBase to a const ProductVectorBase object and thow exception if this fails. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Thyra::VectorBase< Scalar >
virtual void assignVecImpl (const VectorBase< Scalar > &x)=0
 Virtual implementation for NVI assign. More...
 
virtual void randomizeImpl (Scalar l, Scalar u)=0
 Virtual implementation for NVI randomize. More...
 
virtual void absImpl (const VectorBase< Scalar > &x)=0
 Virtual implementation for NVI abs. More...
 
virtual void reciprocalImpl (const VectorBase< Scalar > &x)=0
 Virtual implementation for NVI reciprocal. More...
 
virtual void eleWiseScaleImpl (const VectorBase< Scalar > &x)=0
 Virtual implementation for NVI ele_wise_scale. More...
 
virtual void updateVecImpl (Scalar alpha, const VectorBase< Scalar > &x)=0
 Virtual implementation for NVI update. More...
 
virtual void linearCombinationVecImpl (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const VectorBase< Scalar > > > &x, const Scalar &beta)=0
 Virtual implementation for NVI linear_combination. More...
 
virtual Scalar dotImpl (const VectorBase< Scalar > &x) const =0
 Virtual implementation for NVI dot. More...
 
virtual Teuchos::ScalarTraits
< Scalar >::magnitudeType 
norm1Impl () const =0
 Virtual implementation for NVI norm_1. More...
 
virtual Teuchos::ScalarTraits
< Scalar >::magnitudeType 
norm2Impl () const =0
 Virtual implementation for NVI norm_2. More...
 
virtual Teuchos::ScalarTraits
< Scalar >::magnitudeType 
norm2WeightedImpl (const VectorBase< Scalar > &x) const =0
 Virtual implementation for NVI norm_2 (weighted). More...
 
virtual Teuchos::ScalarTraits
< Scalar >::magnitudeType 
normInfImpl () const =0
 Virtual implementation for NVI norm_inf. More...
 
virtual void applyOpImpl (const RTOpPack::RTOpT< Scalar > &op, const ArrayView< const Ptr< const VectorBase< Scalar > > > &vecs, const ArrayView< const Ptr< VectorBase< Scalar > > > &targ_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal global_offset) const =0
 Apply a reduction/transformation operator over a set of vectors: op(op(v[0]...v[nv-1],z[0]...z[nz-1]),(*reduct_obj)) -> z[0]...z[nz-1],(*reduct_obj). More...
 
virtual void acquireDetachedVectorViewImpl (const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const =0
 Get a non-mutable explicit view of a sub-vector. More...
 
virtual void releaseDetachedVectorViewImpl (RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const =0
 Free an explicit view of a sub-vector. More...
 
virtual void acquireNonconstDetachedVectorViewImpl (const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)=0
 Get a mutable explicit view of a sub-vector. More...
 
virtual void commitNonconstDetachedVectorViewImpl (RTOpPack::SubVectorView< Scalar > *sub_vec)=0
 Commit changes for a mutable explicit view of a sub-vector. More...
 
virtual void setSubVectorImpl (const RTOpPack::SparseSubVectorT< Scalar > &sub_vec)=0
 Set a specific sub-vector. More...
 
- Protected Member Functions inherited from Thyra::MultiVectorBase< Scalar >
void absRowSum (const Teuchos::Ptr< Thyra::VectorBase< Scalar > > &output) const
 
void absColSum (const Teuchos::Ptr< Thyra::VectorBase< Scalar > > &output) const
 
virtual void assignImpl (Scalar alpha)=0
 Virtual implementation for NVI assign(Scalar). More...
 
virtual void assignMultiVecImpl (const MultiVectorBase< Scalar > &mv)=0
 Virtual implementation for NVI assign(MV). More...
 
virtual void scaleImpl (Scalar alpha)=0
 Virtual implementation for NVI scale(). More...
 
virtual void updateImpl (Scalar alpha, const MultiVectorBase< Scalar > &mv)=0
 Virtual implementation for NVI update(). More...
 
virtual void linearCombinationImpl (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)=0
 Virtual implementation for NVI linear_combination(). More...
 
virtual void dotsImpl (const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const =0
 Virtual implementation for NVI dots(). More...
 
virtual void norms1Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const =0
 Virtual implementation for NVI norms_1(). More...
 
virtual void norms2Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const =0
 Virtual implementation for NVI norms_2(). More...
 
virtual void normsInfImpl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const =0
 Virtual implementation for NVI norms_inf(). More...
 
virtual RCP< const VectorBase
< Scalar > > 
colImpl (Ordinal j) const
 Return a non-changeable view of a constituent column vector. More...
 
virtual RCP< VectorBase< Scalar > > nonconstColImpl (Ordinal j)=0
 Return a changeable view of a constituent column vector. More...
 
virtual RCP< const
MultiVectorBase< Scalar > > 
contigSubViewImpl (const Range1D &colRng) const =0
 Return a non-changeable sub-view of a contiguous set of columns of the this multi-vector. More...
 
virtual RCP< MultiVectorBase
< Scalar > > 
nonconstContigSubViewImpl (const Range1D &colRng)=0
 Return a changeable sub-view of a contiguous set of columns of the this multi-vector. More...
 
virtual RCP< const
MultiVectorBase< Scalar > > 
nonContigSubViewImpl (const ArrayView< const int > &cols) const =0
 Return a non-changeable sub-view of a non-contiguous set of columns of this multi-vector. More...
 
virtual RCP< MultiVectorBase
< Scalar > > 
nonconstNonContigSubViewImpl (const ArrayView< const int > &cols)=0
 Return a changeable sub-view of a non-contiguous set of columns of this multi-vector. More...
 
virtual void mvMultiReductApplyOpImpl (const RTOpPack::RTOpT< Scalar > &primary_op, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &multi_vecs, const ArrayView< const Ptr< MultiVectorBase< Scalar > > > &targ_multi_vecs, const ArrayView< const Ptr< RTOpPack::ReductTarget > > &reduct_objs, const Ordinal primary_global_offset) const =0
 Apply a reduction/transformation operator column by column and return an array of the reduction objects. More...
 
virtual void mvSingleReductApplyOpImpl (const RTOpPack::RTOpT< Scalar > &primary_op, const RTOpPack::RTOpT< Scalar > &secondary_op, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &multi_vecs, const ArrayView< const Ptr< MultiVectorBase< Scalar > > > &targ_multi_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal primary_global_offset) const =0
 Apply a reduction/transformation operator column by column and reduce the intermediate reduction objects into a single reduction object. More...
 
virtual void acquireDetachedMultiVectorViewImpl (const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const =0
 Get a non-changeable explicit view of a sub-multi-vector. More...
 
virtual void releaseDetachedMultiVectorViewImpl (RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const =0
 Free a non-changeable explicit view of a sub-multi-vector. More...
 
virtual void acquireNonconstDetachedMultiVectorViewImpl (const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)=0
 Get a changeable explicit view of a sub-multi-vector. More...
 
virtual void commitNonconstDetachedMultiVectorViewImpl (RTOpPack::SubMultiVectorView< Scalar > *sub_mv)=0
 Commit changes for a changeable explicit view of a sub-multi-vector. More...
 
virtual bool rowStatIsSupportedImpl (const RowStatLinearOpBaseUtils::ERowStat rowStat) const
 
virtual void getRowStatImpl (const RowStatLinearOpBaseUtils::ERowStat rowStat, const Ptr< VectorBase< Scalar > > &rowStatVec) const
 
virtual bool supportsScaleLeftImpl () const
 
virtual bool supportsScaleRightImpl () const
 
virtual void scaleLeftImpl (const VectorBase< Scalar > &row_scaling)
 
virtual void scaleRightImpl (const VectorBase< Scalar > &col_scaling)
 
- Protected Member Functions inherited from Thyra::LinearOpBase< Scalar >
virtual bool opSupportedImpl (EOpTransp M_trans) const =0
 Override in subclass. More...
 
virtual void applyImpl (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const =0
 Override in subclass. More...
 
- Protected Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar >
- Protected Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar >

Detailed Description

template<class Scalar>
class Thyra::ProductVectorBase< Scalar >

Base interface for product vectors.

This class defines an abstract interface for a vector that is built out of the one or more other vectors to form what mathematicians like to call a "product vector".

A product vector is simply the concatenation of two or more vectors to form a larger "composite" vector. Specifically, a product vector with numBlock constituent block vectors represents the blocked vector

        [ v[0]           ]
        [ v[1]           ]
this =  [ .              ]
        [ v[numBlocks-1] ]

The constituent vectors v[k] can be accessed through the const and non-const access functions getBlock().

A product vector knows its product space which is returned by the productSpace() function. A ProductVectorBase object is created by a ProductVectorSpaceBase object and never directly created by clients.

This class is only an interface. A standard implementation of this interface that should be sufficient for 99% or so of use cases is provided in the concrete subclass DefaultProductVector.

Definition at line 87 of file Thyra_ProductVectorBase.hpp.

Member Function Documentation

template<class Scalar>
virtual RCP<VectorBase<Scalar> > Thyra::ProductVectorBase< Scalar >::getNonconstVectorBlock ( const int  k)
pure virtual

Returns a non-persisting non-const view of the (zero-based) kth block vector.

Parameters
k[in] The (zero-based) kth block index specifying which vector block to access.

Preconditions:

Note that *this is not guaranteed to be modified until the smart pointer returned from this function, as well as any other smart pointers created from this smart pointer, are destroyed. This requirement allows more flexibility in how this function is implemented.

Also note that no further interactions with *this should be performed until the view returned from this function is released as described above.

Implemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, and Thyra::DefaultClusteredSpmdProductVector< Scalar >.

template<class Scalar>
virtual RCP<const VectorBase<Scalar> > Thyra::ProductVectorBase< Scalar >::getVectorBlock ( const int  k) const
pure virtual

Returns a non-persisting const view of the (zero-based) kth block vector.

Parameters
k[in] The (zero-based) kth block index specifying which vectorblock to access.

Preconditions:

Implemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, and Thyra::DefaultClusteredSpmdProductVector< Scalar >.

Friends And Related Function Documentation

template<class Scalar >
RCP< Thyra::ProductVectorBase< Scalar > > nonconstProductVectorBase ( const RCP< Thyra::VectorBase< Scalar > > &  v)
related

Dynamic cast from a VectorBase to a ProductVectorBase object and thow exception if this fails.

Definition at line 148 of file Thyra_ProductVectorBase.hpp.

template<class Scalar >
RCP< const Thyra::ProductVectorBase< Scalar > > productVectorBase ( const RCP< const Thyra::VectorBase< Scalar > > &  v)
related

Dynamic cast from a const VectorBase to a const ProductVectorBase object and thow exception if this fails.

Definition at line 164 of file Thyra_ProductVectorBase.hpp.


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