Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
List of all members
Epetra_SerialDenseVector Class Reference

Epetra_SerialDenseVector: A class for constructing and using dense vectors. More...

#include <Epetra_SerialDenseVector.h>

Inheritance diagram for Epetra_SerialDenseVector:
Inheritance graph
[legend]
Collaboration diagram for Epetra_SerialDenseVector:
Collaboration graph
[legend]

Public Member Functions

Constructors/destructors
 Epetra_SerialDenseVector ()
 Default constructor; defines a zero size object. More...
 
 Epetra_SerialDenseVector (int Length)
 Sized constructor; defines a variable-sized object. More...
 
 Epetra_SerialDenseVector (Epetra_DataAccess CV, double *Values, int Length)
 Set object values from one-dimensional array. More...
 
 Epetra_SerialDenseVector (const Epetra_SerialDenseVector &Source)
 Epetra_SerialDenseVector copy constructor.
 
virtual ~Epetra_SerialDenseVector ()
 Epetra_SerialDenseVector destructor.
 
Post-construction modification routines
int Size (int Length_in)
 Set length of a Epetra_SerialDenseVector object; init values to zero. More...
 
int Resize (int Length_in)
 Resize a Epetra_SerialDenseVector object. More...
 
Element access methods
Epetra_SerialDenseVectoroperator= (const Epetra_SerialDenseVector &Source)
 Value copy from one vector to another. More...
 
double & operator() (int Index)
 Element access function. More...
 
const double & operator() (int Index) const
 Element access function. More...
 
double & operator[] (int Index)
 Element access function. More...
 
const double & operator[] (int Index) const
 Column access function. More...
 
Mathematical methods
int Random ()
 Set vector values to random numbers. More...
 
double Dot (const Epetra_SerialDenseVector &x) const
 Compute 1-norm of each vector in multi-vector. More...
 
double Norm1 () const
 Compute 1-norm of each vector in multi-vector. More...
 
double Norm2 () const
 Compute 2-norm of each vector in multi-vector. More...
 
double NormInf () const
 Compute Inf-norm of each vector in multi-vector. More...
 
Attribute access methods
int Length () const
 Returns length of vector.
 
double * Values () const
 Returns pointer to the values in vector.
 
Epetra_DataAccess CV () const
 Returns the data access mode of the this vector.
 
I/O methods
virtual void Print (std::ostream &os) const
 Print service methods; defines behavior of ostream << operator.
 
- Public Member Functions inherited from Epetra_SerialDenseMatrix
 Epetra_SerialDenseMatrix (bool set_object_label=true)
 Default constructor; defines a zero size object. More...
 
 Epetra_SerialDenseMatrix (int NumRows, int NumCols, bool set_object_label=true)
 Shaped constructor; defines a variable-sized object. More...
 
 Epetra_SerialDenseMatrix (Epetra_DataAccess CV, double *A_in, int LDA_in, int NumRows, int NumCols, bool set_object_label=true)
 Set object values from two-dimensional array. More...
 
 Epetra_SerialDenseMatrix (const Epetra_SerialDenseMatrix &Source)
 Epetra_SerialDenseMatrix copy constructor.
 
virtual ~Epetra_SerialDenseMatrix ()
 Epetra_SerialDenseMatrix destructor.
 
int Shape (int NumRows, int NumCols)
 Set dimensions of a Epetra_SerialDenseMatrix object; init values to zero. More...
 
int Reshape (int NumRows, int NumCols)
 Reshape a Epetra_SerialDenseMatrix object. More...
 
int Multiply (char TransA, char TransB, double ScalarAB, const Epetra_SerialDenseMatrix &A, const Epetra_SerialDenseMatrix &B, double ScalarThis)
 Matrix-Matrix multiplication, this = ScalarThis*this + ScalarAB*A*B. More...
 
int Multiply (bool transA, const Epetra_SerialDenseMatrix &x, Epetra_SerialDenseMatrix &y)
 Matrix-Vector multiplication, y = A*x, where 'this' == A.
 
int Multiply (char SideA, double ScalarAB, const Epetra_SerialSymDenseMatrix &A, const Epetra_SerialDenseMatrix &B, double ScalarThis)
 Matrix-Matrix multiplication with a symmetric matrix A. More...
 
int Scale (double ScalarA)
 Inplace scalar-matrix product A = a A. More...
 
virtual double NormOne () const
 Computes the 1-Norm of the this matrix. More...
 
Epetra_SerialDenseMatrixoperator= (const Epetra_SerialDenseMatrix &Source)
 Value copy from one matrix to another. More...
 
