Thyra
Version of the Day
|
Base node implementation class for SPMD multi-vectors. More...
#include <Thyra_SpmdMultiVectorDefaultBase_decl.hpp>
Constructors / initializers / accessors | |
SpmdMultiVectorDefaultBase () | |
Overridden public functions from MultiVectorAdapterBase | |
RCP< const ScalarProdVectorSpaceBase < Scalar > > | rangeScalarProdVecSpc () const |
Returns spmdSpace . More... | |
Protected funtions overridden from SpmdMultiVectorBase. | |
RTOpPack::SubMultiVectorView < Scalar > | getNonconstLocalSubMultiVectorImpl () |
RTOpPack::ConstSubMultiVectorView < Scalar > | getLocalSubMultiVectorImpl () const |
Protected functions overridden from MultiVectorBase | |
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 |
void | acquireDetachedMultiVectorViewImpl (const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const |
void | releaseDetachedMultiVectorViewImpl (RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const |
void | acquireNonconstDetachedMultiVectorViewImpl (const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv) |
void | commitNonconstDetachedMultiVectorViewImpl (RTOpPack::SubMultiVectorView< Scalar > *sub_mv) |
Protected functions overridden from MultiVectorAdapterBase | |
void | euclideanApply (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const |
Uses GEMM() and Teuchos::reduceAll() to implement. More... | |
Protected functions for subclasses to call. | |
virtual void | updateSpmdSpace () |
Subclasses should call whenever the structure of the VectorSpaceBase changes. More... | |
Range1D | validateRowRange (const Range1D &rowRng) const |
Validate and resize the row range. More... | |
Range1D | validateColRange (const Range1D &rowCol) const |
Validate and resize the column range. More... | |
Additional Inherited Members | |
Public Member Functions inherited from Thyra::SpmdMultiVectorBase< Scalar > | |
RCP< const SpmdVectorSpaceBase < Scalar > > | spmdSpace () const |
Returns the SPMD vector space object for the range of *this multi-vector. More... | |
RTOpPack::SubMultiVectorView < Scalar > | getNonconstLocalSubMultiVector () |
Get a non-const generalized view of local multi-vector data. More... | |
RTOpPack::ConstSubMultiVectorView < Scalar > | getLocalSubMultiVector () const |
Get a const generalized view of local multi-vector data. More... | |
void | getNonconstLocalData (const Ptr< ArrayRCP< Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) |
Returns a non-const pointer to a Fortran-style view of the local multi-vector data. More... | |
void | getLocalData (const Ptr< ArrayRCP< const Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) const |
Returns a const pointer to a Fortran-style view of the local multi-vector data. 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... | |
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 > | |
bool | opSupported (EOpTransp M_trans) const |
Return if the M_trans operation of apply() is supported or not. More... | |
void | apply (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const |
Apply the linear operator to a multi-vector : Y = alpha*op(M)*X + beta*Y . More... | |
Public Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar > | |
bool | rowStatIsSupported (const RowStatLinearOpBaseUtils::ERowStat rowStat) const |
Determine if a given row stat is supported. More... | |
void | getRowStat (const RowStatLinearOpBaseUtils::ERowStat rowStat, const Ptr< VectorBase< Scalar > > &rowStatVec) const |
Get some statistics about a supported row. More... | |
Public Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar > | |
bool | supportsScaleLeft () const |
Determines if this objects supports left scaling. More... | |
bool | supportsScaleRight () const |
Determines if this objects supports right scaling. More... | |
void | scaleLeft (const VectorBase< Scalar > &row_scaling) |
Left scales operator with diagonal scaling operator. More... | |
void | scaleRight (const VectorBase< Scalar > &col_scaling) |
Right scales operator with diagonal scaling operator. More... | |
Public Member Functions inherited from Thyra::MultiVectorAdapterBase< Scalar > | |
virtual RCP< const ScalarProdVectorSpaceBase < Scalar > > | domainScalarProdVecSpc () const =0 |
RCP< const VectorSpaceBase < Scalar > > | range () const |
Returns this->rangeScalarProdVecSpc() More... | |
RCP< const VectorSpaceBase < Scalar > > | domain () const |
Returns this->domainScalarProdVecSpc() More... | |
Public Member Functions inherited from Thyra::MultiVectorDefaultBase< Scalar > | |
virtual RCP< MultiVectorBase < Scalar > > | clone_mv () const |
Protected Member Functions inherited from Thyra::SpmdMultiVectorBase< Scalar > | |
virtual RCP< const SpmdVectorSpaceBase< Scalar > > | spmdSpaceImpl () const =0 |
Virtual implementation for spmdSpace(). More... | |
virtual void | getNonconstLocalMultiVectorDataImpl (const Ptr< ArrayRCP< Scalar > > &localValues, const Ptr< Ordinal > &leadingDim)=0 |
Virtual implementation for getNonconstLocalData(). More... | |
virtual void | getLocalMultiVectorDataImpl (const Ptr< ArrayRCP< const Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) const =0 |
Virtual implementation for getLocalData(). 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 > | |
Protected Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar > | |
Protected Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar > | |
Protected Member Functions inherited from Thyra::MultiVectorAdapterBase< Scalar > | |
bool | opSupportedImpl (EOpTransp M_trans) const |
void | applyImpl (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const |
Protected Member Functions inherited from Thyra::MultiVectorDefaultBase< Scalar > | |
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 | 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 |
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... | |
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 node implementation class for SPMD multi-vectors.
By inheriting from this base class, multi-vector implementations allow their multi-vector objects to be seamlessly combined with other SPMD multi-vector objects (of different concrete types) in applyOp()
and apply()
. A big part of this protocol is that every multi-vector object can expose an SpmdVectorSpaceBase
object through the virtual function spmdSpace()
.
Notes to subclass developers
Concrete subclasses must override only five functions: spmdSpace()
, getLocalData(const Scalar**,Ordinal*)
, freeLocalData(const Scalar**,Ordinal*)
, getLocalData(Scalar**,Ordinal*)
, and commitLocalData(Scalar**,Ordinal*)
. Note that overriding the spmdSpace()
function requires implementing or using a pre-implemented concrete SpmdVectorSpaceBase
subclass.
If the acquireDetachedView()
functions are ever called with index ranges outside of those of the local process, then the default implementations in MultiVectorBase
of all of the functions (const
) MultiVectorBase::acquireDetachedView()
, MultiVectorBase::releaseDetachedView()
, (non-const
) MultiVectorBase::acquireDetachedView()
and MultiVectorBase::commitDetachedView()
are called in instead. Alternatively, a subclass could provide more specialized implementations of these functions (for more efficient gather/scatter operations) if desired but this should not be needed for most use cases.
It is interesting to note that in the above use case that the explicit subvector access functions call on its default implementation defined in MultiVectorBase
(which calls on applyOp()
) and the operator will be properly applied since the version of applyOp()
implemented in this class will only request local vector data and hence there will only be two levels of recursion for any call to an explicit subvector access function. This is a truly elegant result.
Note that a multi-vector subclass derived from this base class must only be directly used in SPMD mode for this to work properly.
Definition at line 66 of file Thyra_SpmdMultiVectorDefaultBase_decl.hpp.
Thyra::SpmdMultiVectorDefaultBase< Scalar >::SpmdMultiVectorDefaultBase | ( | ) |
Definition at line 45 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
virtual |
Returns spmdSpace
.
Implements Thyra::MultiVectorAdapterBase< Scalar >.
Definition at line 59 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::SpmdMultiVectorBase< Scalar >.
Definition at line 72 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::SpmdMultiVectorBase< Scalar >.
Definition at line 91 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 112 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 200 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::MultiVectorBase< Scalar >.
Definition at line 233 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 252 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 290 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Uses GEMM()
and Teuchos::reduceAll()
to implement.
ToDo: Finish documentation!
Implements Thyra::MultiVectorAdapterBase< Scalar >.
Reimplemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 312 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protectedvirtual |
Subclasses should call whenever the structure of the VectorSpaceBase changes.
WARNING! This function can be overridden by subclasses but this particular function implementation must be called back from within any override (i.e. call SpmdMultiVectorDefaultBase<Scalar>::updateSpmdSpace();
).
Definition at line 633 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protected |
Validate and resize the row range.
This function throws an exception if the input range is invalid
Definition at line 654 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.
|
protected |
Validate and resize the column range.
This function throws an exception if the input range is invalid
Definition at line 670 of file Thyra_SpmdMultiVectorDefaultBase_def.hpp.