| Thyra
    Version of the Day
    | 
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>

| 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 thisoperator.  More... | |
| virtual RCP< const VectorSpaceBase< Scalar > > | domain () const =0 | 
| Return a smart pointer for the domain space for thisoperator.  More... | |
| bool | opSupported (EOpTransp M_trans) const | 
| Return if the M_transoperation ofapply()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... | |
|  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 ), fori = 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... | |
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.
| 
 | 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.
| 
 | protectedvirtual | 
Default implementation of assign(scalar) using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 94 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of assign(MV) using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 104 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of scale using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 114 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of update using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 134 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of linear_combination using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 147 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of dots using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 221 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of norms_1 using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 173 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of norms_2 using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 197 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Default implementation of norms_inf using RTOps.
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductVector< Scalar >, Thyra::DefaultMultiVectorProductVector< Scalar >, Thyra::DefaultColumnwiseMultiVector< Scalar >, Thyra::TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, Thyra::DefaultProductMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 246 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 271 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 299 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 327 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::VectorDefaultBase< Scalar >, Thyra::DefaultSpmdMultiVector< Scalar >, and Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 363 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | protectedvirtual | 
This implementation calls VectorBase::applyOp() on each column this->col(j) for j = 0 ... this->range()->dim()-1. 
Implements Thyra::MultiVectorBase< Scalar >.
Reimplemented in Thyra::DefaultProductMultiVector< Scalar >, Thyra::TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >, and Thyra::SpmdMultiVectorDefaultBase< Scalar >.
Definition at line 398 of file Thyra_MultiVectorDefaultBase_def.hpp.
| 
 | 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.
| 
 | 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.
| 
 | 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.
| 
 | 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.
| 
 | 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.
 1.8.5
 1.8.5