bool operator== (const Epetra_SerialDenseMatrix &rhs) const
 Comparison operator. More...
 
bool operator!= (const Epetra_SerialDenseMatrix &rhs) const
 Inequality operator. More...
 
Epetra_SerialDenseMatrixoperator+= (const Epetra_SerialDenseMatrix &Source)
 Add one matrix to another. More...
 
double & operator() (int RowIndex, int ColIndex)
 Element access function. More...
 
const double & operator() (int RowIndex, int ColIndex) const
 Element access function. More...
 
double * operator[] (int ColIndex)
 Column access function. More...
 
const double * operator[] (int ColIndex) const
 Column access function. More...
 
int Random ()
 Set matrix values to random numbers. More...
 
int M () const
 Returns row dimension of system.
 
int N () const
 Returns column dimension of system.
 
double * A () const
 Returns pointer to the this matrix.
 
double * A ()
 Returns pointer to the this matrix.
 
int LDA () const
 Returns the leading dimension of the this matrix.
 
Epetra_DataAccess CV () const
 Returns the data access mode of the this matrix.
 
virtual double OneNorm () const
 Computes the 1-Norm of the this matrix (identical to NormOne() method). More...
 
virtual double InfNorm () const
 Computes the Infinity-Norm of the this matrix (identical to NormInf() method).
 
virtual int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied. More...
 
virtual int Apply (const Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &Y)
 Returns the result of a Epetra_SerialDenseOperator applied to a Epetra_SerialDenseMatrix X in Y. More...
 
virtual int ApplyInverse (const Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &Y)
 Returns the result of a Epetra_SerialDenseOperator inverse applied to an Epetra_SerialDenseMatrix X in Y. More...
 
virtual const char * Label () const
 Returns a character string describing the operator.
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting.
 
virtual bool HasNormInf () const
 Returns true if the this object can provide an approximate Inf-norm, false otherwise.
 
virtual int RowDim () const
 Returns the row dimension of operator.
 
virtual int ColDim () const
 Returns the column dimension of operator.
 
- Public Member Functions inherited from Epetra_CompObject
Epetra_CompObjectoperator= (const Epetra_CompObject &src)
 
 Epetra_CompObject ()
 Basic Epetra_CompObject constuctor.
 
 Epetra_CompObject (const Epetra_CompObject &Source)
 Epetra_CompObject copy constructor.
 
virtual ~Epetra_CompObject ()
 Epetra_CompObject destructor.
 
void SetFlopCounter (const Epetra_Flops &FlopCounter_in)
 Set the internal Epetra_Flops() pointer.
 
void SetFlopCounter (const Epetra_CompObject &CompObject)
 Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject.
 
void UnsetFlopCounter ()
 Set the internal Epetra_Flops() pointer to 0 (no flops counted).
 
Epetra_FlopsGetFlopCounter () const
 Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none.
 
void ResetFlops () const
 Resets the number of floating point operations to zero for this multi-vector.
 
double Flops () const
 Returns the number of floating point operations with this multi-vector.
 
