NOX
Development
|
Implementation of NOX::Abstract::MultiVector for Thyra multi-vectors. More...
#include <NOX_Thyra_MultiVector.H>
Public Member Functions | |
MultiVector (const Teuchos::RCP< ::Thyra::MultiVectorBase< double > > &source) | |
Constructor that creates a VIEW of the Thyra multivector. | |
MultiVector (const ::Thyra::MultiVectorBase< double > &source) | |
Construct from a given Thyra multivector. | |
MultiVector (const NOX::Thyra::MultiVector &source, NOX::CopyType type=NOX::DeepCopy) | |
Copy constructor. | |
~MultiVector () | |
Destruct MultiVector. | |
virtual NOX::size_type | length () const |
Return the length of multi-vector. | |
virtual int | numVectors () const |
Return the number of vectors in the multi-vector. | |
virtual void | print (std::ostream &stream) const |
Print the vector. This is meant for debugging purposes only. | |
void | setWeightVector (const Teuchos::RCP< const ::Thyra::VectorBase< double > > &weightVec) |
Set the weighting vector used for inner products and norms. | |
bool | hasWeightVector () const |
Teuchos::RCP< const ::Thyra::VectorBase< double > > | getWeightVector () const |
Returns a weighting vector if one was set, otherwise throws. | |
virtual Teuchos::RCP < ::Thyra::MultiVectorBase < double > > | getThyraMultiVector () |
Get RCP to underlying Thyra vector. | |
virtual Teuchos::RCP< const ::Thyra::MultiVectorBase < double > > | getThyraMultiVector () const |
Get const RCP to underlying Thyra vector. | |
virtual NOX::Abstract::MultiVector & | init (double value) |
Initialize every element of this multi-vector with gamma . | |
virtual NOX::Abstract::MultiVector & | random (bool useSeed=false, int seed=1) |
Initialize each element of this multi-vector with a random value. | |
virtual NOX::Abstract::MultiVector & | operator= (const NOX::Abstract::MultiVector &source) |
Copy source multi-vector source into this multi-vector. | |
NOX::Abstract::MultiVector & | operator= (const NOX::Thyra::MultiVector &source) |
virtual NOX::Abstract::MultiVector & | setBlock (const NOX::Abstract::MultiVector &source, const std::vector< int > &index) |
Copy the vectors in source to a set of vectors in *this . The index.size() vectors in source are copied to a subset of vectors in *this indicated by the indices given in index . | |
virtual NOX::Abstract::MultiVector & | augment (const NOX::Abstract::MultiVector &source) |
Append the vectors in source to *this . | |
virtual NOX::Abstract::Vector & | operator[] (int i) |
Return a reference to the i-th column of the multivector as an abstract vector. | |
virtual const NOX::Abstract::Vector & | operator[] (int i) const |
Return a const reference to the i-th column of the multivector as an abstract vector. | |
virtual NOX::Abstract::MultiVector & | scale (double gamma) |
Scale each element of this multivector by gamma . | |
virtual NOX::Abstract::MultiVector & | update (double alpha, const NOX::Abstract::MultiVector &a, double gamma=0.0) |
Compute x = (alpha * a) + (gamma * x) where a is a multi-vector and x = *this . | |
virtual NOX::Abstract::MultiVector & | update (double alpha, const NOX::Abstract::MultiVector &a, double beta, const NOX::Abstract::MultiVector &b, double gamma=0.0) |
Compute x = (alpha * a) + (beta * b) + (gamma * x) where a and b are multi-vectors and x = *this . | |
virtual NOX::Abstract::MultiVector & | update (Teuchos::ETransp transb, double alpha, const NOX::Abstract::MultiVector &a, const NOX::Abstract::MultiVector::DenseMatrix &b, double gamma=0.0) |
Compute x = (alpha * a * b) + (gamma * x) where a is a multivector, b is a dense matrix, x = *this , and op(b) = b if transb = Teuchos::NO_TRANS and op(b) is b transpose if transb = Teuchos::TRANS. | |
virtual Teuchos::RCP < NOX::Abstract::MultiVector > | clone (CopyType type=DeepCopy) const |
Create a new Vector of the same underlying type by cloning "this", and return a pointer to the new vector. More... | |
virtual Teuchos::RCP < NOX::Abstract::MultiVector > | clone (int numvecs) const |
Creates a new multi-vector with numvecs columns. | |
virtual Teuchos::RCP < NOX::Abstract::MultiVector > | subCopy (const std::vector< int > &index) const |
Creates a new multi-vector with index.size() columns whose columns are copies of the columns of *this given by index . | |
virtual Teuchos::RCP < NOX::Abstract::MultiVector > | subView (const std::vector< int > &index) const |
Creates a new multi-vector with ndex.size() columns that shares the columns of *this given by index . | |
virtual void | norm (std::vector< double > &result, NOX::Abstract::Vector::NormType type=NOX::Abstract::Vector::TwoNorm) const |
Norm. | |
virtual void | multiply (double alpha, const NOX::Abstract::MultiVector &y, NOX::Abstract::MultiVector::DenseMatrix &b) const |
Computes the matrix-matrix product . | |
bool | getImplicitWeighting () const |
Return true if implicit weighting is currently enabled. More... | |
void | setImplicitWeighting (bool do_implicit_weighting) |
Set to true to enable implicit weighting, false disables. More... | |
Public Member Functions inherited from NOX::Abstract::MultiVector | |
MultiVector () | |
Default constructor. Does nothing. | |
Protected Member Functions | |
bool | isContiguous (const std::vector< int > &index) const |
Check whether an index array is contiguous. | |
Protected Attributes | |
Teuchos::RCP < ::Thyra::MultiVectorBase < double > > | thyraMultiVec |
Pointer to petra vector owned by this object. | |
std::vector< Teuchos::RCP < NOX::Thyra::Vector > > | noxThyraVectors |
NOX::Thyra::Vector's for each column of the multivector. More... | |
Teuchos::RCP< const ::Thyra::VectorBase< double > > | weightVec_ |
Thyra vector used for weighting inner product and norms. | |
Teuchos::RCP < ::Thyra::MultiVectorBase < double > > | tmpMultiVec_ |
Thyra vector used for weighting inner product and norms. | |
bool | do_implicit_weighting_ |
True if implicit weighting is enabled (i.e. a nonnull wieghtVec_) | |
Additional Inherited Members | |
Public Types inherited from NOX::Abstract::MultiVector | |
typedef Teuchos::SerialDenseMatrix < int, double > | DenseMatrix |
Typename of dense matrices. | |
Implementation of NOX::Abstract::MultiVector for Thyra multi-vectors.
|
virtual |
Create a new Vector of the same underlying type by cloning "this", and return a pointer to the new vector.
If type is NOX::DeepCopy, then we need to create an exact replica of "this". Otherwise, if type is NOX::ShapeCopy, we need only replicate the shape of "this". Note that there is no assumption that a vector created by ShapeCopy is initialized to zeros.
Implements NOX::Abstract::MultiVector.
References Teuchos::rcp().
|
virtual |
Return true if implicit weighting is currently enabled.
CAUTION: This is a power user feature and should only be used in concert with specialized NOX::Abstract::Vector implementations.
Implements NOX::Abstract::ImplicitWeighting.
|
virtual |
Set to true to enable implicit weighting, false disables.
CAUTION: This is a power user feature and should only be used in concert with specialized NOX::Abstract::Vector implementations.
Implements NOX::Abstract::ImplicitWeighting.
Referenced by NOX::Thyra::Vector::createMultiVector().
|
mutableprotected |
NOX::Thyra::Vector's for each column of the multivector.
Each Thyra_Vector in the NOX::Thyra::Vector has a view into a column of the multivector and get filled in as needed by operator[].