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

Node subclass that uses a default MultiVectorBase implementation to provide default implementations for as many other functions in MultiVectorBase interface the as is reasonable. More...

#include <Thyra_MultiVectorDefaultBase_decl.hpp>

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

Overridden public member functions from MultiVectorBase

virtual RCP< MultiVectorBase
< Scalar > > 
clone_mv () const
 

Overridden protected member functions from MultiVectorBase

virtual void assignImpl (Scalar alpha)
 Default implementation of assign(scalar) using RTOps. More...
 
virtual void assignMultiVecImpl (const MultiVectorBase< Scalar > &mv)
 Default implementation of assign(MV) using RTOps. More...
 
virtual void scaleImpl (Scalar alpha)
 Default implementation of scale using RTOps. More...
 
virtual void updateImpl (Scalar alpha, const MultiVectorBase< Scalar > &mv)
 Default implementation of update using RTOps. More...
 
virtual void linearCombinationImpl (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
 Default implementation of linear_combination using RTOps. More...
 
virtual void dotsImpl (const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
 Default implementation of dots using RTOps. More...
 
virtual void norms1Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
 Default implementation of norms_1 using RTOps. More...
 
virtual void norms2Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
 Default implementation of norms_2 using RTOps. More...
 
virtual void normsInfImpl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
 Default implementation of norms_inf using RTOps. More...
 
RCP< const MultiVectorBase
< Scalar > > 
contigSubViewImpl (const Range1D &colRng) const
 
RCP< MultiVectorBase< Scalar > > nonconstContigSubViewImpl (const Range1D &colRng)
 
RCP< const MultiVectorBase
< Scalar > > 
nonContigSubViewImpl (const ArrayView< const int > &cols) const
 
RCP< MultiVectorBase< Scalar > > nonconstNonContigSubViewImpl (const ArrayView< const int > &cols)
 
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
 
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
 
virtual void acquireDetachedMultiVectorViewImpl (const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
 
virtual void releaseDetachedMultiVectorViewImpl (RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
 
virtual void acquireNonconstDetachedMultiVectorViewImpl (const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
 
virtual void commitNonconstDetachedMultiVectorViewImpl (RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
 

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...
 
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 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 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 >
- Protected Member Functions inherited from Thyra::LinearOpDefaultBase< Scalar >
std::string description () const
 Default description that gives the label, type, and dimenstion . More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
 Generates a default outputting for all linear operators. More...
 

Detailed Description

template<class Scalar>
class Thyra::MultiVectorDefaultBase< Scalar >

Node subclass that uses a default MultiVectorBase implementation to provide default implementations for as many other functions in MultiVectorBase interface the as is reasonable.

Notes to subclass developers

Only three function overrides are required in order to create a concrete MultiVectorBase subclass: range(), domain() and the non-const version of col(). All of the other functions have default implementations. However, a good implementation will provide optimized overrides of at least the functions apply() and applyTranspose(). The non-const versions of subView() should be overridden if subviews are important. The default implementation will not achieve near-optimal performance in many cases.

Definition at line 70 of file Thyra_MultiVectorDefaultBase_decl.hpp.

Member Function Documentation

template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::MultiVectorDefaultBase< Scalar >::clone_mv ( ) const
virtual

This implementation uses the vector space to create a new multi-vector object and then uses a transformation operator to assign the vector elements. A subclass should only override this function if it can do something more sophisticated (i.e. lazy evaluation) but in general, this is not needed.

Implements Thyra::MultiVectorBase< Scalar >.

Reimplemented in Thyra::DefaultProductMultiVector< Scalar >, and Thyra::VectorDefaultBase< Scalar >.

Definition at line 76 of file Thyra_MultiVectorDefaultBase_def.hpp.

template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::assignImpl ( Scalar  alpha)
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::assignMultiVecImpl ( const MultiVectorBase< Scalar > &  mv)
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::scaleImpl ( Scalar  alpha)
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::updateImpl ( Scalar  alpha,
const MultiVectorBase< Scalar > &  mv 
)
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::linearCombinationImpl ( const ArrayView< const Scalar > &  alpha,
const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &  mv,
const Scalar &  beta 
)
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::dotsImpl ( const MultiVectorBase< Scalar > &  mv,
const ArrayView< Scalar > &  prods 
) const
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::norms1Impl ( const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &  norms) const
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::norms2Impl ( const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &  norms) const
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::normsInfImpl ( const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &  norms) const
protectedvirtual
template<class Scalar >
RCP< const MultiVectorBase< Scalar > > Thyra::MultiVectorDefaultBase< Scalar >::contigSubViewImpl ( const Range1D colRng) const
protectedvirtual
template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::MultiVectorDefaultBase< Scalar >::nonconstContigSubViewImpl ( const Range1D colRng)
protectedvirtual
template<class Scalar >
RCP< const MultiVectorBase< Scalar > > Thyra::MultiVectorDefaultBase< Scalar >::nonContigSubViewImpl ( const ArrayView< const int > &  cols) const
protectedvirtual
template<class Scalar >
RCP< MultiVectorBase< Scalar > > Thyra::MultiVectorDefaultBase< Scalar >::nonconstNonContigSubViewImpl ( const ArrayView< const int > &  cols)
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::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
protectedvirtual
template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::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
protectedvirtual

This implementation calls applyOp() where an array of reduction objects is taken.

Implements Thyra::MultiVectorBase< Scalar >.

Definition at line 457 of file Thyra_MultiVectorDefaultBase_def.hpp.

template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::acquireDetachedMultiVectorViewImpl ( const Range1D rowRng,
const Range1D colRng,
RTOpPack::ConstSubMultiVectorView< Scalar > *  sub_mv 
) const
protectedvirtual

This implementation is based on the vector operation VectorBase::acquireDetachedView() called on the non-changeable vector objects returned from col(). Note that the footprint of the reduction object (both internal and external state) will be O(rowRng.size()*colRng.size()). For serial applications this is fairly reasonable and will not be a major performance penalty. For parallel applications, however, this is a terrible implementation and must be overridden if rowRng.size() is large at all. Although, this function should not even be used in cases where the multi-vector is very large. If a subclass does override this function, it must also override releaseDetachedView() which has an implementation which is a companion to this function's implementation.

Implements Thyra::MultiVectorBase< Scalar >.

Reimplemented in Thyra::DefaultProductMultiVector< Scalar >, Thyra::VectorDefaultBase< Scalar >, Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Thyra::SpmdMultiVectorDefaultBase< Scalar >.

Definition at line 505 of file Thyra_MultiVectorDefaultBase_def.hpp.

template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::releaseDetachedMultiVectorViewImpl ( RTOpPack::ConstSubMultiVectorView< Scalar > *  sub_mv) const
protectedvirtual

This implementation is a companion to the implementation for acquireDetachedView(). If acquireDetachedView() is overridden by a subclass then this function must be overridden also!

Implements Thyra::MultiVectorBase< Scalar >.

Reimplemented in Thyra::DefaultProductMultiVector< Scalar >, Thyra::VectorDefaultBase< Scalar >, and Thyra::SpmdMultiVectorDefaultBase< Scalar >.

Definition at line 555 of file Thyra_MultiVectorDefaultBase_def.hpp.

template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::acquireNonconstDetachedMultiVectorViewImpl ( const Range1D rowRng,
const Range1D colRng,
RTOpPack::SubMultiVectorView< Scalar > *  sub_mv 
)
protectedvirtual

This implementation is based on the vector operation VectorBase::acquireDetachedView() called on the changeable vector objects returned from col(). Note that the footprint of the reduction object (both internal and external state) will be O(rowRng.size()*colRng.size()). For serial applications this is fairly reasonable and will not be a major performance penalty. For parallel applications, however, this is a terrible implementation and must be overridden if rowRng.size() is large at all. Although, this function should not even be used in case where the multi-vector is very large. If a subclass does override this function, it must also override commitDetachedView() which has an implementation which is a companion to this function's implementation.

Implements Thyra::MultiVectorBase< Scalar >.

Reimplemented in Thyra::DefaultProductMultiVector< Scalar >, Thyra::VectorDefaultBase< Scalar >, Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Thyra::SpmdMultiVectorDefaultBase< Scalar >.

Definition at line 565 of file Thyra_MultiVectorDefaultBase_def.hpp.

template<class Scalar >
void Thyra::MultiVectorDefaultBase< Scalar >::commitNonconstDetachedMultiVectorViewImpl ( RTOpPack::SubMultiVectorView< Scalar > *  sub_mv)
protectedvirtual

This implementation is a companion to the default implementation for acquireDetachedView(). If acquireDetachedView() is overridden by a subclass then this function must be overridden also!

Implements Thyra::MultiVectorBase< Scalar >.

Reimplemented in Thyra::DefaultProductMultiVector< Scalar >, Thyra::VectorDefaultBase< Scalar >, Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Thyra::SpmdMultiVectorDefaultBase< Scalar >.

Definition at line 586 of file Thyra_MultiVectorDefaultBase_def.hpp.


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