| 
    Thyra Version of the Day
    
   | 
 
Base interface class for SPMD multi-vectors. More...
#include <Thyra_SpmdMultiVectorBase.hpp>

Public non-virtual interface functions | |
| RCP< const SpmdVectorSpaceBase< Scalar > > | spmdSpace () const | 
Returns the SPMD vector space object for the range of *this multi-vector.   | |
| RTOpPack::SubMultiVectorView< Scalar > | getNonconstLocalSubMultiVector () | 
| Get a non-const generalized view of local multi-vector data.   | |
| RTOpPack::ConstSubMultiVectorView< Scalar > | getLocalSubMultiVector () const | 
| Get a const generalized view of local multi-vector data.   | |
| 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.   | |
| 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.   | |
Virtual functions to be overridden by sublcasses. | |
| virtual RCP< const SpmdVectorSpaceBase< Scalar > > | spmdSpaceImpl () const =0 | 
| Virtual implementation for spmdSpace().   | |
| virtual RTOpPack::SubMultiVectorView< Scalar > | getNonconstLocalSubMultiVectorImpl ()=0 | 
| Virtual implementation for getNonconstLocalSubMultiVector().   | |
| virtual RTOpPack::ConstSubMultiVectorView< Scalar > | getLocalSubMultiVectorImpl () const =0 | 
| Virtual implementation for getLocalSubMultiVector().   | |
| virtual void | getNonconstLocalMultiVectorDataImpl (const Ptr< ArrayRCP< Scalar > > &localValues, const Ptr< Ordinal > &leadingDim)=0 | 
| Virtual implementation for getNonconstLocalData().   | |
| virtual void | getLocalMultiVectorDataImpl (const Ptr< ArrayRCP< const Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) const =0 | 
| Virtual implementation for getLocalData().   | |
Additional Inherited Members | |
  Public Member Functions inherited from Thyra::MultiVectorBase< Scalar > | |
