| 
    NOX
    Development
    
   | 
 
Implemenatation of the NOX::Abstract::MultiVector class for extended multi-vectors comprised of an arbitrary number of multi-vectors and scalars. More...
#include <LOCA_Extended_MultiVector.H>


Public Member Functions | |
| MultiVector (const MultiVector &source, NOX::CopyType type=NOX::DeepCopy) | |
| Copy constructor.  | |
| MultiVector (const MultiVector &source, int nColumns) | |
| Copy constructor that creates a new multivector with nColumns columns.  | |
| MultiVector (const MultiVector &source, const std::vector< int > &index, bool view) | |
| Copy constructor that creates a sub copy or view of the given multivector.  | |
| virtual | ~MultiVector () | 
| Vector destructor.  | |
| 
virtual  NOX::Abstract::MultiVector &  | init (double gamma) | 
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.  | |
| virtual MultiVector & | operator= (const MultiVector &y) | 
Copy source multi-vector source into this multi-vector.  | |
| 
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 &  | setBlock (const MultiVector &source, const std::vector< int > &index) | 
| 
virtual  NOX::Abstract::MultiVector &  | augment (const NOX::Abstract::MultiVector &source) | 
Append the vectors in source to *this.  | |
| 
virtual  NOX::Abstract::MultiVector &  | augment (const 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 a is a multi-vector and x = *this.  | |
| 
virtual  NOX::Abstract::MultiVector &  | update (double alpha, const 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 a and b are multi-vectors and x = *this.  | |
| 
virtual  NOX::Abstract::MultiVector &  | update (double alpha, const MultiVector &a, double beta, const 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 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  NOX::Abstract::MultiVector &  | update (Teuchos::ETransp transb, double alpha, const MultiVector &a, const NOX::Abstract::MultiVector::DenseMatrix &b, double gamma=0.0) | 
| virtual Teuchos::RCP < NOX::Abstract::MultiVector >  | clone (NOX::CopyType type=NOX::DeepCopy) const | 
| Create a new MultiVector 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 index.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  .  | |
| virtual void | multiply (double alpha, const MultiVector &y, NOX::Abstract::MultiVector::DenseMatrix &b) const | 
| 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 Teuchos::RCP< const  NOX::Abstract::MultiVector >  | getMultiVector (int i) const | 
| Returns const ref-count pointer to the ith multi-vector.  | |
| 
virtual Teuchos::RCP < NOX::Abstract::MultiVector >  | getMultiVector (int i) | 
| Returns ref-count pointer to the ith multi-vector.  | |
| 
virtual Teuchos::RCP< const  NOX::Abstract::MultiVector::DenseMatrix >  | getScalars () const | 
| Returns const ref-count pointer to scalar matrix.  | |
| 
virtual Teuchos::RCP < NOX::Abstract::MultiVector::DenseMatrix >  | getScalars () | 
| Returns ref-count pointer to scalar matrix.  | |
| 
virtual Teuchos::RCP< const  NOX::Abstract::MultiVector::DenseMatrix >  | getScalarRows (int num_rows, int row) const | 
Returns const ref-count pointer to num_rows rows of scalar matrix starting at row row.  | |
| 
virtual Teuchos::RCP < NOX::Abstract::MultiVector::DenseMatrix >  | getScalarRows (int num_rows, int row) | 
Returns ref-count pointer to num_rows rows of scalar matrix starting at row row.  | |
| virtual const double & | getScalar (int i, int j) const | 
| Returns const reference to the scalar for row i, column j.  | |
| virtual double & | getScalar (int i, int j) | 
| Returns reference to the scalar for row i, column j.  | |
| 
virtual Teuchos::RCP < LOCA::Extended::Vector >  | getVector (int i) | 
| Return a ref-count pointer to the i-th column of the multivector as an abstract vector.  | |
| 
virtual Teuchos::RCP< const  LOCA::Extended::Vector >  | getVector (int i) const | 
| Return a const ref-count pointer to the i-th column of the multivector as an abstract vector.  | |
| virtual int | getNumScalarRows () const | 
| Returns number of scalars rows.  | |
| virtual int | getNumMultiVectors () const | 
| Returns number of multi vectors.  | |
  Public Member Functions inherited from NOX::Abstract::MultiVector | |
| MultiVector () | |
| Default constructor. Does nothing.  | |
Protected Member Functions | |
| MultiVector (const Teuchos::RCP< LOCA::GlobalData > &global_data, int nColumns, int nVectorRows, int nScalarRows) | |
| Constructor that creates an empty multivector to be filled in later.  | |
| virtual Teuchos::RCP < LOCA::Extended::Vector >  | generateVector (int nVecs, int nScalarRows) const | 
| Generate a derived extended vector.  More... | |
| void | setMultiVectorPtr (int i, Teuchos::RCP< NOX::Abstract::MultiVector > v) | 
| Sets the pointer to the ith multivector.  | |
| void | checkDimensions (const std::string &callingFunction, const LOCA::Extended::MultiVector &a) const | 
| Checks multi-vec argument dimensions are consistent.  | |
| void | checkIndex (const std::string &callingFunction, int i) const | 
| Checks validity of column index.  | |
| void | checkVectorRowIndex (const std::string &callingFunction, int i) const | 
| Checks validity of vector row index.  | |
| void | checkIndex (const std::string &callingFunction, int i, int j) const | 
| Checks validity of column and row index for scalars.  | |
| bool | isContiguous (const std::vector< int > &index) const | 
| Checks is index array is contiguous.  | |
Protected Attributes | |
| Teuchos::RCP< LOCA::GlobalData > | globalData | 
| Global data.  | |
| int | numColumns | 
| Number of columns in each multivec and number of scalar vector columns.  | |
| int | numMultiVecRows | 
| Number of multivec block rows.  | |
| int | numScalarRows | 
| Number of scalar rows.  | |
| 
std::vector< Teuchos::RCP < NOX::Abstract::MultiVector > >  | multiVectorPtrs | 
| Array of multi-vector pointers, one for each block ROW.  | |
| 
Teuchos::RCP < NOX::Abstract::MultiVector::DenseMatrix >  | scalarsPtr | 
| Dense matrix of scalars.  | |
| 
std::vector< Teuchos::RCP < LOCA::Extended::Vector > >  | extendedVectorPtrs | 
| Pointers to each column as a LOCA::Extended::Vector.  | |
| bool | isView | 
| Flag indicating whether this vector is a view.  | |
Friends | |
| class | Vector | 
| Declare LOCA::Extended::Vector as a friend class.  | |
Additional Inherited Members | |
  Public Types inherited from NOX::Abstract::MultiVector | |
| 
typedef  Teuchos::SerialDenseMatrix < int, double >  | DenseMatrix | 
| Typename of dense matrices.  | |
Implemenatation of the NOX::Abstract::MultiVector class for extended multi-vectors comprised of an arbitrary number of multi-vectors and scalars.
      
  | 
  virtual | 
Create a new MultiVector 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.
Reimplemented in LOCA::MultiContinuation::ExtendedMultiVector, LOCA::Pitchfork::MooreSpence::ExtendedMultiVector, LOCA::TurningPoint::MooreSpence::ExtendedMultiVector, LOCA::Hopf::MooreSpence::ExtendedMultiVector, LOCA::PhaseTransition::ExtendedMultiVector, and LOCA::Hopf::ComplexMultiVector.
References Teuchos::rcp().
      
  | 
  protectedvirtual | 
Generate a derived extended vector.
Extended multi-vectors derived from this class should implement this method and return a vector of the appropriate type. This allows the operator[] methods to work correctly for derived classes.
Reimplemented in LOCA::Hopf::MooreSpence::ExtendedMultiVector, LOCA::Pitchfork::MooreSpence::ExtendedMultiVector, LOCA::TurningPoint::MooreSpence::ExtendedMultiVector, LOCA::PhaseTransition::ExtendedMultiVector, LOCA::Hopf::ComplexMultiVector, and LOCA::MultiContinuation::ExtendedMultiVector.
References Teuchos::rcp().
 1.8.5