NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
NOX::Abstract::MultiVector Class Referenceabstract

Abstract interface for multi-vectors used by NOX. More...

#include <NOX_Abstract_MultiVector.H>

Inheritance diagram for NOX::Abstract::MultiVector:
Inheritance graph
[legend]

Public Types

typedef
Teuchos::SerialDenseMatrix
< int, double > 
DenseMatrix
 Typename of dense matrices.
 

Public Member Functions

 MultiVector ()
 Default constructor. Does nothing.
 
virtual ~MultiVector ()
 Destructor. Does nothing.
 
virtual NOX::size_type length () const =0
 Return the length of multi-vector.
 
virtual int numVectors () const =0
 Return the number of vectors in the multi-vector.
 
virtual void print (std::ostream &stream) const =0
 Print the vector. This is meant for debugging purposes only.
 
virtual
NOX::Abstract::MultiVector
init (double gamma)=0
 Initialize every element of this multi-vector with gamma.
 
virtual
NOX::Abstract::MultiVector
random (bool useSeed=false, int seed=1)=0
 Initialize each element of this multi-vector with a random value.
 
virtual
NOX::Abstract::MultiVector
operator= (const NOX::Abstract::MultiVector &source)=0
 Copy source multi-vector source into this multi-vector.
 
virtual
NOX::Abstract::MultiVector
setBlock (const NOX::Abstract::MultiVector &source, const std::vector< int > &index)=0
 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)=0
 Append the vectors in source to *this.
 
virtual NOX::Abstract::Vectoroperator[] (int i)=0
 Return a reference to the i-th column of the multivector as an abstract vector. More...
 
virtual const
NOX::Abstract::Vector
operator[] (int i) const =0
 Return a const reference to the i-th column of the multivector as an abstract vector. More...
 
virtual
NOX::Abstract::MultiVector
scale (double gamma)=0
 Scale each element of this multivector by gamma.
 