void UpdateFlops (int Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (long int Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (long long Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (double Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (float Flops_in) const
 Increment Flop count for this object.
 
- Public Member Functions inherited from Epetra_Object
 Epetra_Object (int TracebackModeIn=-1, bool set_label=true)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const char *const Label, int TracebackModeIn=-1)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const Epetra_Object &Object)
 Epetra_Object Copy Constructor. More...
 
virtual ~Epetra_Object ()
 Epetra_Object Destructor. More...
 
virtual int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method.
 
virtual void SetLabel (const char *const Label)
 Epetra_Object Label definition using char *. More...
 
- Public Member Functions inherited from Epetra_SerialDenseOperator
virtual ~Epetra_SerialDenseOperator ()
 Destructor.
 
- Public Member Functions inherited from Epetra_BLAS
 Epetra_BLAS (void)
 Epetra_BLAS Constructor. More...
 
 Epetra_BLAS (const Epetra_BLAS &BLAS)
 Epetra_BLAS Copy Constructor. More...
 
virtual ~Epetra_BLAS (void)
 Epetra_BLAS Destructor.
 
float ASUM (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS one norm function (SASUM).
 
double ASUM (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS one norm function (DASUM).
 
float DOT (const int N, const float *X, const float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS dot product function (SDOT).
 
double DOT (const int N, const double *X, const double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS dot product function (DDOT).
 
float NRM2 (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS norm function (SNRM2).
 
double NRM2 (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS norm function (DNRM2).
 
void SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const
 Epetra_BLAS vector scale function (SSCAL)
 
void SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const
 Epetra_BLAS vector scale function (DSCAL)
 
void COPY (const int N, const float *X, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector copy function (SCOPY)
 
void COPY (const int N, const double *X, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector scale function (DCOPY)
 
int IAMAX (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS arg maximum of absolute value function (ISAMAX)
 
int IAMAX (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS arg maximum of absolute value function (IDAMAX)
 
void AXPY (const int N, const float ALPHA, const float *X, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector update function (SAXPY)
 
void AXPY (const int N, const double ALPHA, const double *X, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector update function (DAXPY)
 
void GEMV (const char TRANS, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *X, const float BETA, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS matrix-vector multiply function (SGEMV)
 
void GEMV (const char TRANS, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *X, const double BETA, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS matrix-vector multiply function (DGEMV)
 
void GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (SGEMM)
 
void GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (DGEMM)
 
void SYMM (const char SIDE, const char UPLO, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
 Epetra_BLAS symmetric matrix-matrix multiply function (SSYMM)
 
void SYMM (const char SIDE, const char UPLO, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (DSYMM)
 
void TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const float ALPHA, const float *A, const int LDA, float *B, const int LDB) const
 Epetra_BLAS triangular matrix-matrix multiply function (STRMM)
 
void TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const double ALPHA, const double *A, const int LDA, double *B, const int LDB) const
 Epetra_BLAS triangular matrix-matrix multiply function (DTRMM)
 
void SYRK (const char UPLO, const char TRANS, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float BETA, float *C, const int LDC) const
 Eperta_BLAS symetric rank k funtion (ssyrk)
 
void SYRK (const char UPLO, const char TRANS, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double BETA, double *C, const int LDC) const
 Eperta_BLAS symetric rank k funtion (dsyrk)
 

Additional Inherited Members

- Static Public Member Functions inherited from Epetra_Object
static void SetTracebackMode (int TracebackModeValue)
 Set the value of the Epetra_Object error traceback report mode. More...
 
static int GetTracebackMode ()
 Get the value of the Epetra_Object error report mode.
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting.
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 
- Protected Member Functions inherited from Epetra_SerialDenseMatrix
void CopyMat (const double *Source, int Source_LDA, int NumRows, int NumCols, double *Target, int Target_LDA, bool add=false)
 
void CleanupData ()
 
- Protected Member Functions inherited from Epetra_Object
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 
- Protected Attributes inherited from Epetra_SerialDenseMatrix
int M_
 
int N_
 
bool A_Copied_
 
Epetra_DataAccess CV_
 
int LDA_
 
double * A_
 
bool UseTranspose_
 
- Protected Attributes inherited from Epetra_CompObject
Epetra_FlopsFlopCounter_
 

Detailed Description

Epetra_SerialDenseVector: A class for constructing and using dense vectors.

The Epetra_SerialDenseVector class enables the construction and use of real-valued,
double-precision dense vectors.  It is built on the BLAS and LAPACK and derives from the Epetra_SerialDenseMatrix class.

The Epetra_SerialDenseVector class is intended to provide convenient vector notation but derives all signficant functionality from Epetra_SerialDenseMatrix.

Constructing Epetra_SerialDenseVector Objects

There are four Epetra_SerialDenseVector constructors. The first constructs a zero-length object which should be made to appropriate length using the Size() or Resize() functions and then filled with the [] or () operators. The second constructs an object sized to the dimension specified, which should be filled with the [] or () operators. The third is a constructor that accepts user data as a 1D array, and the fourth is a copy constructor. The third constructor has two data access modes (specified by the Epetra_DataAccess argument):

  1. Copy mode - Allocates memory and makes a copy of the user-provided data. In this case, the user data is not needed after construction.
  2. View mode - Creates a "view" of the user data. In this case, the user data is required to remain intact for the life of the object.
Warning
View mode is extremely dangerous from a data hiding perspective. Therefore, we strongly encourage users to develop code using Copy mode first and only use the View mode in a secondary optimization phase.

Extracting Data from Epetra_SerialDenseVector Objects

Once a Epetra_SerialDenseVector is constructed, it is possible to view the data via access functions.

Warning
Use of these access functions cam be extremely dangerous from a data hiding perspective.

The final useful function is Flops(). Each Epetra_SerialDenseVector object keep track of the number of serial floating point operations performed using the specified object as the this argument to the function. The Flops() function returns this number as a double precision number. Using this information, in conjunction with the Epetra_Time class, one can get accurate parallel performance numbers.

Constructor & Destructor Documentation

Epetra_SerialDenseVector::Epetra_SerialDenseVector ( )

Default constructor; defines a zero size object.

Epetra_SerialDenseVector objects defined by the default constructor should be sized with the Size() or Resize functions. Values should be defined by using the [] or () operators.

Epetra_SerialDenseVector::Epetra_SerialDenseVector ( int  Length)

Sized constructor; defines a variable-sized object.

Parameters
InLength - Length of vector.

Epetra_SerialDenseVector objects defined by the sized constructor are already sized to the dimension given as a parameter. All values are initialized to 0. Calling this constructor is equivalent to using the default constructor, and then calling the Size function on it. Values should be defined by using the [] or () operators.

Epetra_SerialDenseVector::Epetra_SerialDenseVector ( Epetra_DataAccess  CV,
double *  Values,
int  Length 
)

Set object values from one-dimensional array.

Parameters
InEpetra_DataAccess - Enumerated type set to Copy or View.
InValues - Pointer to an array of double precision numbers containing the values.
InLength - Length of vector.

See Detailed Description section for further discussion.

Member Function Documentation

double Epetra_SerialDenseVector::Dot ( const Epetra_SerialDenseVector x) const

Compute 1-norm of each vector in multi-vector.

Parameters
x(In) Input vector x.
Returns
Dot-product of the this vector and x.
double Epetra_SerialDenseVector::Norm1 ( ) const

Compute 1-norm of each vector in multi-vector.

Returns
1-norm of the vector.
double Epetra_SerialDenseVector::Norm2 ( ) const

Compute 2-norm of each vector in multi-vector.

Parameters
Out
Returns
2-norm of the vector.
double Epetra_SerialDenseVector::NormInf ( ) const
virtual

Compute Inf-norm of each vector in multi-vector.

Returns
Infinity-norm of the vector.

Reimplemented from Epetra_SerialDenseMatrix.

double & Epetra_SerialDenseVector::operator() ( int  Index)
inline

Element access function.

Returns the specified element of the vector. Bounds checking is enforced.

Returns
Specified element in vector.
Warning
No bounds checking is done unless Epetra is compiled with HAVE_EPETRA_ARRAY_BOUNDS_CHECK.

References Epetra_Object::ReportError().

const double & Epetra_SerialDenseVector::operator() ( int  Index) const
inline

Element access function.

Returns the specified element of the vector. Bounds checking is enforced.

Returns
Specified element in vector.
Warning
No bounds checking is done unless Epetra is compiled with HAVE_EPETRA_ARRAY_BOUNDS_CHECK.

References Epetra_Object::ReportError().

Epetra_SerialDenseVector& Epetra_SerialDenseVector::operator= ( const Epetra_SerialDenseVector Source)

Value copy from one vector to another.

The operator= allows one to copy the values from one existing SerialDenseVector to another, as long as there is enough room in the target to hold the source.

Returns
Values of the left hand side vector are modified by the values of the right hand side vector.
double & Epetra_SerialDenseVector::operator[] ( int  Index)
inline

Element access function.

Returns the specified element of the vector.

Returns
Specified element in vector.
Warning
No bounds checking is done unless Epetra is compiled with HAVE_EPETRA_ARRAY_BOUNDS_CHECK.

References Epetra_Object::ReportError().

const double & Epetra_SerialDenseVector::operator[] ( int  Index) const
inline

Column access function.

Returns the specified element of the vector.

Returns
Specified element in vector.
Warning
No bounds checking is done unless Epetra is compiled with HAVE_EPETRA_ARRAY_BOUNDS_CHECK.

References Epetra_Object::ReportError().

int Epetra_SerialDenseVector::Random ( )

Set vector values to random numbers.

SerialDenseVector uses the random number generator provided by Epetra_Util.
The vector values will be set to random values on the interval (-1.0, 1.0).
Returns
Integer error code, set to 0 if successful.
int Epetra_SerialDenseVector::Resize ( int  Length_in)
inline

Resize a Epetra_SerialDenseVector object.

Parameters
InLength - Length of vector object.

Allows user to define the dimension of a Epetra_SerialDenseVector. This function can be called at any point after construction. Any values that were previously in this object are copied into the new size. If the new shape is smaller than the original, the first Length values are copied to the new vector.

Returns
Integer error code, set to 0 if successful.

References Epetra_SerialDenseMatrix::Reshape().

int Epetra_SerialDenseVector::Size ( int  Length_in)
inline

Set length of a Epetra_SerialDenseVector object; init values to zero.

Parameters
InLength - Length of vector object.

Allows user to define the dimension of a Epetra_SerialDenseVector. This function can be called at any point after construction. Any values that were previously in this object are destroyed and the resized vector starts off with all zero values.

Returns
Integer error code, set to 0 if successful.

References Epetra_SerialDenseMatrix::Shape().


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