Thyra
Version of the Day
|
Base interface for product multi-vectors. More...
#include <Thyra_ProductMultiVectorBase.hpp>
Public Member Functions | |
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... | |
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... | |
Additional Inherited Members | |
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 > | |
Related Functions inherited from Thyra::MultiVectorBase< Scalar > | |
template<class Scalar > | |
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=0) |
Apply a reduction/transformation operator column by column and return an array of the reduction objects. More... | |
template<class Scalar > | |
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=0) |
Apply a reduction/transformation operator column by column and reduce the intermediate reduction objects into one reduction object. More... | |
template<class Scalar > | |
void | norms (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) |
Column-wise multi-vector natural norm. More... | |
template<class Scalar , class NormOp > | |
void | reductions (const MultiVectorBase< Scalar > &V, const NormOp &op, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) |
Column-wise multi-vector reductions. More... | |
template<class Scalar > | |
void | norms_1 (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) |
Column-wise multi-vector one norm. More... | |
template<class Scalar > | |
void | norms_2 (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) |
Column-wise multi-vector 2 (Euclidean) norm. More... | |
template<class Scalar > | |
void | norms_inf (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) |
Column-wise multi-vector infinity norm. More... | |
template<class Scalar > | |
Array< typename ScalarTraits < Scalar >::magnitudeType > | norms_inf (const MultiVectorBase< Scalar > &V) |
Column-wise multi-vector infinity norm. More... | |
template<class Scalar > | |
void | dots (const MultiVectorBase< Scalar > &V1, const MultiVectorBase< Scalar > &V2, const ArrayView< Scalar > &dots) |
Multi-vector dot product. More... | |
template<class Scalar > | |
void | sums (const MultiVectorBase< Scalar > &V, const ArrayView< Scalar > &sums) |
Multi-vector column sum. More... | |
template<class Scalar > | |
ScalarTraits< Scalar > ::magnitudeType | norm_1 (const MultiVectorBase< Scalar > &V) |
Take the induced matrix one norm of a multi-vector. More... | |
template<class Scalar > | |
void | scale (Scalar alpha, const Ptr< MultiVectorBase< Scalar > > &V) |
V = alpha*V. More... | |
template<class Scalar > | |
void | scaleUpdate (const VectorBase< Scalar > &a, const MultiVectorBase< Scalar > &U, const Ptr< MultiVectorBase< Scalar > > &V) |
A*U + V -> V (where A is a diagonal matrix with diagonal a). More... | |
template<class Scalar > | |
void | assign (const Ptr< MultiVectorBase< Scalar > > &V, Scalar alpha) |
V = alpha. More... | |
template<class Scalar > | |
void | assign (const Ptr< MultiVectorBase< Scalar > > &V, const MultiVectorBase< Scalar > &U) |
V = U. More... | |
template<class Scalar > | |
void | update (Scalar alpha, const MultiVectorBase< Scalar > &U, const Ptr< MultiVectorBase< Scalar > > &V) |
alpha*U + V -> V. More... | |
template<class Scalar > | |
void | update (const ArrayView< const Scalar > &alpha, Scalar beta, const MultiVectorBase< Scalar > &U, const Ptr< MultiVectorBase< Scalar > > &V) |
alpha[j]*beta*U(j) + V(j) - > V(j), for j = 0 ,,, More... | |
template<class Scalar > | |
void | update (const MultiVectorBase< Scalar > &U, const ArrayView< const Scalar > &alpha, Scalar beta, const Ptr< MultiVectorBase< Scalar > > &V) |
U(j) + alpha[j]*beta*V(j) - > V(j), for j = 0 ,,, U.domain()->dim()-1. More... | |
template<class Scalar > | |
void | linear_combination (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &X, const Scalar &beta, const Ptr< MultiVectorBase< Scalar > > &Y) |
Y.col(j)(i) = beta*Y.col(j)(i) + sum( alpha[k]*X[k].col(j)(i), k=0...m-1 ) , for i = 0...Y->range()->dim()-1 , j = 0...Y->domain()->dim()-1 . More... | |
template<class Scalar > | |
void | randomize (Scalar l, Scalar u, const Ptr< MultiVectorBase< Scalar > > &V) |
Generate a random multi-vector with elements uniformly distributed elements. More... | |
template<class Scalar > | |
void | Vt_S (const Ptr< MultiVectorBase< Scalar > > &Z, const Scalar &alpha) |
Z(i,j) *= alpha, i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1 . More... | |
template<class Scalar > | |
void | Vp_S (const Ptr< MultiVectorBase< Scalar > > &Z, const Scalar &alpha) |
Z(i,j) += alpha, i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1 . More... | |
template<class Scalar > | |
void | Vp_V (const Ptr< MultiVectorBase< Scalar > > &Z, const MultiVectorBase< Scalar > &X) |
Z(i,j) += X(i,j), i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1 . More... | |
template<class Scalar > | |
void | V_VpV (const Ptr< MultiVectorBase< Scalar > > &Z, const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y) |
Z(i,j) = X(i,j) + Y(i,j), i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1 . More... | |
template<class Scalar > | |
void | V_VmV (const Ptr< MultiVectorBase< Scalar > > &Z, const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y) |
Z(i,j) = X(i,j) - Y(i,j), i = 0...Z->range()->dim()-1, j = 0...Z->domain()->dim()-1 . More... | |
template<class Scalar > | |
void | V_StVpV (const Ptr< MultiVectorBase< Scalar > > &Z, const Scalar &alpha, const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y) |
Z(i,j) = alpha*X(i,j) + Y(i), i = 0...z->space()->dim()-1 , , j = 0...Z->domain()->dim()-1. More... | |
Related Functions inherited from Thyra::LinearOpBase< Scalar > | |
template<class Scalar > | |
bool | isFullyUninitialized (const LinearOpBase< Scalar > &M) |
Determines if a linear operator is in the "Fully Uninitialized" state or not. More... | |
template<class Scalar > | |
bool | isPartiallyInitialized (const LinearOpBase< Scalar > &M) |
Determines if a linear operator is in the "Partially Initialized" state or not. More... | |
template<class Scalar > | |
bool | isFullyInitialized (const LinearOpBase< Scalar > &M) |
Determines if a linear operator is in the "Fully Initialized" state or not. More... | |
template<class Scalar > | |
bool | opSupported (const LinearOpBase< Scalar > &M, EOpTransp M_trans) |
Determines if an operation is supported for a single scalar type. More... | |
template<class Scalar > | |
void | apply (const LinearOpBase< Scalar > &M, const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha=static_cast< Scalar >(1.0), const Scalar beta=static_cast< Scalar >(0.0)) |
Non-member function call for M.apply(...) . More... | |
void | apply (const LinearOpBase< double > &M, const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha=1.0, const double beta=0.0) |
Calls apply<double>(...) . More... | |
Base interface for product multi-vectors.
This class defines an abstract interface for a multi-vector that is built out of the one or more other multi-vectors to form a product multi-vector. 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 DefaultProductMultiVector
.
ToDo: Finish documentation!
Definition at line 34 of file Thyra_ProductMultiVectorBase.hpp.
|
pure virtual |
Returns the associated product vector space that represents the range.
If *this
is uninitialized then return.get()==NULL
.
Implemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::DefaultClusteredSpmdProductVector< Scalar >.
|
pure virtual |
Return if the kth
multi-vector block is const-only.
k | [in] The (zero-based) kth block index specifying which multi-vector block to access. |
Preconditions:
productSpace().get()!=NULL
0 <= k && k < productSpace()->numBlocks()
Implemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::DefaultClusteredSpmdProductVector< Scalar >.
|
pure virtual |
Returns a non-persisting non-const
view of the zero-based kth
block multi-vector.
k | [in] The (zero-based) kth block index specifying which multi-vector block to access. |
Preconditions:
productSpace().get()!=NULL
0 <= k && k < productSpace()->numBlocks()
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 >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::DefaultClusteredSpmdProductVector< Scalar >.
|
pure virtual |
Returns a non-persisting const
view of the (zero-based) kth
block multi-vector.
k | [in] The (zero-based) kth block index specifying which multi-vector block to access. |
Preconditions:
productSpace().get()!=NULL
0 <= k && k < productSpace()->numBlocks()
Implemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::DefaultClusteredSpmdProductVector< Scalar >.