| void | assign (Scalar alpha) | 
| V = alpha.   | |
| void | assign (const MultiVectorBase< Scalar > &mv) | 
| V = mv.   | |
| 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),   | |
| void | dots (const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const | 
| Column-wise Euclidean dot product.   | |
| void | norms_1 (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const | 
| Column-wise 1-norms.   | |
| void | norms_2 (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const | 
| Column-wise 2-norms.   | |
| void | norms_inf (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const | 
| Column-wise infinity-norms.   | |
| RCP< const VectorBase< Scalar > > | col (Ordinal j) const | 
| Calls colImpl().   | |
| RCP< VectorBase< Scalar > > | col (Ordinal j) | 
| Calls nonconstColImpl().   | |
| RCP< const MultiVectorBase< Scalar > > | subView (const Range1D &colRng) const | 
| Calls contigSubViewImpl().   | |
| RCP< MultiVectorBase< Scalar > > | subView (const Range1D &colRng) | 
| Calls nonconstContigSubViewImpl().   | |
| RCP< const MultiVectorBase< Scalar > > | subView (const ArrayView< const int > &cols) const | 
| nonContigSubViewImpl().   | |
| RCP< MultiVectorBase< Scalar > > | subView (const ArrayView< const int > &cols) | 
| nonconstNonContigSubViewImpl().   | |
| 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().   | |
| 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().   | |
| void | acquireDetachedView (const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const | 
| Calls acquireDetachedMultiVectorViewImpl().   | |
| void | releaseDetachedView (RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const | 
| Calls releaseDetachedMultiVectorViewImpl().   | |
| void | acquireDetachedView (const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv) | 
| Calls acquireNonconstDetachedMultiVectorViewImpl().   | |
| void | commitDetachedView (RTOpPack::SubMultiVectorView< Scalar > *sub_mv) | 
| Calls commitNonconstDetachedMultiVectorViewImpl().   | |
| virtual RCP< MultiVectorBase< Scalar > > | clone_mv () const =0 | 
| Clone the multi-vector object (if supported).   | |
| RCP< const LinearOpBase< Scalar > > | clone () const | 
This function is simply overridden to return this->clone_mv().   | |
  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.   | |
| virtual RCP< const VectorSpaceBase< Scalar > > | domain () const =0 | 
Return a smart pointer for the domain space for this operator.   | |
| bool | opSupported (EOpTransp M_trans) const | 
Return if the M_trans operation of apply() is supported or not.   | |
| 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.   | |
  Public Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar > | |
| bool | rowStatIsSupported (const RowStatLinearOpBaseUtils::ERowStat rowStat) const | 
| Determine if a given row stat is supported.   | |
| void | getRowStat (const RowStatLinearOpBaseUtils::ERowStat rowStat, const Ptr< VectorBase< Scalar > > &rowStatVec) const | 
| Get some statistics about a supported row.   | |
  Public Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar > | |
| bool | supportsScaleLeft () const | 
| Determines if this objects supports left scaling.   | |
| bool | supportsScaleRight () const | 
| Determines if this objects supports right scaling.   | |
| void | scaleLeft (const VectorBase< Scalar > &row_scaling) | 
| Left scales operator with diagonal scaling operator.   | |
| void | scaleRight (const VectorBase< Scalar > &col_scaling) | 
| Right scales operator with diagonal scaling operator.   | |
  Protected Member Functions inherited from Thyra::MultiVectorBase< Scalar > | |
| virtual void | assignImpl (Scalar alpha)=0 | 
| Virtual implementation for NVI assign(Scalar).   | |
| virtual void | assignMultiVecImpl (const MultiVectorBase< Scalar > &mv)=0 | 
| Virtual implementation for NVI assign(MV).   | |
| virtual void | scaleImpl (Scalar alpha)=0 | 
| Virtual implementation for NVI scale().   | |
| virtual void | updateImpl (Scalar alpha, const MultiVectorBase< Scalar > &mv)=0 | 
| Virtual implementation for NVI update().   | |
| 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().   | |
| virtual void | dotsImpl (const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const =0 | 
| Virtual implementation for NVI dots().   | |
| virtual void | norms1Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const =0 | 
| Virtual implementation for NVI norms_1().   | |
| virtual void | norms2Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const =0 | 
| Virtual implementation for NVI norms_2().   | |
| virtual void | normsInfImpl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const =0 | 
| Virtual implementation for NVI norms_inf().   | |
| virtual RCP< const VectorBase< Scalar > > | colImpl (Ordinal j) const | 
| Return a non-changeable view of a constituent column vector.   | |
| virtual RCP< VectorBase< Scalar > > | nonconstColImpl (Ordinal j)=0 | 
| Return a changeable view of a constituent column vector.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| virtual void | releaseDetachedMultiVectorViewImpl (RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const =0 | 
| Free a non-changeable explicit view of a sub-multi-vector.   | |
| 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.   | |
| virtual void | commitNonconstDetachedMultiVectorViewImpl (RTOpPack::SubMultiVectorView< Scalar > *sub_mv)=0 | 
| Commit changes for a changeable explicit view of a sub-multi-vector.   | |
| 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) | 
| void | absRowSum (const Teuchos::Ptr< Thyra::VectorBase< Scalar > > &output) const | 
| void | absColSum (const Teuchos::Ptr< Thyra::VectorBase< Scalar > > &output) const | 
  Protected Member Functions inherited from Thyra::LinearOpBase< Scalar > | |
| virtual bool | opSupportedImpl (EOpTransp M_trans) const =0 | 
| Override in subclass.   | |
| 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.   | |
  Protected Member Functions inherited from Thyra::RowStatLinearOpBase< Scalar > | |
  Protected Member Functions inherited from Thyra::ScaledLinearOpBase< Scalar > | |
  Related Symbols 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.   | |
| 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.   | |
| template<class Scalar > | |
| void | norms (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) | 
| Column-wise multi-vector natural norm.   | |
| 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.   | |
| template<class Scalar > | |
| void | norms_1 (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) | 
| Column-wise multi-vector one norm.   | |
| template<class Scalar > | |
| void | norms_2 (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) | 
| Column-wise multi-vector 2 (Euclidean) norm.   | |
| template<class Scalar > | |
| void | norms_inf (const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) | 
| Column-wise multi-vector infinity norm.   | |
| template<class Scalar > | |
| Array< typename ScalarTraits< Scalar >::magnitudeType > | norms_inf (const MultiVectorBase< Scalar > &V) | 
| Column-wise multi-vector infinity norm.   | |
| template<class Scalar > | |
| void | dots (const MultiVectorBase< Scalar > &V1, const MultiVectorBase< Scalar > &V2, const ArrayView< Scalar > &dots) | 
| Multi-vector dot product.   | |
| template<class Scalar > | |
| void | sums (const MultiVectorBase< Scalar > &V, const ArrayView< Scalar > &sums) | 
| Multi-vector column sum.   | |
| template<class Scalar > | |
| ScalarTraits< Scalar >::magnitudeType | norm_1 (const MultiVectorBase< Scalar > &V) | 
| Take the induced matrix one norm of a multi-vector.   | |
| template<class Scalar > | |
| void | scale (Scalar alpha, const Ptr< MultiVectorBase< Scalar > > &V) | 
| V = alpha*V.   | |
| 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).   | |
| template<class Scalar > | |
| void | assign (const Ptr< MultiVectorBase< Scalar > > &V, Scalar alpha) | 
| V = alpha.   | |
| template<class Scalar > | |
| void | assign (const Ptr< MultiVectorBase< Scalar > > &V, const MultiVectorBase< Scalar > &U) | 
| V = U.   | |
| template<class Scalar > | |
| void | update (Scalar alpha, const MultiVectorBase< Scalar > &U, const Ptr< MultiVectorBase< Scalar > > &V) | 
| alpha*U + V -> V.   | |
| 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 ,,,   | |
| 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.   | |
| 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.   | |
| template<class Scalar > | |
| void | randomize (Scalar l, Scalar u, const Ptr< MultiVectorBase< Scalar > > &V) | 
| Generate a random multi-vector with elements uniformly distributed elements.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
| 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.   | |
  Related Symbols 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.   | |
| template<class Scalar > | |
| bool | isPartiallyInitialized (const LinearOpBase< Scalar > &M) | 
| Determines if a linear operator is in the "Partially Initialized" state or not.   | |
| template<class Scalar > | |
| bool | isFullyInitialized (const LinearOpBase< Scalar > &M) | 
| Determines if a linear operator is in the "Fully Initialized" state or not.   | |
| template<class Scalar > | |
| bool | opSupported (const LinearOpBase< Scalar > &M, EOpTransp M_trans) | 
| Determines if an operation is supported for a single scalar type.   | |
| 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(...).   | |
| 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>(...).   | |
Base interface 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 function spmdSpace(). 
Definition at line 35 of file Thyra_SpmdMultiVectorBase.hpp.
      
  | 
  inline | 
Returns the SPMD vector space object for the range of *this multi-vector. 
Definition at line 45 of file Thyra_SpmdMultiVectorBase.hpp.
      
  | 
  inline | 
Get a non-const generalized view of local multi-vector data.
Definition at line 50 of file Thyra_SpmdMultiVectorBase.hpp.
      
  | 
  inline | 
Get a const generalized view of local multi-vector data.
Definition at line 55 of file Thyra_SpmdMultiVectorBase.hpp.
      
  | 
  inline | 
Returns a non-const pointer to a Fortran-style view of the local multi-vector data. 
| localValues | [out] On output *localValues will point to the first element in the first column of the local multi-vector stored as a column-major dense Fortran-style matrix. | 
| leadingDim | [out] On output *leadingDim gives the leading dimension of the Fortran-style local multi-vector. | 
Preconditions:
localValues!=NULL leadingDim!=NULL Preconditions:
*localValues!=NULL *leadingDim!=0 Definition at line 78 of file Thyra_SpmdMultiVectorBase.hpp.
      
  | 
  inline | 
Returns a const pointer to a Fortran-style view of the local multi-vector data. 
| localValues | [out] On output *localValues will point to the first element in the first column of the local multi-vector stored as a column-major dense Fortran-style matrix. | 
| leadingDim | [out] On output *leadingDim gives the leading dimension of the Fortran-style local multi-vector. | 
Preconditions:
localValues!=NULL leadingDim!=NULL Preconditions:
*localValues!=NULL *leadingDim!=0 Definition at line 103 of file Thyra_SpmdMultiVectorBase.hpp.
      
  | 
  protectedpure virtual | 
Virtual implementation for spmdSpace().
Implemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::DefaultSpmdVector< Scalar >.
      
  | 
  protectedpure virtual | 
Virtual implementation for getNonconstLocalSubMultiVector().
Implemented in Thyra::SpmdMultiVectorDefaultBase< Scalar >, and Thyra::SpmdVectorDefaultBase< Scalar >.
      
  | 
  protectedpure virtual | 
Virtual implementation for getLocalSubMultiVector().
Implemented in Thyra::SpmdMultiVectorDefaultBase< Scalar >, and Thyra::SpmdVectorDefaultBase< Scalar >.
      
  | 
  protectedpure virtual | 
Virtual implementation for getNonconstLocalData().
Implemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::SpmdVectorDefaultBase< Scalar >.
      
  | 
  protectedpure virtual | 
Virtual implementation for getLocalData().
Implemented in Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::SpmdVectorDefaultBase< Scalar >.