| NOX
    Development
    | 
Implementation of NOX::Abstract::MultiVector for Epetra multi-vectors. More...
#include <NOX_Epetra_MultiVector.H>


| Public Types | |
| enum | MemoryType { CreateView, CreateCopy } | 
| Type of memory management to use when constructing the vector.  More... | |
|  Public Types inherited from NOX::Abstract::MultiVector | |
| typedef Teuchos::SerialDenseMatrix < int, double > | DenseMatrix | 
| Typename of dense matrices. | |
| Public Member Functions | |
| MultiVector (const Teuchos::RCP< Epetra_MultiVector > &source, NOX::CopyType type=NOX::DeepCopy, NOX::Epetra::MultiVector::MemoryType memoryType=NOX::Epetra::MultiVector::CreateCopy) | |
| Constructor that creates a COPY or VIEW of the Epetra_MultiVector.  More... | |
| MultiVector (const Epetra_MultiVector &source, NOX::CopyType type=NOX::DeepCopy) | |
| Construct by copying map and/or elements of an Epetra_MultiVector. | |
| MultiVector (const NOX::Epetra::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. | |
| virtual Epetra_MultiVector & | getEpetraMultiVector () | 
| Get reference to underlying Epetra vector. | |
| virtual const Epetra_MultiVector & | getEpetraMultiVector () const | 
| Get const reference to underlying Epetra 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 Epetra_MultiVector &source) | 
| Copy source multi-vector sourceinto this multi-vector. | |
| virtual NOX::Abstract::MultiVector & | operator= (const NOX::Epetra::MultiVector &source) | 
| Copy source multi-vector sourceinto this multi-vector. | |
| virtual NOX::Abstract::MultiVector & | operator= (const NOX::Abstract::MultiVector &source) | 
| Copy source multi-vector sourceinto this multi-vector. | |
| virtual NOX::Abstract::MultiVector & | setBlock (const NOX::Abstract::MultiVector &source, const std::vector< int > &index) | 
| Copy the vectors in sourceto a set of vectors in*this. Theindex.size()vectors insourceare copied to a subset of vectors in*thisindicated by the indices given inindex. | |
| virtual NOX::Abstract::MultiVector & | setBlock (const NOX::Epetra::MultiVector &source, const std::vector< int > &index) | 
| virtual NOX::Abstract::MultiVector & | augment (const NOX::Abstract::MultiVector &source) | 
| Append the vectors in sourceto*this. | |
| virtual NOX::Abstract::MultiVector & | augment (const NOX::Epetra::MultiVector &source) | 
| 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 ais a multi-vector andx=*this. | |
| virtual NOX::Abstract::MultiVector & | update (double alpha, const NOX::Epetra::MultiVector &a, double gamma=0.0) | 
| 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 aandbare multi-vectors andx=*this. | |
| virtual NOX::Abstract::MultiVector & | update (double alpha, const NOX::Epetra::MultiVector &a, double beta, const NOX::Epetra::MultiVector &b, double gamma=0.0) | 
| 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 ais a multivector,bis 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 NOX::Abstract::MultiVector & | update (Teuchos::ETransp transb, double alpha, const NOX::Epetra::MultiVector &a, const NOX::Abstract::MultiVector::DenseMatrix &b, double gamma=0.0) | 
| 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 numvecscolumns. | |
| 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*thisgiven byindex. | |
| 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*thisgiven byindex. | |
| 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  . | |
| virtual void | multiply (double alpha, const NOX::Epetra::MultiVector &y, NOX::Abstract::MultiVector::DenseMatrix &b) const | 
|  Public Member Functions inherited from NOX::Abstract::MultiVector | |
| MultiVector () | |
| Default constructor. Does nothing. | |
| Protected Member Functions | |
| MultiVector (int numvecs) | |
| Constructor (Protected) | |
| void | checkIndex (int idx) const | 
| Checks whether an index is valid. Throws an error if invalid. | |
| Protected Attributes | |
| Teuchos::RCP< Epetra_MultiVector > | epetraMultiVec | 
| Pointer to petra vector owned by this object. | |
| std::vector < NOX::Epetra::Vector * > | noxEpetraVectors | 
| NOX::Epetra::Vector's for each column of the multivector.  More... | |
Implementation of NOX::Abstract::MultiVector for Epetra multi-vectors.
Type of memory management to use when constructing the vector.
| Enumerator | |
|---|---|
| CreateView | Keeps a pointer to and uses the actual Epetra_Vector passed in. | 
| CreateCopy | Allocates a new underlying Epetra_Vector object. | 
| NOX::Epetra::MultiVector::MultiVector | ( | const Teuchos::RCP< Epetra_MultiVector > & | source, | 
| NOX::CopyType | type = NOX::DeepCopy, | ||
| NOX::Epetra::MultiVector::MemoryType | memoryType = NOX::Epetra::MultiVector::CreateCopy | ||
| ) | 
Constructor that creates a COPY or VIEW of the Epetra_MultiVector.
NOTE: This ctor should just always create a view. It should be implicit from the fact that a RCP object is being passed in that a persisting relationship is present. However, since this could cause confusion, the default is to make a copy and if a user wants a view, they must pass in an explicit flag.
A VIEW of a vector uses the same underlying memory. WARNING: A View can be dangerous since multiple objects can access the same memory locations.
References CreateView, NOX::DeepCopy, epetraMultiVec, noxEpetraVectors, Teuchos::rcp(), and NOX::ShapeCopy.
| 
 | 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().
| 
 | mutableprotected | 
NOX::Epetra::Vector's for each column of the multivector.
Each Epetra_Vector in the NOX::Epetra::Vector has a view into a column of the multivector and get filled in as needed by operator[].
Referenced by MultiVector().
 1.8.5
 1.8.5