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

Abstract interface for finite-dimensional dense vectors. More...

#include <Thyra_VectorBase.hpp>

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

Related Functions

(Note that these are not member functions.)

template<class Scalar >
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=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...
 
template<class Scalar >
Scalar sum (const VectorBase< Scalar > &v)
 Sum of vector elements: result = sum( v(i), i = 0...v.space()->dim()-1 ). More...
 
template<class Scalar >
Scalar scalarProd (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
 Scalar product result = <x,y>. More...
 
template<class Scalar >
Scalar inner (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
 Inner/Scalar product result = <x,y>. More...
 
template<class Scalar >
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm (const VectorBase< Scalar > &v)
 Natural norm: result = sqrt(<v,v>). More...
 
template<class Scalar >
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_1 (const VectorBase< Scalar > &v)
 One (1) norm: result = ||v||1. More...
 
template<class Scalar >
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_2 (const VectorBase< Scalar > &v)
 Euclidean (2) norm: result = ||v||2. More...
 
template<class Scalar >
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_2 (const VectorBase< Scalar > &w, const VectorBase< Scalar > &v)
 Weighted Euclidean (2) norm: result = sqrt( sum( w(i)*conj(v(i))*v(i)) ). More...
 
template<class Scalar >
Teuchos::ScalarTraits< Scalar >
::magnitudeType 
norm_inf (const VectorBase< Scalar > &v_rhs)
 Infinity norm: result = ||v||inf. More...
 
template<class Scalar >
Scalar dot (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
 Dot product: result = conj(x)'*y. More...
 
template<class Scalar >
Scalar get_ele (const VectorBase< Scalar > &v, Ordinal i)
 Get single element: result = v(i). More...
 
template<class Scalar >
void set_ele (Ordinal i, Scalar alpha, const Ptr< VectorBase< Scalar > > &v)
 Set single element: v(i) = alpha. More...
 
template<class Scalar >
void put_scalar (const Scalar &alpha, const Ptr< VectorBase< Scalar > > &y)
 Assign all elements to a scalar: y(i) = alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void copy (const VectorBase< Scalar > &x, const Ptr< VectorBase< Scalar > > &y)
 Vector assignment: y(i) = x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void add_scalar (const Scalar &alpha, const Ptr< VectorBase< Scalar > > &y)
 Add a scalar to all elements: y(i) += alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void scale (const Scalar &alpha, const Ptr< VectorBase< Scalar > > &y)
 Scale all elements by a scalar: y(i) *= alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void abs (const VectorBase< Scalar > &x, const Ptr< VectorBase< Scalar > > &y)
 Element-wise absolute value: y(i) = abs(x(i)), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void reciprocal (const VectorBase< Scalar > &x, const Ptr< VectorBase< Scalar > > &y)
 Element-wise reciprocal: y(i) = 1/x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void ele_wise_prod (const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &v, const Ptr< VectorBase< Scalar > > &y)
 Element-wise product update: y(i) += alpha * x(i) * v(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void pair_wise_max (const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &v, const Ptr< VectorBase< Scalar > > &y)
 Element-wise maximum: y(i) = alpha * max(x(i), v(i)), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void ele_wise_conj_prod (const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &v, const Ptr< VectorBase< Scalar > > &y)
 Element-wise conjugate product update: y(i) += alpha * conj(x(i)) * v(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void ele_wise_scale (const VectorBase< Scalar > &x, const Ptr< VectorBase< Scalar > > &y)
 Element-wise scaling: y(i) *= x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void Vp_StVtV (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &v)
 Element-wise product update: y(i) += alpha * x(i) * v(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void ele_wise_prod_update (const Scalar &alpha, const VectorBase< Scalar > &x, const Ptr< VectorBase< Scalar > > &y)
 Element-wise product update: y(i) *= alpha * x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void pair_wise_max_update (const Scalar &alpha, const VectorBase< Scalar > &x, const Ptr< VectorBase< Scalar > > &y)
 Element-wise maximum update: y(i) = alpha * max(x(i), y(i)), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void Vt_StV (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha, const VectorBase< Scalar > &x)
 Element-wise product update: y(i) *= alpha * x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void ele_wise_divide (const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &v, const Ptr< VectorBase< Scalar > > &y)
 Element-wise division update: y(i) += alpha * x(i) / v(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void linear_combination (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const VectorBase< Scalar > > > &x, const Scalar &beta, const Ptr< VectorBase< Scalar > > &y)
 Linear combination: y(i) = beta*y(i) + sum( alpha[k]*x[k](i), k=0...m-1 ), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void seed_randomize (unsigned int s)
 Seed the random number generator used in randomize(). More...
 
template<class Scalar >
void randomize (Scalar l, Scalar u, const Ptr< VectorBase< Scalar > > &v)
 Random vector generation: v(i) = rand(l,u), , i = 1...v->space()->dim(). More...
 
template<class Scalar >
void assign (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha)
 Assign all elements to a scalar: y(i) = alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void assign (const Ptr< VectorBase< Scalar > > &y, const VectorBase< Scalar > &x)
 Vector assignment: y(i) = x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void Vp_S (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha)
 Add a scalar to all elements: y(i) += alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void Vt_S (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha)
 Scale all elements by a scalar: y(i) *= alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void V_StV (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha, const VectorBase< Scalar > &x)
 Assign scaled vector: y(i) = alpha * x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void Vp_StV (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha, const VectorBase< Scalar > &x)
 AXPY: y(i) = alpha * x(i) + y(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void Vp_V (const Ptr< VectorBase< Scalar > > &y, const VectorBase< Scalar > &x, const Scalar &beta=static_cast< Scalar >(1.0))
 y(i) = x(i) + beta*y(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void V_V (const Ptr< VectorBase< Scalar > > &y, const VectorBase< Scalar > &x)
 y(i) = x(i), i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void V_S (const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha)
 y(i) = alpha, i = 0...y->space()->dim()-1. More...
 
template<class Scalar >
void V_VpV (const Ptr< VectorBase< Scalar > > &z, const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
 z(i) = x(i) + y(i), i = 0...z->space()->dim()-1. More...
 
template<class Scalar >
void V_VmV (const Ptr< VectorBase< Scalar > > &z, const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
 z(i) = x(i) - y(i), i = 0...z->space()->dim()-1. More...
 
template<class Scalar >
void V_StVpV (const Ptr< VectorBase< Scalar > > &z, const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &y)
 z(i) = alpha*x(i) + y(i), i = 0...z->space()->dim()-1. More...
 
template<class Scalar >
void V_VpStV (const Ptr< VectorBase< Scalar > > &z, const VectorBase< Scalar > &x, const Scalar &alpha, const VectorBase< Scalar > &y)
 z(i) = x(i) + alpha*y(i), i = 0...z->space()->dim()-1. More...
 
template<class Scalar >
void V_StVpStV (const Ptr< VectorBase< Scalar > > &z, const Scalar &alpha, const VectorBase< Scalar > &x, const Scalar &beta, const VectorBase< Scalar > &y)
 z(i) = alpha*x(i) + beta*y(i), i = 0...z->space()->dim()-1. More...
 
template<class Scalar >
Scalar min (const VectorBase< Scalar > &x)
 Min element: result = min{ x(i), i = 0...x.space()->dim()-1 } . More...
 
template<class Scalar >
void min (const VectorBase< Scalar > &x, const Ptr< Scalar > &maxEle, const Ptr< Ordinal > &maxIndex)
 Min element and its index: Returns maxEle = x(k) and maxIndex = k such that x(k) <= x(i) for all i = 0...x.space()->dim()-1. More...
 
template<class Scalar >
void minGreaterThanBound (const VectorBase< Scalar > &x, const Scalar &bound, const Ptr< Scalar > &minEle, const Ptr< Ordinal > &minIndex)
 Minimum element greater than some bound and its index: Returns minEle = x(k) and minIndex = k such that x(k) <= x(i) for all i where x(i) > bound. More...
 
template<class Scalar >
Scalar max (const VectorBase< Scalar > &x)
 Max element: result = max{ x(i), i = 1...n } . More...
 
template<class Scalar >
void max (const VectorBase< Scalar > &x, const Ptr< Scalar > &maxEle, const Ptr< Ordinal > &maxIndex)
 Max element and its index: Returns maxEle = x(k) and maxIndex = k such that x(k) >= x(i) for i = 0...x.space()->dim()-1. More...
 
template<class Scalar >
void maxLessThanBound (const VectorBase< Scalar > &x, const Scalar &bound, const Ptr< Scalar > &maxEle, const Ptr< Ordinal > &maxIndex)
 Max element less than bound and its index: Returns maxEle = x(k) and maxIndex = k such that x(k) >= x(i) for all i where x(i) < bound. More...
 

Minimal mathematical functions

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...
 

Space membership

virtual RCP< const
VectorSpaceBase< Scalar > > 
space () const =0
 Return a smart pointer to the vector space that this vector belongs to. More...
 

Reduction/Transformation operator support

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...
 

Vector Cloning

virtual RCP< VectorBase< Scalar > > clone_v () const =0
 Returns a seprate cloned copy of *this vector with the same values but different storage. More...
 

Explicit sub-vector access

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...
 

Protected virtual functions to be overridden by subclasses

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...
 

Additional Inherited Members

- 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...
 
- 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::VectorBase< Scalar >

Abstract interface for finite-dimensional dense vectors.

This interface contains the minimal set of operations needed to define an abstract vector.

Outline

Reduction/transformation operator (RTOp) support

The main feature of this interface is the function applyOp() which is used to implement all types of vector reduction and transformation operations (RTOp) through RTOp operators . Every standard (i.e. BLAS) and nearly every non-standard element-wise operation that can be performed on a set of vectors can be performed efficiently through reduction/transformation operators. More standard vector operations could be included in this interface and allow for specialized implementations but, in general, assuming the sub-vectors are large enough, such implementations would not be significantly faster than those implemented through reduction/transformation operators. There are some operations however that can not always be efficiently implemented with reduction/transformation operators and a few of these important operations are included in this interface. The applyOp() function allows to client to specify a sub-set of the vector elements to include in reduction/transformation operation. This greatly increases the generality of this vector interface as vector objects can be used as sub objects in larger composite vectors and sub-views of a vector can be created.

Collection of pre-written RTOps and wrapper functions

There already exists RTOp-based implementations of several standard vector operations and some convenience functions that wrap these operators and call applyOp(). These wrapper functions can be found here

Explicit vector coefficient access

This interface also allows a client to extract a sub-set of vector coefficients in an explicit form as non-mutable RTOpPack::ConstSubVectorView or mutable RTOpPack::SubVectorView objects using the acquireDetachedView() functions. In general, this is a very inefficient thing to do and should be avoided. However, there are some situations where getting explicit access to the coefficients of a vector is a very reasonable and efficient thing to do (i.e. for vectors in the domain of a multi-vector for instance) and therefore this functionality is supported. These views and the parent vector follow the state behavior outlined here.

Explicit vector coefficient access utilities

Note that client code in general should not directly call the above explicit sub-vector access functions but should use the utility classes ConstDetachedVectorView and DetachedVectorView instead since these are easier an safer in the event that an exception is thrown.

Explicit vector coefficient assignment

In addition to being able to extract an explicit non-mutable and mutable views of some (small?) sub-set of elements, this interface allows a client to set sub-vectors using setSubVector().

Vector is a MultiVectorBase is a LinearOpBase

It is also worth mentioning that that this VectorBase interface class also inherits from MultiVectorBase so that every VectorBase object is also a MultiVectorBase object. This allows any piece of code that accepts MultiVectorBase objects to automatically accept VectorBase objects as well. In addition, since MultiVectorBase inherits from LinearOpBase, then this means that every vector is also a linear operator.

Notes for subclass developers

The support subclass VectorDefaultBase provides default implementations for as many functions as possible and should be considered a first choice for creating concrete subclasses.

Definition at line 370 of file Thyra_OperatorVectorTypes.hpp.

Member Function Documentation

template<class Scalar>
void Thyra::VectorBase< Scalar >::assign ( const VectorBase< Scalar > &  x)
inline

Vector assignment:

y(i) = x(i), i = 0...y->space()->dim()-1.

NVI function.

Definition at line 166 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::randomize ( Scalar  l,
Scalar  u 
)
inline

Random vector generation:

v(i) = rand(l,u), , i = 1...v->space()->dim().

The elements v(i) are randomly generated between [l,u].

NVI function.

Definition at line 178 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::update ( Scalar  alpha,
const VectorBase< Scalar > &  x 
)
inline

AXPY:

y(i) = alpha * x(i) + y(i), i = 0...y->space()->dim()-1.

NVI function.

Definition at line 190 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::linear_combination ( const ArrayView< const Scalar > &  alpha,
const ArrayView< const Ptr< const VectorBase< Scalar > > > &  x,
const Scalar &  beta 
)
inline

Linear combination:

y(i) = beta*y(i) + sum( alpha[k]*x[k](i), k=0...m-1 ), i = 0...y->space()->dim()-1.

Parameters
m[in] Number of vectors x[]
alpha[in] Array (length m) of input scalars.
x[in] Array (length m) of input vectors.
beta[in] Scalar multiplier for y
y[in/out] Target vector that is the result of the linear combination.

This function implements a general linear combination:

  y(i) = beta*y(i) + alpha[0]*x[0](i) + alpha[1]*x[1](i)
         + ... + alpha[m-1]*x[m-1](i), i = 0...y->space()->dim()-1

NVI function.

Definition at line 221 of file Thyra_VectorBase.hpp.

template<class Scalar>
Scalar Thyra::VectorBase< Scalar >::dot ( const VectorBase< Scalar > &  x) const
inline

Euclidean dot product: result = x^H * this.

Definition at line 231 of file Thyra_VectorBase.hpp.

template<class Scalar>
Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm_1 ( ) const
inline

One (1) norm: result = ||v||1.

Definition at line 238 of file Thyra_VectorBase.hpp.

template<class Scalar>
Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm_2 ( ) const
inline

Euclidean (2) norm: result = ||v||2.

NVI function.

Definition at line 246 of file Thyra_VectorBase.hpp.

template<class Scalar>
Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm_2 ( const VectorBase< Scalar > &  x) const
inline

Weighted Euclidean (2) norm: result = ||v||2.

NVI function.

Definition at line 254 of file Thyra_VectorBase.hpp.

template<class Scalar>
Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm_inf ( ) const
inline

Infinity norm: result = ||v||inf.

NVI function.

Definition at line 262 of file Thyra_VectorBase.hpp.

template<class Scalar>
virtual RCP< const VectorSpaceBase<Scalar> > Thyra::VectorBase< Scalar >::space ( ) const
pure virtual

Return a smart pointer to the vector space that this vector belongs to.

A return value of space().get()==NULL is a flag that *this is not fully initialized.

If return.get()!=NULL, then it is required that the object referenced by *return.get() must have lifetime that extends past the lifetime of the returned smart pointer object. However, the object referenced by *return.get() may change if *this is modified so this reference should not be maintained for too long.

New Behavior! It is required that the VectorSpaceBase object embedded in return must be valid past the lifetime of *this vector object.

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

template<class Scalar>
void Thyra::VectorBase< Scalar >::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
inline

Calls applyOpImpl().

Temporary NVI function.

Definition at line 297 of file Thyra_VectorBase.hpp.

template<class Scalar>
virtual RCP<VectorBase<Scalar> > Thyra::VectorBase< Scalar >::clone_v ( ) const
pure virtual

Returns a seprate cloned copy of *this vector with the same values but different storage.

This function performs a deep copy and assigns the values to the returned copy. It basically performs:

auto copy = Thyra::createMember(this->space());
Thyra::assign<Scalar>(copy.ptr(), *this);
return copy;

This function exists to be consistent with the clone functions clone() which creates a LinearOpBase object and clone_mv() which creates a MultiVectorBase object. However, this function is not really necessary because the capability to create a new vector object is already present in the VectorSpaceBase object returned from this->space() followed by the assignment function as shown above.

Subclasses should only consider overriding this function if there they want to be very sophisticated and implement some form of lazy evaluation in case the created object might not actually be modified before it is destroyed. However, this is not advised.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
void Thyra::VectorBase< Scalar >::acquireDetachedView ( const Range1D rng,
RTOpPack::ConstSubVectorView< Scalar > *  sub_vec 
) const
inline

Calls acquireDetachedVectorViewImpl().

Temporary NVI function.

Definition at line 349 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::releaseDetachedView ( RTOpPack::ConstSubVectorView< Scalar > *  sub_vec) const
inline

Calls releaseDetachedVectorViewImpl().

Temporary NVI function.

Definition at line 358 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::acquireDetachedView ( const Range1D rng,
RTOpPack::SubVectorView< Scalar > *  sub_vec 
)
inline

Calls acquireNonconstDetachedVectorViewImpl().

Temporary NVI function.

Definition at line 367 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::commitDetachedView ( RTOpPack::SubVectorView< Scalar > *  sub_vec)
inline

Calls commitDetachedView().

Temporary NVI function.

Definition at line 376 of file Thyra_VectorBase.hpp.

template<class Scalar>
void Thyra::VectorBase< Scalar >::setSubVector ( const RTOpPack::SparseSubVectorT< Scalar > &  sub_vec)
inline

Calls setSubVectorImpl().

Temporary NVI function.

Definition at line 385 of file Thyra_VectorBase.hpp.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::assignVecImpl ( const VectorBase< Scalar > &  x)
protectedpure virtual

Virtual implementation for NVI assign.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::randomizeImpl ( Scalar  l,
Scalar  u 
)
protectedpure virtual
template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::absImpl ( const VectorBase< Scalar > &  x)
protectedpure virtual
template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::reciprocalImpl ( const VectorBase< Scalar > &  x)
protectedpure virtual
template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::eleWiseScaleImpl ( const VectorBase< Scalar > &  x)
protectedpure virtual
template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::updateVecImpl ( Scalar  alpha,
const VectorBase< Scalar > &  x 
)
protectedpure virtual

Virtual implementation for NVI update.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::linearCombinationVecImpl ( const ArrayView< const Scalar > &  alpha,
const ArrayView< const Ptr< const VectorBase< Scalar > > > &  x,
const Scalar &  beta 
)
protectedpure virtual

Virtual implementation for NVI linear_combination.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual Scalar Thyra::VectorBase< Scalar >::dotImpl ( const VectorBase< Scalar > &  x) const
protectedpure virtual

Virtual implementation for NVI dot.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm1Impl ( ) const
protectedpure virtual

Virtual implementation for NVI norm_1.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm2Impl ( ) const
protectedpure virtual

Virtual implementation for NVI norm_2.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::norm2WeightedImpl ( const VectorBase< Scalar > &  x) const
protectedpure virtual
template<class Scalar>
virtual Teuchos::ScalarTraits<Scalar>::magnitudeType Thyra::VectorBase< Scalar >::normInfImpl ( ) const
protectedpure virtual

Virtual implementation for NVI norm_inf.

Implemented in Thyra::VectorDefaultBase< Scalar >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::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
protectedpure virtual

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).

Preconditions:

  • this->space().get()!=NULL (throw std::logic_error)

The vector *this that this function is called on is assumed to be one of the vectors in v[0]...v[nv-1],z[0]...z[nz-1]. This function is generally should not called directly by a client but instead the client should call the nonmember function Thyra::applyOp().

See the documentation for the nonmember function Thyra::applyOp() for a description of what this function does.

Implemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultClusteredSpmdProductVector< Scalar >, Thyra::SpmdVectorDefaultBase< Scalar >, and Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::acquireDetachedVectorViewImpl ( const Range1D rng,
RTOpPack::ConstSubVectorView< Scalar > *  sub_vec 
) const
protectedpure virtual

Get a non-mutable explicit view of a sub-vector.

Parameters
rng[in] The range of the elements to extract the sub-vector view.
sub_vec[in/out] View of the sub-vector. Prior to the first call to this function, sub_vec->set_uninitialized() must be called. Technically *sub_vec owns the memory but this memory can be freed only by calling this->releaseDetachedView(sub_vec).

Preconditions:

  • this->space().get()!=NULL (throw std::logic_error)
  • [!rng.full_range()] (rng.ubound() < this->space()->dim()) == true (throw std::out_of_range)

Postconditions:

  • *sub_vec contains an explicit non-mutable view to the elements in the range full_range(rng,0,this->space()->dim()-1)

This is only a transient view of a sub-vector that is to be immediately used and then released with a call to releaseDetachedView().

Note that calling this function might require some dynamic memory allocations and temporary memory. Therefore, it is critical that this->releaseDetachedView(sub_vec) is called to clean up memory and avoid memory leaks after the sub-vector is used.

Heads Up! Note that client code in general should not directly call this function but should instead use the utility class ConstDetachedVectorView which will also take care of calling releaseDetachedView().

If this->acquireDetachedView(...,sub_vec) was previously called on sub_vec then it may be possible to reuse this memory if it is sufficiently sized. The user is encouraged to make multiple calls to this->acquireDetachedView(...,sub_vec) before this->releaseDetachedView(sub_vec) to finally clean up all of the memory. Of course, the same sub_vec object must be passed to the same vector object for this to work correctly.

Implemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::SpmdVectorDefaultBase< Scalar >, and Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::releaseDetachedVectorViewImpl ( RTOpPack::ConstSubVectorView< Scalar > *  sub_vec) const
protectedpure virtual

Free an explicit view of a sub-vector.

Parameters
sub_vec[in/out] The memory referred to by sub_vec->values() will be released if it was allocated and *sub_vec will be zeroed out using sub_vec->set_uninitialized().

Preconditions:

  • this->space().get()!=NULL (throw std::logic_error)
  • sub_vec must have been passed through a call to this->acquireDetachedView(...,sub_vec)

Postconditions:

  • See RTOpPack::ConstSubVectorView::set_uninitialized() for sub_vec

The sub-vector view must have been allocated by this->acquireDetachedView() first.

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

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::acquireNonconstDetachedVectorViewImpl ( const Range1D rng,
RTOpPack::SubVectorView< Scalar > *  sub_vec 
)
protectedpure virtual

Get a mutable explicit view of a sub-vector.

Parameters
rng[in] The range of the elements to extract the sub-vector view.
sub_vec[in/out] Mutable view of the sub-vector. Prior to the first call to this function sub_vec->set_uninitialized() must have been called for the correct behavior. Technically *sub_vec owns the memory but this memory must be committed and freed by calling this->commitDetachedView(sub_vec) after the client is finished modifying the view.

Preconditions:

  • this->space().get()!=NULL (throw std::logic_error)
  • [!rng.full_range()] rng.ubound() < this->space()->dim() (throw std::out_of_range)

Postconditions:

  • *sub_vec contains an explicit mutable view to the elements in the range full_range(rng,0,this->space()->dim()-1)

This is only a transient view of a sub-vector that is to be immediately used and then committed back with a call to commitDetachedView().

Note that calling this function might require some internal allocations and temporary memory. Therefore, it is critical that this->commitDetachedView(sub_vec) is called to commit the changed entries, clean up memory, and avoid memory leaks after the sub-vector is modified.

Heads Up! Note that client code in general should not directly call this function but should instead use the utility class DetachedVectorView which will also take care of calling commitDetachedView().

If this->acquireDetachedView(...,sub_vec) was previously called on sub_vec then it may be possible to reuse this memory if it is sufficiently sized. The user is encouraged to make multiple calls to this->acquireDetachedView(...,sub_vec) before this->commitDetachedView(sub_vec) is called to finally clean up all of the memory. Of course the same sub_vec object must be passed to the same vector object for this to work correctly.

Changes to the underlying sub-vector are not guaranteed to become permanent until this->acquireDetachedView(...,sub_vec) is called again, or this->commitDetachedView(sub_vec) is called.

Implemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::SpmdVectorDefaultBase< Scalar >, and Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::commitNonconstDetachedVectorViewImpl ( RTOpPack::SubVectorView< Scalar > *  sub_vec)
protectedpure virtual

Commit changes for a mutable explicit view of a sub-vector.

Parameters
sub_vec[in/out] The data in sub_vec->values() will be written back internal storage and the memory referred to by sub_vec->values() will be released if it was allocated and *sub_vec will be zeroed out using sub_vec->set_uninitialized().

Preconditions:

  • this->space().get()!=NULL (throw std::logic_error)
  • sub_vec must have been passed through a call to this->acquireDetachedView(...,sub_vec)

Postconditions:

  • See RTOpPack::SubVectorView::set_uninitialized() for sub_vec
  • *this will be updated according the the changes made to sub_vec

The sub-vector view must have been allocated by this->acquireDetachedView() first.

Implemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::SpmdVectorDefaultBase< Scalar >, and Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar>
virtual void Thyra::VectorBase< Scalar >::setSubVectorImpl ( const RTOpPack::SparseSubVectorT< Scalar > &  sub_vec)
protectedpure virtual

Set a specific sub-vector.

Parameters
sub_vec[in] Represents the elements in the sub-vector to be set.

Preconditions:

  • this->space().get()!=NULL (throw std::logic_error)
  • sub_vec.globalOffset() + sub_vec.subDim() < this->space()->dim() (throw std::out_of_range)

Postconditions:

  • All of the elements in the range [sub_vec.globalOffset(),sub_vec.globalOffset()+sub_vec.subDim()-1] in *this are set to 0.0 except for those that have that have entries in sub_vec which are set to the values specified by (*this)(sub_vec.globalOffset()+vec.localOffset()+sub_vec.indices()[sub_vec.indicesStride()*k]) = vec.values[vec.valueStride()*k], for k = 0..sub_vec.subNz()-1

After this function returns, the corresponding elements in *this vector object will be set equal to those in the input view sub_vec.

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

Friends And Related Function Documentation

template<class Scalar >
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 = 0 
)
related

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).

Parameters
op[in] Reduction/transformation operator to apply over each sub-vector and assemble the intermediate targets into reduct_obj (if reduct_obj != RTOp_REDUCT_OBJ_NULL).
vecs[in] Array (length num_vecs) of a set of pointers to non-mutable vectors to include in the operation. The order of these vectors is significant to op. If vecs.size()==0, then op is called with no non-mutable vectors.
targ_vecs[in] Array (length num_targ_vecs) of a set of pointers to mutable vectors to include in the operation. The order of these vectors is significant to op. If targ_vecs.size()==0, then op is called with no mutable vectors.
reduct_obj[in/out] Target object of the reduction operation. This object must have been created by the op.reduct_obj_create_raw() function first. The reduction operation will be added to *reduct_obj if *reduct_obj has already been through a reduction. By allowing the info in *reduct_obj to be added to the reduction over all of these vectors, the reduction operation can be accumulated over a set of abstract vectors which can be useful for implementing composite vectors for instance. If op.get_reduct_type_num_entries(...) returns num_values == 0, num_indexes == 0 and num_chars == 0 then reduct_obj must be set to null and no reduction will be performed.
global_offset[in] (default = 0) The offset applied to the included vector elements.

Preconditions:

  • [vecs.size() > 0] vecs[k]->space()->isCompatible(*this->space()) == true, for k = 0...vecs.size()-1 (throw Exceptions::IncompatibleVectorSpaces)

  • [targ_vecs.size() > 0] targ_vecs[k]->space()->isCompatible(*this->space()) == true, for k = 0...targ_vecs.size()-1 (throw Exceptions::IncompatibleVectorSpaces)

  • [targ_vecs.size() > 0] The vectors pointed to by targ_vecs[k], for k = 0...targ_vecs.size()-1 must not alias each other or any of the vectors vecs[k], for k = 0...vecs.size()-1. You have be warned!!!!

  • global_offset >= 0 (throw std::invalid_argument)

Postconditions:

  • The vectors in targ_vecs[] may be modified as determined by the definition of op.

  • [reduct_obj!=null] The reduction object reduct_obj contains the combined reduction from the input state of reduct_obj and the reductions that where accumulated during this this function invocation.

Definition at line 770 of file Thyra_VectorBase.hpp.

template<class Scalar >
Scalar sum ( const VectorBase< Scalar > &  v)
related

Sum of vector elements: result = sum( v(i), i = 0...v.space()->dim()-1 ).

template<class Scalar >
Scalar scalarProd ( const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
)
related

Scalar product result = <x,y>.

Returns x.space()->scalarProd(x,y).

template<class Scalar >
Scalar inner ( const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
)
related

Inner/Scalar product result = <x,y>.

Returns x.space()->scalarProd(x,y).

template<class Scalar >
Teuchos::ScalarTraits< Scalar >::magnitudeType norm ( const VectorBase< Scalar > &  v)
related

Natural norm: result = sqrt(<v,v>).

Returns Teuchos::ScalarTraits<Scalar>::squareroot(v.space()->scalarProd(v,v)).

template<class Scalar >
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_1 ( const VectorBase< Scalar > &  v)
related

One (1) norm: result = ||v||1.

template<class Scalar >
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_2 ( const VectorBase< Scalar > &  v)
related

Euclidean (2) norm: result = ||v||2.

template<class Scalar >
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_2 ( const VectorBase< Scalar > &  w,
const VectorBase< Scalar > &  v 
)
related

Weighted Euclidean (2) norm: result = sqrt( sum( w(i)*conj(v(i))*v(i)) ).

template<class Scalar >
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_inf ( const VectorBase< Scalar > &  v_rhs)
related

Infinity norm: result = ||v||inf.

template<class Scalar >
Scalar dot ( const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
)
related

Dot product: result = conj(x)'*y.

template<class Scalar >
Scalar get_ele ( const VectorBase< Scalar > &  v,
Ordinal  i 
)
related

Get single element: result = v(i).

template<class Scalar >
void set_ele ( Ordinal  i,
Scalar  alpha,
const Ptr< VectorBase< Scalar > > &  v 
)
related

Set single element: v(i) = alpha.

template<class Scalar >
void put_scalar ( const Scalar &  alpha,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Assign all elements to a scalar: y(i) = alpha, i = 0...y->space()->dim()-1.

template<class Scalar >
void copy ( const VectorBase< Scalar > &  x,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Vector assignment: y(i) = x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void add_scalar ( const Scalar &  alpha,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Add a scalar to all elements: y(i) += alpha, i = 0...y->space()->dim()-1.

template<class Scalar >
void scale ( const Scalar &  alpha,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Scale all elements by a scalar: y(i) *= alpha, i = 0...y->space()->dim()-1.

This takes care of the special cases of alpha == 0.0 (set y = 0.0) and alpha == 1.0 (don't do anything).

template<class Scalar >
void abs ( const VectorBase< Scalar > &  x,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise absolute value: y(i) = abs(x(i)), i = 0...y->space()->dim()-1.

template<class Scalar >
void reciprocal ( const VectorBase< Scalar > &  x,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise reciprocal: y(i) = 1/x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void ele_wise_prod ( const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  v,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise product update: y(i) += alpha * x(i) * v(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void pair_wise_max ( const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  v,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise maximum: y(i) = alpha * max(x(i), v(i)), i = 0...y->space()->dim()-1.

template<class Scalar >
void ele_wise_conj_prod ( const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  v,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise conjugate product update: y(i) += alpha * conj(x(i)) * v(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void ele_wise_scale ( const VectorBase< Scalar > &  x,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise scaling: y(i) *= x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void Vp_StVtV ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  v 
)
related

Element-wise product update: y(i) += alpha * x(i) * v(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void ele_wise_prod_update ( const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise product update: y(i) *= alpha * x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void pair_wise_max_update ( const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise maximum update: y(i) = alpha * max(x(i), y(i)), i = 0...y->space()->dim()-1.

template<class Scalar >
void Vt_StV ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha,
const VectorBase< Scalar > &  x 
)
related

Element-wise product update: y(i) *= alpha * x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void ele_wise_divide ( const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  v,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Element-wise division update: y(i) += alpha * x(i) / v(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void linear_combination ( const ArrayView< const Scalar > &  alpha,
const ArrayView< const Ptr< const VectorBase< Scalar > > > &  x,
const Scalar &  beta,
const Ptr< VectorBase< Scalar > > &  y 
)
related

Linear combination: y(i) = beta*y(i) + sum( alpha[k]*x[k](i), k=0...m-1 ), i = 0...y->space()->dim()-1.

Parameters
m[in] Number of vectors x[]
alpha[in] Array (length m) of input scalars.
x[in] Array (length m) of input vectors.
beta[in] Scalar multiplier for y
y[in/out] Target vector that is the result of the linear combination.

This function implements a general linear combination:

  y(i) = beta*y(i) + alpha[0]*x[0](i) + alpha[1]*x[1](i)
         + ... + alpha[m-1]*x[m-1](i), i = 0...y->space()->dim()-1
template<class Scalar >
void seed_randomize ( unsigned int  s)
related

Seed the random number generator used in randomize().

Parameters
s[in] The seed for the random number generator.

Note, this just calls Teuchos::TOpRandomize<Scalar>::set_static_seed(s).

template<class Scalar >
void randomize ( Scalar  l,
Scalar  u,
const Ptr< VectorBase< Scalar > > &  v 
)
related

Random vector generation: v(i) = rand(l,u), , i = 1...v->space()->dim().

The elements v->getEle(i) are randomly generated between [l,u].

The seed is set using the above seed_randomize() function.

template<class Scalar >
void assign ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha 
)
related

Assign all elements to a scalar: y(i) = alpha, i = 0...y->space()->dim()-1.

template<class Scalar >
void assign ( const Ptr< VectorBase< Scalar > > &  y,
const VectorBase< Scalar > &  x 
)
related

Vector assignment: y(i) = x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void Vp_S ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha 
)
related

Add a scalar to all elements: y(i) += alpha, i = 0...y->space()->dim()-1.

template<class Scalar >
void Vt_S ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha 
)
related

Scale all elements by a scalar: y(i) *= alpha, i = 0...y->space()->dim()-1.

This takes care of the special cases of alpha == 0.0 (set y = 0.0) and alpha == 1.0 (don't do anything).

template<class Scalar >
void V_StV ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha,
const VectorBase< Scalar > &  x 
)
related

Assign scaled vector: y(i) = alpha * x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void Vp_StV ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha,
const VectorBase< Scalar > &  x 
)
related

AXPY: y(i) = alpha * x(i) + y(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void Vp_V ( const Ptr< VectorBase< Scalar > > &  y,
const VectorBase< Scalar > &  x,
const Scalar &  beta = static_cast< Scalar >(1.0) 
)
related

y(i) = x(i) + beta*y(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void V_V ( const Ptr< VectorBase< Scalar > > &  y,
const VectorBase< Scalar > &  x 
)
related

y(i) = x(i), i = 0...y->space()->dim()-1.

template<class Scalar >
void V_S ( const Ptr< VectorBase< Scalar > > &  y,
const Scalar &  alpha 
)
related

y(i) = alpha, i = 0...y->space()->dim()-1.

template<class Scalar >
void V_VpV ( const Ptr< VectorBase< Scalar > > &  z,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
)
related

z(i) = x(i) + y(i), i = 0...z->space()->dim()-1.

template<class Scalar >
void V_VmV ( const Ptr< VectorBase< Scalar > > &  z,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
)
related

z(i) = x(i) - y(i), i = 0...z->space()->dim()-1.

template<class Scalar >
void V_StVpV ( const Ptr< VectorBase< Scalar > > &  z,
const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const VectorBase< Scalar > &  y 
)
related

z(i) = alpha*x(i) + y(i), i = 0...z->space()->dim()-1.

template<class Scalar >
void V_VpStV ( const Ptr< VectorBase< Scalar > > &  z,
const VectorBase< Scalar > &  x,
const Scalar &  alpha,
const VectorBase< Scalar > &  y 
)
related

z(i) = x(i) + alpha*y(i), i = 0...z->space()->dim()-1.

template<class Scalar >
void V_StVpStV ( const Ptr< VectorBase< Scalar > > &  z,
const Scalar &  alpha,
const VectorBase< Scalar > &  x,
const Scalar &  beta,
const VectorBase< Scalar > &  y 
)
related

z(i) = alpha*x(i) + beta*y(i), i = 0...z->space()->dim()-1.

template<class Scalar >
Scalar min ( const VectorBase< Scalar > &  x)
related

Min element: result = min{ x(i), i = 0...x.space()->dim()-1 } .

template<class Scalar >
void min ( const VectorBase< Scalar > &  x,
const Ptr< Scalar > &  maxEle,
const Ptr< Ordinal > &  maxIndex 
)
related

Min element and its index: Returns maxEle = x(k) and maxIndex = k such that x(k) <= x(i) for all i = 0...x.space()->dim()-1.

Parameters
x[in] Input vector.
minEle[out] The minimum element value.
maxindex[out] The global index of the minimum element. If there is more than one element with the maximum entry then this returns the lowest index in order to make the output independent of the order of operations.

Preconditions:

  • minEle!=NULL
  • minIndex!=NULL
template<class Scalar >
void minGreaterThanBound ( const VectorBase< Scalar > &  x,
const Scalar &  bound,
const Ptr< Scalar > &  minEle,
const Ptr< Ordinal > &  minIndex 
)
related

Minimum element greater than some bound and its index: Returns minEle = x(k) and minIndex = k such that x(k) <= x(i) for all i where x(i) > bound.

Parameters
x[in] Input vector.
bound[in] The upper bound
minEle[out] The minimum element value as defined above.
minIndex[out] The global index of the maximum element. If there is more than one element with the minimum value then this returns the lowest index in order to make the output independent of the order of operations. If no entries are less than bound then minIndex < 0 on return.

Preconditions:

  • minEle!=NULL
  • minIndex!=NULL

Postconditions:

  • If *minIndex > 0 then such an element was found.
  • If *minIndex < 0 then no such element was found.
template<class Scalar >
Scalar max ( const VectorBase< Scalar > &  x)
related

Max element: result = max{ x(i), i = 1...n } .

template<class Scalar >
void max ( const VectorBase< Scalar > &  x,
const Ptr< Scalar > &  maxEle,
const Ptr< Ordinal > &  maxIndex 
)
related

Max element and its index: Returns maxEle = x(k) and maxIndex = k such that x(k) >= x(i) for i = 0...x.space()->dim()-1.

Parameters
x[in] Input vector.
maxEle[out] The maximum element value.
maxindex[out] The global index of the maximum element. If there is more than one element with the maximum value then this returns the lowest index in order to make the output independent of the order of operations.

Preconditions:

  • maxEle!=NULL
  • maxIndex!=NULL
template<class Scalar >
void maxLessThanBound ( const VectorBase< Scalar > &  x,
const Scalar &  bound,
const Ptr< Scalar > &  maxEle,
const Ptr< Ordinal > &  maxIndex 
)
related

Max element less than bound and its index: Returns maxEle = x(k) and maxIndex = k such that x(k) >= x(i) for all i where x(i) < bound.

Parameters
x[in] Input vector.
bound[in] The upper bound
maxEle[out] The maximum element value as defined above.
maxindex[out] The global index of the maximum element. If there is more than one element with the maximum index then this returns the lowest index in order to make the output independent of the order of operations. If no entries are less than bound then minIndex < 0 on return.

Preconditions:

  • maxEle!=NULL
  • maxIndex!=NULL

Postconditions:

  • If *maxIndex > 0 then such an element was found.
  • If *maxIndex < 0 then no such element was found.

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