virtual
NOX::Abstract::MultiVector
update (double alpha, const NOX::Abstract::MultiVector &a, double gamma=0.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)=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 DenseMatrix &b, double gamma=0.0)=0
 Compute x = (alpha * a * op(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 (NOX::CopyType type=NOX::DeepCopy) const =0
 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 =0
 Creates a new multi-vector with numvecs columns.
 
virtual Teuchos::RCP
< NOX::Abstract::MultiVector
subCopy (const std::vector< int > &index) const =0
 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 =0
 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 =0
 Norm.
 
virtual void multiply (double alpha, const NOX::Abstract::MultiVector &y, DenseMatrix &b) const =0
 Computes the matrix-matrix product $\alpha * y^T * (*this)$.
 

Detailed Description

Abstract interface for multi-vectors used by NOX.

Member Function Documentation

virtual Teuchos::RCP<NOX::Abstract::MultiVector> NOX::Abstract::MultiVector::clone ( NOX::CopyType  type = NOX::DeepCopy) const
pure 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.

Returns
Pointer to newly created vector or NULL if clone is not supported.

Implemented in NOX::Epetra::MultiVector, LOCA::Extended::MultiVector, NOX::Thyra::MultiVector, NOX::MultiVector, LOCA::MultiContinuation::ExtendedMultiVector, LOCA::Pitchfork::MooreSpence::ExtendedMultiVector, LOCA::TurningPoint::MooreSpence::ExtendedMultiVector, LOCA::Hopf::MooreSpence::ExtendedMultiVector, LOCA::PhaseTransition::ExtendedMultiVector, and LOCA::Hopf::ComplexMultiVector.

Referenced by LOCA::AnasaziOperator::ShiftInvert::apply(), LOCA::AnasaziOperator::ShiftInvert2Matrix::apply(), LOCA::AnasaziOperator::Cayley2Matrix::apply(), LOCA::Epetra::AnasaziOperator::Floquet::apply(), LOCA::AnasaziOperator::Cayley::apply(), LOCA::BorderedSolver::Bordering::applyInverse(), LOCA::BorderedSolver::Bordering::applyInverseTranspose(), LOCA::Hopf::MooreSpence::ExtendedGroup::applyJacobianMultiVector(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::applyJacobianMultiVector(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianMultiVector(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianTransposeMultiVector(), Anasazi::MultiVecTraits< double, NOX::Abstract::MultiVector >::Clone(), Anasazi::MultiVecTraits< double, NOX::Abstract::MultiVector >::CloneCopy(), LOCA::Hopf::ComplexMultiVector::ComplexMultiVector(), LOCA::MultiContinuation::CompositeConstraintMVDX::CompositeConstraintMVDX(), LOCA::TurningPoint::MinimallyAugmented::ModifiedConstraint::computeConstraints(), LOCA::Thyra::Group::computeDfDpMulti(), LOCA::Hopf::MinimallyAugmented::Constraint::computeDOmega(), LOCA::DerivUtils::computeDwtJnDx(), LOCA::BorderedSolver::HouseholderQR::computeQR(), LOCA::Hopf::MinimallyAugmented::Constraint::Constraint(), LOCA::MultiContinuation::CompositeConstraintMVDX::copy(), LOCA::Hopf::MooreSpence::ExtendedGroup::ExtendedGroup(), LOCA::PhaseTransition::ExtendedMultiVector::ExtendedMultiVector(), LOCA::MultiContinuation::ExtendedMultiVector::ExtendedMultiVector(), LOCA::Hopf::MooreSpence::ExtendedMultiVector::ExtendedMultiVector(), LOCA::Pitchfork::MooreSpence::ExtendedMultiVector::ExtendedMultiVector(), LOCA::TurningPoint::MooreSpence::ExtendedMultiVector::ExtendedMultiVector(), LOCA::MultiPredictor::Tangent::operator=(), LOCA::AnasaziOperator::Cayley::preProcessSeedVector(), LOCA::AnasaziOperator::Cayley2Matrix::preProcessSeedVector(), LOCA::Epetra::AnasaziOperator::Floquet::rayleighQuotient(), LOCA::BorderedSolver::UpperTriangularBlockElimination::solve(), LOCA::TurningPoint::MooreSpence::SalingerBordering::solve(), LOCA::Pitchfork::MooreSpence::SalingerBordering::solve(), LOCA::TurningPoint::MooreSpence::PhippsBordering::solve(), LOCA::Hopf::MooreSpence::SalingerBordering::solve(), LOCA::Pitchfork::MooreSpence::PhippsBordering::solve(), LOCA::BorderedSolver::TpetraHouseholder::solve(), LOCA::BorderedSolver::EpetraHouseholder::solve(), LOCA::Pitchfork::MooreSpence::SalingerBordering::solveContiguous(), LOCA::TurningPoint::MooreSpence::SalingerBordering::solveContiguous(), LOCA::Hopf::MooreSpence::SalingerBordering::solveContiguous(), LOCA::Pitchfork::MooreSpence::PhippsBordering::solveContiguous(), LOCA::TurningPoint::MooreSpence::PhippsBordering::solveContiguous(), LOCA::BorderedSolver::Bordering::solveFZero(), LOCA::BorderedSolver::Bordering::solveFZeroTrans(), LOCA::BorderedSolver::UpperTriangularBlockElimination::solveTranspose(), LOCA::TurningPoint::MooreSpence::SalingerBordering::solveTranspose(), LOCA::TurningPoint::MooreSpence::PhippsBordering::solveTranspose(), LOCA::BorderedSolver::EpetraHouseholder::solveTranspose(), LOCA::TurningPoint::MooreSpence::SalingerBordering::solveTransposeContiguous(), LOCA::TurningPoint::MooreSpence::PhippsBordering::solveTransposeContiguous(), and LOCA::MultiPredictor::Tangent::Tangent().

virtual NOX::Abstract::Vector& NOX::Abstract::MultiVector::operator[] ( int  i)
pure virtual

Return a reference to the i-th column of the multivector as an abstract vector.

It is assumed that any implementation of this method returns a vector that has a view of the underlying data in the multivector

Implemented in NOX::Epetra::MultiVector, LOCA::Extended::MultiVector, NOX::Thyra::MultiVector, and NOX::MultiVector.

virtual const NOX::Abstract::Vector& NOX::Abstract::MultiVector::operator[] ( int  i) const
pure virtual

Return a const reference to the i-th column of the multivector as an abstract vector.

It is assumed that any implementation of this method returns a vector that has a view of the underlying data in the multivector.

Implemented in NOX::Epetra::MultiVector, LOCA::Extended::MultiVector, NOX::Thyra::MultiVector, and NOX::MultiVector.


The documentation for this class was generated from the following file: