NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | Protected Attributes | List of all members
NOX::Epetra::Vector Class Reference

Implementation of NOX::Abstract::Vector for Epetra vectors. More...

#include <NOX_Epetra_Vector.H>

Inheritance diagram for NOX::Epetra::Vector:
Inheritance graph
[legend]
Collaboration diagram for NOX::Epetra::Vector:
Collaboration graph
[legend]

Public Types

enum  MemoryType { CreateView, CreateCopy }
 Type of memory management to use when constructing the vector. More...
 
- Public Types inherited from NOX::Abstract::Vector
enum  NormType { TwoNorm, OneNorm, MaxNorm }
 Norm types used in norm() calculations. More...
 

Public Member Functions

 Vector (const Teuchos::RCP< Epetra_Vector > &source, NOX::Epetra::Vector::MemoryType memoryType=NOX::Epetra::Vector::CreateCopy, NOX::CopyType type=NOX::DeepCopy, Teuchos::RCP< NOX::Epetra::VectorSpace > vs=Teuchos::null)
 Constructor that creates a COPY or VIEW of the Epetra_Vector. More...
 
 Vector (const Epetra_Vector &source, NOX::CopyType type=NOX::DeepCopy, Teuchos::RCP< NOX::Epetra::VectorSpace > vs=Teuchos::null)
 Construct by copying map and/or elements of an Epetra_Vector. More...
 
 Vector (const NOX::Epetra::Vector &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
 ~Vector ()
 Destruct Vector.
 
virtual NOX::size_type length () const
 Return the length of vector. More...
 
virtual void print (std::ostream &stream) const
 Print the vector. To be used for debugging only.
 
virtual Teuchos::RCP
< NOX::Epetra::VectorSpace
getVectorSpace () const
 Returns the NOX::Epetra::VectorSpace associated with this vector.
 
virtual Epetra_VectorgetEpetraVector ()
 Get reference to underlying Epetra vector.
 
virtual const Epetra_VectorgetEpetraVector () const
 Get const reference to underlying Epetra vector.
 
virtual NOX::Abstract::Vectorinit (double gamma)
 Initialize every element of this vector with gamma. More...
 
virtual NOX::Abstract::Vectorrandom (bool useSeed=false, int seed=1)
 Initialize each element of this vector with a random value. More...
 
virtual NOX::Abstract::Vectoroperator= (const Epetra_Vector &y)
 Copies source vector into "this". More...
 
virtual NOX::Abstract::Vectoroperator= (const NOX::Epetra::Vector &y)
 
virtual NOX::Abstract::Vectoroperator= (const NOX::Abstract::Vector &y)
 Copy source vector y into this vector. More...
 
virtual NOX::Abstract::Vectorabs (const NOX::Epetra::Vector &y)
 
virtual NOX::Abstract::Vectorabs (const NOX::Abstract::Vector &y)
 Put element-wise absolute values of source vector y into this vector. More...
 
virtual NOX::Abstract::Vectorreciprocal (const NOX::Epetra::Vector &y)
 
virtual NOX::Abstract::Vectorreciprocal (const NOX::Abstract::Vector &y)
 Put element-wise reciprocal of source vector y into this vector. More...
 
virtual NOX::Abstract::Vectorscale (double gamma)
 Scale each element of this vector by gamma. More...
 
virtual NOX::Abstract::Vectorscale (const NOX::Epetra::Vector &a)
 
virtual NOX::Abstract::Vectorscale (const NOX::Abstract::Vector &a)
 Scale this vector element-by-element by the vector a. More...
 
virtual NOX::Abstract::Vectorupdate (double alpha, const NOX::Epetra::Vector &a, double gamma=0.0)
 
virtual NOX::Abstract::Vectorupdate (double alpha, const NOX::Abstract::Vector &a, double gamma=0.0)
 Compute x = (alpha * a) + (gamma * x) where x is this vector. More...
 
virtual NOX::Abstract::Vectorupdate (double alpha, const NOX::Epetra::Vector &a, double beta, const NOX::Epetra::Vector &b, double gamma=0.0)
 
virtual NOX::Abstract::Vectorupdate (double alpha, const NOX::Abstract::Vector &a, double beta, const NOX::Abstract::Vector &b, double gamma=0.0)
 Compute x = (alpha * a) + (beta * b) + (gamma * x) where x is this vector. More...
 
virtual Teuchos::RCP
< NOX::Abstract::Vector
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
createMultiVector (const NOX::Abstract::Vector *const *vecs, int numVecs, NOX::CopyType type=NOX::DeepCopy) const
 Create a MultiVector with numVecs+1 columns out of an array of Vectors. The vector stored under this will be the first column with the remaining numVecs columns given by vecs. More...
 
virtual Teuchos::RCP
< NOX::Abstract::MultiVector
createMultiVector (int numVecs, NOX::CopyType type=NOX::DeepCopy) const
 Create a MultiVector with numVecs columns. More...
 
virtual double norm (NOX::Abstract::Vector::NormType type=TwoNorm) const
 Norm. More...
 
virtual double norm (const NOX::Epetra::Vector &weights) const
 
virtual double norm (const NOX::Abstract::Vector &weights) const
 Weighted 2-Norm. More...
 
virtual double innerProduct (const NOX::Epetra::Vector &y) const
 
virtual double innerProduct (const NOX::Abstract::Vector &y) const
 Inner product with y. More...
 
- Public Member Functions inherited from NOX::Abstract::Vector
 Vector ()
 Abstract Vector constructor (does nothing)
 

Protected Attributes

Teuchos::RCP< Epetra_VectorepetraVec
 Pointer to petra vector owned by this object.
 
Teuchos::RCP
< NOX::Epetra::VectorSpace
vectorSpace
 Pointer to the vector space.
 

Detailed Description

Implementation of NOX::Abstract::Vector for Epetra vectors.

Member Enumeration Documentation

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.

Constructor & Destructor Documentation

NOX::Epetra::Vector::Vector ( const Teuchos::RCP< Epetra_Vector > &  source,
NOX::Epetra::Vector::MemoryType  memoryType = NOX::Epetra::Vector::CreateCopy,
NOX::CopyType  type = NOX::DeepCopy,
Teuchos::RCP< NOX::Epetra::VectorSpace vs = Teuchos::null 
)

Constructor that creates a COPY or VIEW of the Epetra_Vector.

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, epetraVec, Teuchos::is_null(), Teuchos::rcp(), NOX::ShapeCopy, and vectorSpace.

NOX::Epetra::Vector::Vector ( const Epetra_Vector source,
NOX::CopyType  type = NOX::DeepCopy,
Teuchos::RCP< NOX::Epetra::VectorSpace vs = Teuchos::null 
)

Construct by copying map and/or elements of an Epetra_Vector.

Allocates an entirely new vector. Does NOT allow for a view.

References NOX::DeepCopy, Teuchos::is_null(), Teuchos::rcp(), and NOX::ShapeCopy.

Member Function Documentation

NOX::Abstract::Vector & NOX::Epetra::Vector::abs ( const NOX::Abstract::Vector y)
virtual

Put element-wise absolute values of source vector y into this vector.

Here x represents this vector, and we update it as

\[ x_i = | y_i | \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

Teuchos::RCP< NOX::Abstract::Vector > NOX::Epetra::Vector::clone ( CopyType  type = DeepCopy) const
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" (the memory is allocated for the objects, but the current values are not copied into the vector). 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.

Implements NOX::Abstract::Vector.

References Teuchos::rcp().

Teuchos::RCP< NOX::Abstract::MultiVector > NOX::Epetra::Vector::createMultiVector ( const NOX::Abstract::Vector *const *  vecs,
int  numVecs,
NOX::CopyType  type = NOX::DeepCopy 
) const
virtual

Create a MultiVector with numVecs+1 columns out of an array of Vectors. The vector stored under this will be the first column with the remaining numVecs columns given by vecs.

The implementation here creates a NOX::Epetra::MultiVector with either Shape or Deep copies of the supplied vectors.

Reimplemented from NOX::Abstract::Vector.

References Epetra_Vector::ExtractView(), getEpetraVector(), Teuchos::rcp(), and View.

Teuchos::RCP< NOX::Abstract::MultiVector > NOX::Epetra::Vector::createMultiVector ( int  numVecs,
NOX::CopyType  type = NOX::DeepCopy 
) const
virtual

Create a MultiVector with numVecs columns.

The implementation here creates a NOX::Epetra::MultiVector with either Shape or Deep copies of the supplied vector.

Reimplemented from NOX::Abstract::Vector.

References Teuchos::rcp(), and NOX::ShapeCopy.

NOX::Abstract::Vector & NOX::Epetra::Vector::init ( double  gamma)
virtual

Initialize every element of this vector with gamma.

Here x represents this vector, and we update it as

\[ x_i = \gamma \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

Referenced by NOX::Epetra::LinearSystemAztecOO::applyJacobianInverse(), NOX::Epetra::LinearSystemAztecOO::applyRightPreconditioning(), NOX::Epetra::Group::computeNewton(), and NOX::Epetra::MatrixFree::MatrixFree().

double NOX::Epetra::Vector::innerProduct ( const NOX::Abstract::Vector y) const
virtual

Inner product with y.

Here x represents this vector, and we compute its inner product with y as follows:

\[ \langle x,y \rangle = \sum_{i=1}^n x_i y_i \]

Returns
$\langle x,y \rangle$

Implements NOX::Abstract::Vector.

NOX::size_type NOX::Epetra::Vector::length ( ) const
virtual

Return the length of vector.

Returns
The length of this vector
Note
Even if the vector is distributed across processors, this should return the global length of the vector.

Implements NOX::Abstract::Vector.

double NOX::Epetra::Vector::norm ( NOX::Abstract::Vector::NormType  type = TwoNorm) const
virtual
double NOX::Epetra::Vector::norm ( const NOX::Abstract::Vector weights) const
virtual

Weighted 2-Norm.

Here x represents this vector, and we compute its weighted norm as follows:

\[ \|x\|_w = \sqrt{\sum_{i=1}^{n} w_i \; x_i^2} \]

Returns
$ \|x\|_w $

Implements NOX::Abstract::Vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::operator= ( const Epetra_Vector y)
virtual

Copies source vector into "this".

NOTE: this will NOT copy the underlying vector space into the new vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::operator= ( const NOX::Abstract::Vector y)
virtual

Copy source vector y into this vector.

Here x represents this vector, and we update it as

\[ x_i = y_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::random ( bool  useSeed = false,
int  seed = 1 
)
virtual

Initialize each element of this vector with a random value.

If useSeed is true, uses the value of seed to seed the random number generator before filling the entries of this vector. So, if two calls are made where useSeed is true and seed is the same, then the vectors returned should be the same.

Default implementation throw an error. Only referenced by LOCA methods.

Returns
Reference to this object

Implements NOX::Abstract::Vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::reciprocal ( const NOX::Abstract::Vector y)
virtual

Put element-wise reciprocal of source vector y into this vector.

Here x represents this vector, and we update it as

\[ x_i = \frac{1}{y_i} \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::scale ( double  gamma)
virtual

Scale each element of this vector by gamma.

Here x represents this vector, and we update it as

\[ x_i = \gamma x_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

Referenced by NOX::Epetra::MatrixFree::Apply(), and NOX::Epetra::Group::computeNewton().

NOX::Abstract::Vector & NOX::Epetra::Vector::scale ( const NOX::Abstract::Vector a)
virtual

Scale this vector element-by-element by the vector a.

Here x represents this vector, and we update it as

\[ x_i = x_i \cdot a_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::update ( double  alpha,
const NOX::Abstract::Vector a,
double  gamma = 0.0 
)
virtual

Compute x = (alpha * a) + (gamma * x) where x is this vector.

Here x represents this vector, and we update it as

\[ x_i = \alpha \; a_i + \gamma \; x_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.

NOX::Abstract::Vector & NOX::Epetra::Vector::update ( double  alpha,
const NOX::Abstract::Vector a,
double  beta,
const NOX::Abstract::Vector b,
double  gamma = 0.0 
)
virtual

Compute x = (alpha * a) + (beta * b) + (gamma * x) where x is this vector.

Here x represents this vector, and we update it as

\[ x_i = \alpha \; a_i + \beta \; b_i + \gamma \; x_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implements NOX::Abstract::Vector.


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