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

Implementation of NOX::Abstract::Vector for STL std::vector<double> (using LAPACK for some computations) More...

#include <NOX_LAPACK_Vector.H>

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

Public Member Functions

 Vector ()
 Construct an empty vector.
 
 Vector (int n)
 Construct a zero vector of length n.
 
 Vector (int n, double *v)
 Construct a vector of length n from given array.
 
 Vector (const NOX::LAPACK::Vector &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
 ~Vector ()
 Destruct Vector.
 
NOX::size_type length () const
 Return the length of vector. More...
 
double & operator() (int i)
 Return the i-th element.
 
const double & operator() (int i) const
 Return the i-th element (const version)
 
std::ostream & leftshift (std::ostream &stream) const
 Prints out the vector to the specified stream. More...
 
void print (std::ostream &stream) const
 Print the vector. To be used for debugging only.
 
NOX::Abstract::Vectorinit (double value)
 Initialize every element of this vector with gamma. More...
 
virtual NOX::Abstract::Vectorrandom (bool useSeed=false, int seed=1)
 Initialize every element of this vector with random values.
 
NOX::Abstract::Vectoroperator= (const std::vector< double > &y)
 Replace this vector with STL double vector y.
 
NOX::Abstract::Vectoroperator= (const NOX::LAPACK::Vector &y)
 
NOX::Abstract::Vectoroperator= (const NOX::Abstract::Vector &y)
 Copy source vector y into this vector. More...
 
NOX::Abstract::Vectorabs (const NOX::LAPACK::Vector &y)
 
NOX::Abstract::Vectorabs (const NOX::Abstract::Vector &y)
 Put element-wise absolute values of source vector y into this vector. More...
 
NOX::Abstract::Vectorreciprocal (const NOX::LAPACK::Vector &y)
 
NOX::Abstract::Vectorreciprocal (const NOX::Abstract::Vector &y)
 Put element-wise reciprocal of source vector y into this vector. More...
 
NOX::Abstract::Vectorscale (double gamma)
 Scale each element of this vector by gamma. More...
 
NOX::Abstract::Vectorscale (const NOX::LAPACK::Vector &a)
 
NOX::Abstract::Vectorscale (const NOX::Abstract::Vector &a)
 Scale this vector element-by-element by the vector a. More...
 
NOX::Abstract::Vectorupdate (double alpha, const NOX::LAPACK::Vector &a, double gamma=0.0)
 
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...
 
NOX::Abstract::Vectorupdate (double alpha, const NOX::LAPACK::Vector &a, double beta, const NOX::LAPACK::Vector &b, double gamma=0.0)
 
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...
 
Teuchos::RCP
< NOX::Abstract::Vector
clone (NOX::CopyType type=NOX::DeepCopy) const
 Create a new Vector of the same underlying type by cloning "this", and return a pointer to the new vector. More...
 
double norm (NOX::Abstract::Vector::NormType type=NOX::Abstract::Vector::TwoNorm) const
 Norm. More...
 
double norm (const NOX::LAPACK::Vector &weights) const
 
double norm (const NOX::Abstract::Vector &weights) const
 Weighted 2-Norm. More...
 
double innerProduct (const NOX::LAPACK::Vector &y) const
 
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)
 
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...
 

Additional Inherited Members

- Public Types inherited from NOX::Abstract::Vector
enum  NormType { TwoNorm, OneNorm, MaxNorm }
 Norm types used in norm() calculations. More...
 

Detailed Description

Implementation of NOX::Abstract::Vector for STL std::vector<double> (using LAPACK for some computations)

Member Function Documentation

NOX::Abstract::Vector & NOX::LAPACK::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::LAPACK::Vector::clone ( NOX::CopyType  type = NOX::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().

NOX::Abstract::Vector & NOX::LAPACK::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.

double NOX::LAPACK::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.

std::ostream & NOX::LAPACK::Vector::leftshift ( std::ostream &  stream) const

Prints out the vector to the specified stream.

For example, a vector would appear as

\[ \left[ \; 0.1 \; 2.34 \; 5 \; \right] \]

It will be all on one line, with a single space between each entry, bracketed on either side.

Referenced by NOX::LAPACK::operator<<().

NOX::size_type NOX::LAPACK::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::LAPACK::Vector::norm ( NOX::Abstract::Vector::NormType  type = NOX::Abstract::Vector::TwoNorm) const
virtual

Norm.

Here x represents this vector, and we compute its norm as follows: for each NOX::Abstract::Vector::NormType:

Returns
$\|x\|$

Implements NOX::Abstract::Vector.

References NOX::LAPACK::i_one.

Referenced by LOCA::LAPACK::Interface::projectToDraw().

double NOX::LAPACK::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::LAPACK::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::LAPACK::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::LAPACK::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.

NOX::Abstract::Vector & NOX::LAPACK::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::LAPACK::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::LAPACK::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: