Thyra
Version of the Day
|
Default subclass for MultiVectorBase
implemented using columns of separate abstract vectors.
More...
#include <Thyra_DefaultColumnwiseMultiVector_decl.hpp>
Constructors/Initializers | |
DefaultColumnwiseMultiVector () | |
Construct to initialized. More... | |
DefaultColumnwiseMultiVector (const RCP< VectorBase< Scalar > > &col_vec) | |
Calls initialize() . More... | |
DefaultColumnwiseMultiVector (const RCP< const VectorSpaceBase< Scalar > > &range, const RCP< const VectorSpaceBase< Scalar > > &domain, const ArrayView< const RCP< VectorBase< Scalar > > > &col_vecs=Teuchos::null) | |
Calls initialize() . More... | |
void | initialize (const RCP< VectorBase< Scalar > > &col_vec) |
Initialize given a single vector object. More... | |
void | initialize (const RCP< const VectorSpaceBase< Scalar > > &range, const RCP< const VectorSpaceBase< Scalar > > &domain, const ArrayView< const RCP< VectorBase< Scalar > > > &col_vecs=Teuchos::null) |
Initialize given the spaces for the columns and rows and possibly the column vectors. More... | |
void | uninitialize () |
Set uninitialized. More... | |
Overridden from LinearOpBAse | |
RCP< const VectorSpaceBase < Scalar > > | range () const |
RCP< const VectorSpaceBase < Scalar > > | domain () const |
Overridden from MultiVectorBase | |
RCP< VectorBase< Scalar > > | nonconstColImpl (Ordinal j) |
RCP< MultiVectorBase< Scalar > > | nonconstContigSubViewImpl (const Range1D &col_rng) |
Overridden protected functions from MultiVectorBase | |
virtual void | assignImpl (Scalar alpha) |
virtual void | assignMultiVecImpl (const MultiVectorBase< Scalar > &mv) |
virtual void | scaleImpl (Scalar alpha) |
virtual void | updateImpl (Scalar alpha, const MultiVectorBase< Scalar > &mv) |
virtual void | linearCombinationImpl (const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta) |
virtual void | dotsImpl (const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const |
virtual void | norms1Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const |
virtual void | norms2Impl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const |
virtual void | normsInfImpl (const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const |
Overridden protected functions from LinearOpBase | |
bool | opSupportedImpl (EOpTransp M_trans) const |
For complex Scalar types returns true for NOTRANS and CONJTRANS and for real types returns true for all values of M_trans . More... | |
void | applyImpl (const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const |
This function is implemented in terms of the multi-vector applyOp() function. More... | |
Additional Inherited Members | |
Public Member Functions inherited from Thyra::MultiVectorDefaultBase< Scalar > | |
virtual RCP< MultiVectorBase < Scalar > > | clone_mv () const |
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... | |
Protected Member Functions inherited from Thyra::MultiVectorDefaultBase< Scalar > | |
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) |
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 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::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... | |
Default subclass for MultiVectorBase
implemented using columns of separate abstract vectors.
This is a very bad implementation of a multi-vector but this will work in situations where you need a multi-vector but some underlying linear algebra library does not directly support them.
This subclass can be used to represent a MultiVectorBase
wrapper around a single VectorBase
object so that a single vector can be passed to a method that expects a MultiVectorBase
object.
Definition at line 34 of file Thyra_DefaultColumnwiseMultiVector_decl.hpp.
Thyra::DefaultColumnwiseMultiVector< Scalar >::DefaultColumnwiseMultiVector | ( | ) |
Construct to initialized.
Postconditions:
this->range().get() == NULL
this->domain().get() == NULL
Definition at line 31 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
Thyra::DefaultColumnwiseMultiVector< Scalar >::DefaultColumnwiseMultiVector | ( | const RCP< VectorBase< Scalar > > & | col_vec | ) |
Calls initialize()
.
Definition at line 36 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
Thyra::DefaultColumnwiseMultiVector< Scalar >::DefaultColumnwiseMultiVector | ( | const RCP< const VectorSpaceBase< Scalar > > & | range, |
const RCP< const VectorSpaceBase< Scalar > > & | domain, | ||
const ArrayView< const RCP< VectorBase< Scalar > > > & | col_vecs = Teuchos::null |
||
) |
Calls initialize()
.
Definition at line 45 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
void Thyra::DefaultColumnwiseMultiVector< Scalar >::initialize | ( | const RCP< VectorBase< Scalar > > & | col_vec | ) |
Initialize given a single vector object.
col_vec | [in] A single column vector. It is not allowed for col_vecs==NULL . |
Preconditions:
col_vec.get() != NULL
(throw std::invalid_argument
) col_vec->dim() > 0
(throw std::invalid_argument
) Postconditions:
this->range().get() == col_vec.space().get()
this->domain()->dim() == 1
this->col(0).get() == col_vec.get()
Definition at line 56 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
void Thyra::DefaultColumnwiseMultiVector< Scalar >::initialize | ( | const RCP< const VectorSpaceBase< Scalar > > & | range, |
const RCP< const VectorSpaceBase< Scalar > > & | domain, | ||
const ArrayView< const RCP< VectorBase< Scalar > > > & | col_vecs = Teuchos::null |
||
) |
Initialize given the spaces for the columns and rows and possibly the column vectors.
range | [in] The space that the columns must lie in. The underlying vector space must not be changed while this is in use. |
domain | [in] The space that the rows must lie in. The underlying vector space must not be changed while this is in use. What this argument really specifies is what vector type will be compatible with the vectors that the client may try to use to interact with the rows of this multivector. |
col_vecs | [in] Array (size domain->dim() ) of column vectors to use for the columns of this . It is allowed for col_vecs==NULL in which case range->createMember() will be used to create the columns of this . |
Preconditions:
range.get() != NULL
(throw std::invalid_argument
) domain.get() != NULL
(throw std::invalid_argument
) range->dim() > 0
(throw std::invalid_argument
) domain->dim() > 0
(throw std::invalid_argument
) col_vecs != NULL
] col_vecs[j].get() != NULL && col_vecs[j]->space()->is_compatible(*range) == true
, for j=0..domain->dim()-1
Postconditions:
this->range().get() == range.get()
this->domain().get() == domain.get()
[col_vecs != NULL
] this->col(j).get() == col_vecs[j].get()
, for j=0..domain->dim()-1
Definition at line 74 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
void Thyra::DefaultColumnwiseMultiVector< Scalar >::uninitialize | ( | ) |
Set uninitialized.
Definition at line 106 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 119 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 127 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
virtual |
Implements Thyra::MultiVectorBase< Scalar >.
Definition at line 342 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
virtual |
Implements Thyra::MultiVectorBase< Scalar >.
Definition at line 356 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 137 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 147 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 162 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 171 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 187 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 213 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 230 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 244 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
Reimplemented from Thyra::MultiVectorDefaultBase< Scalar >.
Definition at line 258 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
For complex Scalar
types returns true
for NOTRANS
and CONJTRANS
and for real types returns true for all values of M_trans
.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 276 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.
|
protectedvirtual |
This function is implemented in terms of the multi-vector applyOp()
function.
The implementation takes care of two types of operations. One (M_trans==TRANS
) is the block dot product of two vectors to form scalar (stored as the vector y
which as one component). The other (M_trans==NOTRANS
) is essentially an axpy operation where x
is a vector with one element. Both of these operations are performed using reduction/transformation operators.
This implementation is near optimal but the default implementation of the multi-vector version of apply()
as implemented in the base class LinearOpBase
will not be a near optimal implementation in its current form do to multiple, sequential reductions but it could be made to be with a little work.
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 284 of file Thyra_DefaultColumnwiseMultiVector_def.hpp.