Teuchos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Teuchos::SerialDenseVector< OrdinalType, ScalarType > Class Template Reference

This class creates and provides basic support for dense vectors of templated type as a specialization of Teuchos::SerialDenseMatrix. Additional methods for the SerialDenseVector class, like mathematical methods, can be found documented in SerialDenseMatrix. More...

#include <Teuchos_SerialDenseVector.hpp>

Inheritance diagram for Teuchos::SerialDenseVector< OrdinalType, ScalarType >:
Inheritance graph
[legend]

Constructor/Destructor methods.

 SerialDenseVector ()
 Default Constructor. More...
 
 SerialDenseVector (OrdinalType length, bool zeroOut=true)
 Shaped Constructor. More...
 
 SerialDenseVector (DataAccess CV, ScalarType *values, OrdinalType length)
 Shaped Constructor with Values. More...
 
 SerialDenseVector (const SerialDenseVector< OrdinalType, ScalarType > &Source)
 Copy Constructor. More...
 
 SerialDenseVector (DataAccess CV, const SerialDenseVector< OrdinalType, ScalarType > &Source)
 Copy Constructor. More...
 
virtual ~SerialDenseVector ()
 Destructor. More...
 

Sizing methods.

int size (OrdinalType length_in)
 Size method for changing the size of a SerialDenseVector, initializing entries to zero. More...
 
int sizeUninitialized (OrdinalType length_in)
 Same as size() except leaves values uninitialized. More...
 
int resize (OrdinalType length_in)
 Resizing method for changing the size of a SerialDenseVector, keeping the entries. More...
 

Set methods.

SerialDenseVector< OrdinalType,
ScalarType > & 
operator= (const ScalarType value)
 Set all values in the matrix to a constant value. More...
 
SerialDenseVector< OrdinalType,
ScalarType > & 
operator= (const SerialDenseVector< OrdinalType, ScalarType > &Source)
 Copies values from one vector to another. More...
 

Comparison methods.

bool operator== (const SerialDenseVector< OrdinalType, ScalarType > &Operand) const
 Equality of two matrices. More...
 
bool operator!= (const SerialDenseVector< OrdinalType, ScalarType > &Operand) const
 Inequality of two matrices. More...
 

Accessor methods.

ScalarType & operator() (OrdinalType index)
 Element access method (non-const). More...
 
const ScalarType & operator() (OrdinalType index) const
 Element access method (const). More...
 
ScalarType & operator[] (OrdinalType index)
 Element access method (non-const). More...
 
const ScalarType & operator[] (OrdinalType index) const
 Element access method (const). More...
 

Mathematical methods.

ScalarType dot (const SerialDenseVector< OrdinalType, ScalarType > &x) const
 Compute the dot product of this vector and x. More...
 

Attribute methods.

OrdinalType length () const
 Returns the length of this vector. More...
 

I/O methods.

std::ostream & print (std::ostream &os) const
 Print method. Define the behavior of the std::ostream << operator inherited from the Object class. More...
 

Additional Inherited Members

- Public Types inherited from Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >
typedef OrdinalType ordinalType
 Typedef for ordinal type. More...
 
typedef ScalarType scalarType
 Typedef for scalar type. More...
 
- Public Types inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType >
typedef details::GivensRotator
< ScalarType >::c_type 
rotg_c_type
 The type used for c in ROTG. More...
 
- Public Member Functions inherited from Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >
 SerialDenseMatrix ()=default
 Default Constructor. More...
 
 SerialDenseMatrix (OrdinalType numRows, OrdinalType numCols, bool zeroOut=true)
 Shaped Constructor. More...
 
 SerialDenseMatrix (DataAccess CV, ScalarType *values, OrdinalType stride, OrdinalType numRows, OrdinalType numCols)
 Shaped Constructor with Values. More...
 
 SerialDenseMatrix (const SerialDenseMatrix< OrdinalType, ScalarType > &Source, ETransp trans=Teuchos::NO_TRANS)
 Copy Constructor. More...
 
 SerialDenseMatrix (DataAccess CV, const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
 Copy Constructor. More...
 
 SerialDenseMatrix (DataAccess CV, const SerialDenseMatrix< OrdinalType, ScalarType > &Source, OrdinalType numRows, OrdinalType numCols, OrdinalType startRow=0, OrdinalType startCol=0)
 Submatrix Copy Constructor. More...
 
virtual ~SerialDenseMatrix ()
 Destructor. More...
 
int shape (OrdinalType numRows, OrdinalType numCols)
 Shape method for changing the size of a SerialDenseMatrix, initializing entries to zero. More...
 
int shapeUninitialized (OrdinalType numRows, OrdinalType numCols)
 Same as shape() except leaves uninitialized. More...
 
int reshape (OrdinalType numRows, OrdinalType numCols)
 Reshaping method for changing the size of a SerialDenseMatrix, keeping the entries. More...
 
SerialDenseMatrix< OrdinalType,
ScalarType > & 
operator= (const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
 Copies values from one matrix to another. More...
 
SerialDenseMatrix< OrdinalType,
ScalarType > & 
assign (const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
 Copies values from one matrix to another. More...
 
SerialDenseMatrix< OrdinalType,
ScalarType > & 
operator= (const ScalarType value)
 Set all values in the matrix to a constant value. More...
 
int putScalar (const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
 Set all values in the matrix to a constant value. More...
 
void swap (SerialDenseMatrix< OrdinalType, ScalarType > &B)
 Swap values between this matrix and incoming matrix. More...
 
int random ()
 Set all values in the matrix to be random numbers. More...
 
ScalarType & operator() (OrdinalType rowIndex, OrdinalType colIndex)
 Element access method (non-const). More...
 
const ScalarType & operator() (OrdinalType rowIndex, OrdinalType colIndex) const
 Element access method (const). More...
 
ScalarType * operator[] (OrdinalType colIndex)
 Column access method (non-const). More...
 
const ScalarType * operator[] (OrdinalType colIndex) const
 Column access method (const). More...
 
ScalarType * values () const
 Data array access method. More...
 
SerialDenseMatrix< OrdinalType,
ScalarType > & 
operator+= (const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
 Add another matrix to this matrix. More...
 
SerialDenseMatrix< OrdinalType,
ScalarType > & 
operator-= (const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
 Subtract another matrix from this matrix. More...
 
SerialDenseMatrix< OrdinalType,
ScalarType > & 
operator*= (const ScalarType alpha)
 Scale this matrix by alpha; *this = alpha**this. More...
 
int scale (const ScalarType alpha)
 Scale this matrix by alpha; *this = alpha**this. More...
 
int scale (const SerialDenseMatrix< OrdinalType, ScalarType > &A)
 Point-wise scale this matrix by A; i.e. *this(i,j) *= A(i,j) More...
 
int multiply (ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
 Multiply A * B and add them to this; this = beta * this + alpha*A*B. More...
 
int multiply (ESide sideA, ScalarType alpha, const SerialSymDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
 Multiply A and B and add them to this; this = beta * this + alpha*A*B or this = beta * this + alpha*B*A. More...
 
bool operator== (const SerialDenseMatrix< OrdinalType, ScalarType > &Operand) const
 Equality of two matrices. More...
 
bool operator!= (const SerialDenseMatrix< OrdinalType, ScalarType > &Operand) const
 Inequality of two matrices. More...
 
OrdinalType numRows () const
 Returns the row dimension of this matrix. More...
 
OrdinalType numCols () const
 Returns the column dimension of this matrix. More...
 
OrdinalType stride () const
 Returns the stride between the columns of this matrix in memory. More...
 
bool empty () const
 Returns whether this matrix is empty. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
normOne () const
 Returns the 1-norm of the matrix. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
normInf () const
 Returns the Infinity-norm of the matrix. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
normFrobenius () const
 Returns the Frobenius-norm of the matrix. More...
 
- Public Member Functions inherited from Teuchos::CompObject
 CompObject ()
 Default constructor. More...
 
 CompObject (const CompObject &source)
 Copy Constructor. More...
 
virtual ~CompObject ()
 Destructor. More...
 
void setFlopCounter (const Flops &FlopCounter)
 Set the internal Teuchos::Flops() pointer. More...
 
void setFlopCounter (const CompObject &compObject)
 Set the internal Teuchos::Flops() pointer to the flop counter of another Teuchos::CompObject. More...
 
void unsetFlopCounter ()
 Set the internal Teuchos::Flops() pointer to 0 (no flops counted). More...
 
FlopsgetFlopCounter () const
 Get the pointer to the Teuchos::Flops() object associated with this object, returns 0 if none. More...
 
void resetFlops () const
 Resets the number of floating point operations to zero for this multi-std::vector. More...
 
double getFlops () const
 Returns the number of floating point operations with this multi-std::vector. More...
 
void updateFlops (int addflops) const
 Increment Flop count for this object. More...
 
void updateFlops (long int addflops) const
 Increment Flop count for this object. More...
 
void updateFlops (double addflops) const
 Increment Flop count for this object. More...
 
void updateFlops (float addflops) const
 Increment Flop count for this object. More...
 
- Public Member Functions inherited from Teuchos::BLAS< OrdinalType, ScalarType >
 BLAS (void)
 Default constructor. More...
 
 BLAS (const BLAS< OrdinalType, ScalarType > &)
 Copy constructor. More...
 
virtual ~BLAS (void)
 Destructor. More...
 
- Public Member Functions inherited from Teuchos::DefaultBLASImpl< OrdinalType, ScalarType >
 DefaultBLASImpl (void)
 Default constructor. More...
 
 DefaultBLASImpl (const DefaultBLASImpl< OrdinalType, ScalarType > &)
 Copy constructor. More...
 
virtual ~DefaultBLASImpl (void)
 Destructor. More...
 
void ROTG (ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
 Computes a Givens plane rotation. More...
 
void ROT (const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
 Applies a Givens plane rotation. More...
 
void SCAL (const OrdinalType &n, const ScalarType &alpha, ScalarType *x, const OrdinalType &incx) const
 Scale the vector x by the constant alpha. More...
 
void COPY (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
 Copy the vector x to the vector y. More...
 
template<typename alpha_type , typename x_type >
void AXPY (const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
 Perform the operation: y <- y+alpha*x. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
ASUM (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Sum the absolute values of the entries of x. More...
 
template<typename x_type , typename y_type >
ScalarType DOT (const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const
 Form the dot product of the vectors x and y. More...
 
ScalarTraits< ScalarType >
::magnitudeType 
NRM2 (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Compute the 2-norm of the vector x. More...
 
OrdinalType IAMAX (const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
 Return the index of the element of x with the maximum magnitude. More...
 
template<typename alpha_type , typename A_type , typename x_type , typename beta_type >
void GEMV (ETransp trans, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const x_type *x, const OrdinalType &incx, const beta_type beta, ScalarType *y, const OrdinalType &incy) const
 Performs the matrix-vector operation: y <- alpha*A*x+beta*y or y <- alpha*A'*x+beta*y where A is a general m by n matrix. More...
 
template<typename A_type >
void TRMV (EUplo uplo, ETransp trans, EDiag diag, const OrdinalType &n, const A_type *A, const OrdinalType &lda, ScalarType *x, const OrdinalType &incx) const
 Performs the matrix-vector operation: x <- A*x or x <- A'*x where A is a unit/non-unit n by n upper/lower triangular matrix. More...
 
template<typename alpha_type , typename x_type , typename y_type >
void GER (const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy, ScalarType *A, const OrdinalType &lda) const
 Performs the rank 1 operation: A <- alpha*x*y'+A. More...
 
template<typename alpha_type , typename A_type , typename B_type , typename beta_type >
void GEMM (ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 General matrix-matrix multiply. More...
 
void SWAP (const OrdinalType &n, ScalarType *const x, const OrdinalType &incx, ScalarType *const y, const OrdinalType &incy) const
 Swap the entries of x and y. More...
 
template<typename alpha_type , typename A_type , typename B_type , typename beta_type >
void SYMM (ESide side, EUplo uplo, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 Performs the matrix-matrix operation: C <- alpha*A*B+beta*C or C <- alpha*B*A+beta*C where A is an m by m or n by n symmetric matrix and B is a general matrix. More...
 
template<typename alpha_type , typename A_type , typename beta_type >
void SYRK (EUplo uplo, ETransp trans, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
 Performs the symmetric rank k operation: C <- alpha*A*A'+beta*C or C <- alpha*A'*A+beta*C, where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case or k by n matrix in the second case. More...
 
template<typename alpha_type , typename A_type >
void TRMM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
 Performs the matrix-matrix operation: B <- alpha*op(A)*B or B <- alpha*B*op(A) where op(A) is an unit/non-unit, upper/lower triangular matrix and B is a general matrix. More...
 
template<typename alpha_type , typename A_type >
void TRSM (ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
 Solves the matrix equations: op(A)*X=alpha*B or X*op(A)=alpha*B where X and B are m by n matrices, A is a unit/non-unit, upper/lower triangular matrix and op(A) is A or A'. The matrix X is overwritten on B. More...
 
- Protected Member Functions inherited from Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >
void copyMat (ScalarType *inputMatrix, OrdinalType strideInput, OrdinalType numRows, OrdinalType numCols, ScalarType *outputMatrix, OrdinalType strideOutput, OrdinalType startRow, OrdinalType startCol, ScalarType alpha=ScalarTraits< ScalarType >::zero())
 
void deleteArrays ()
 
void checkIndex (OrdinalType rowIndex, OrdinalType colIndex=0) const
 
- Static Protected Member Functions inherited from Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >
static ScalarType * allocateValues (const OrdinalType numRows, const OrdinalType numCols)
 
- Protected Attributes inherited from Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >
OrdinalType numRows_ = 0
 
OrdinalType numCols_ = 0
 
OrdinalType stride_ = 0
 
bool valuesCopied_ = false
 
ScalarType * values_ = nullptr
 
- Protected Attributes inherited from Teuchos::CompObject
FlopsflopCounter_
 

Detailed Description

template<typename OrdinalType, typename ScalarType>
class Teuchos::SerialDenseVector< OrdinalType, ScalarType >

This class creates and provides basic support for dense vectors of templated type as a specialization of Teuchos::SerialDenseMatrix. Additional methods for the SerialDenseVector class, like mathematical methods, can be found documented in SerialDenseMatrix.

Definition at line 60 of file Teuchos_SerialDenseVector.hpp.

Constructor & Destructor Documentation

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialDenseVector< OrdinalType, ScalarType >::SerialDenseVector ( )

Default Constructor.

Creates an empty vector of no length. The Sizing methods should be used to size this matrix. Values of this matrix should be set using the [] or the () operators.

Definition at line 218 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialDenseVector< OrdinalType, ScalarType >::SerialDenseVector ( OrdinalType  length,
bool  zeroOut = true 
)

Shaped Constructor.

Parameters
length- Number of elements in this vector.
zeroOut- Initializes values to 0 if true (default)

Creates a shaped vector of length length. All values are initialized to 0 when zeroOut is true. Values of this matrix should be set using the [] or the () operators.

Definition at line 221 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialDenseVector< OrdinalType, ScalarType >::SerialDenseVector ( DataAccess  CV,
ScalarType *  values,
OrdinalType  length 
)

Shaped Constructor with Values.

Parameters
CV- Enumerated type set to Teuchos::Copy or Teuchos::View.
values- Pointer to an array of ScalarType of the given length.
length- Length of vector to be constructed.

Definition at line 224 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialDenseVector< OrdinalType, ScalarType >::SerialDenseVector ( const SerialDenseVector< OrdinalType, ScalarType > &  Source)

Copy Constructor.

Definition at line 228 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialDenseVector< OrdinalType, ScalarType >::SerialDenseVector ( DataAccess  CV,
const SerialDenseVector< OrdinalType, ScalarType > &  Source 
)

Copy Constructor.

Note
Allow explicit control of whether a deep copy or view of Source is made with this copy constructor.

Definition at line 232 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
Teuchos::SerialDenseVector< OrdinalType, ScalarType >::~SerialDenseVector ( )
virtual

Destructor.

Definition at line 236 of file Teuchos_SerialDenseVector.hpp.

Member Function Documentation

template<typename OrdinalType, typename ScalarType>
int Teuchos::SerialDenseVector< OrdinalType, ScalarType >::size ( OrdinalType  length_in)
inline

Size method for changing the size of a SerialDenseVector, initializing entries to zero.

Parameters
length- The length of the new vector.

This allows the user to define the length of a SerialDenseVector at any point. This method can be called at any point after construction. Any values previously in this object will be destroyed and the resized vector starts with all zero values.

Definition at line 112 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType, typename ScalarType>
int Teuchos::SerialDenseVector< OrdinalType, ScalarType >::sizeUninitialized ( OrdinalType  length_in)
inline

Same as size() except leaves values uninitialized.

Definition at line 116 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType, typename ScalarType>
int Teuchos::SerialDenseVector< OrdinalType, ScalarType >::resize ( OrdinalType  length_in)
inline

Resizing method for changing the size of a SerialDenseVector, keeping the entries.

Parameters
length- The length of the new vector. This allows the user to redefine the length of a SerialDenseVector at any point. This method can be called at any point after construction. Any values previously in this object will be copied to the resized vector.

Definition at line 126 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType, typename ScalarType>
SerialDenseVector<OrdinalType, ScalarType>& Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator= ( const ScalarType  value)
inline

Set all values in the matrix to a constant value.

Parameters
value- Value to use;

Definition at line 137 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
bool Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator== ( const SerialDenseVector< OrdinalType, ScalarType > &  Operand) const

Equality of two matrices.

Returns
True if this vector and Operand are of the same length and have the same entries, else False will be returned.

Definition at line 246 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
bool Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator!= ( const SerialDenseVector< OrdinalType, ScalarType > &  Operand) const

Inequality of two matrices.

Returns
True if this vector and Operand are not of the same length or do not have the same entries, else False will be returned.

Definition at line 267 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
SerialDenseVector< OrdinalType, ScalarType > & Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator= ( const SerialDenseVector< OrdinalType, ScalarType > &  Source)

Copies values from one vector to another.

The operator= copies the values from one existing SerialDenseVector to another. If Source is a view (i.e. CV = Teuchos::View), then this method will return a view. Otherwise, it will return a copy of Source. this will be resized if it is not large enough to copy Source into.

Definition at line 239 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarType & Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator() ( OrdinalType  index)
inline

Element access method (non-const).

  Returns the ith element if x(i) is specified, the expression x[i] will return the same element.
Returns
(*this)(index)
Warning
The validity of index will only be checked if Teuchos is configured with –enable-teuchos-abc.

Definition at line 308 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
const ScalarType & Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator() ( OrdinalType  index) const
inline

Element access method (const).

  Returns the ith element if x(i) is specified, the expression x[i] will return the same element.
Returns
(*this)(index)
Warning
The validity of index will only be checked if Teuchos is configured with –enable-teuchos-abc.

Definition at line 317 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarType & Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator[] ( OrdinalType  index)
inline

Element access method (non-const).

  Returns the ith element if x[i] is specified, the expression x(i) will return the same element.
Returns
(*this)[index]
Warning
The validity of index will only be checked if Teuchos is configured with –enable-teuchos-abc.

Definition at line 335 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
const ScalarType & Teuchos::SerialDenseVector< OrdinalType, ScalarType >::operator[] ( OrdinalType  index) const
inline

Element access method (const).

  Returns the ith element if x[i] is specified, the expression x(i) will return the same element.
\return (*this)[index]
Warning
The validity of index will only be checked if Teuchos is configured with –enable-teuchos-abc.

Definition at line 326 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
ScalarType Teuchos::SerialDenseVector< OrdinalType, ScalarType >::dot ( const SerialDenseVector< OrdinalType, ScalarType > &  x) const

Compute the dot product of this vector and x.

Definition at line 273 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType, typename ScalarType>
OrdinalType Teuchos::SerialDenseVector< OrdinalType, ScalarType >::length ( ) const
inline

Returns the length of this vector.

Definition at line 207 of file Teuchos_SerialDenseVector.hpp.

template<typename OrdinalType , typename ScalarType >
std::ostream & Teuchos::SerialDenseVector< OrdinalType, ScalarType >::print ( std::ostream &  os) const
virtual

Print method. Define the behavior of the std::ostream << operator inherited from the Object class.

Reimplemented from Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >.

Definition at line 284 of file Teuchos_SerialDenseVector.hpp.


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