Epetra Package Browser (Single Doxygen Collection)
Development
|
Epetra_SerialSymDenseMatrix: A class for constructing and using symmetric positive definite dense matrices. More...
#include <Epetra_SerialSymDenseMatrix.h>
Public Member Functions | |
void | CopyUPLOMat (bool Upper, double *A, int LDA, int NumRows) |
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. More... | |
virtual | ~Epetra_SerialDenseMatrix () |
Epetra_SerialDenseMatrix destructor. More... | |
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. More... | |
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... | |
Epetra_SerialDenseMatrix & | operator= (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_SerialDenseMatrix & | operator+= (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. More... | |
int | N () const |
Returns column dimension of system. More... | |
double * | A () const |
Returns pointer to the this matrix. More... | |
double * | A () |
Returns pointer to the this matrix. More... | |
int | LDA () const |
Returns the leading dimension of the this matrix. More... | |
Epetra_DataAccess | CV () const |
Returns the data access mode of the this matrix. More... | |
virtual void | Print (std::ostream &os) const |
Print service methods; defines behavior of ostream << operator. More... | |
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. More... | |
virtual bool | UseTranspose () const |
Returns the current UseTranspose setting. More... | |
virtual bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. More... | |
virtual int | RowDim () const |
Returns the row dimension of operator. More... | |
virtual int | ColDim () const |
Returns the column dimension of operator. More... | |
Public Member Functions inherited from Epetra_CompObject | |
Epetra_CompObject & | operator= (const Epetra_CompObject &src) |
Epetra_CompObject () | |
Basic Epetra_CompObject constuctor. More... | |
Epetra_CompObject (const Epetra_CompObject &Source) | |
Epetra_CompObject copy constructor. More... | |
virtual | ~Epetra_CompObject () |
Epetra_CompObject destructor. More... | |
void | SetFlopCounter (const Epetra_Flops &FlopCounter_in) |
Set the internal Epetra_Flops() pointer. More... | |
void | SetFlopCounter (const Epetra_CompObject &CompObject) |
Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject. More... | |
void | UnsetFlopCounter () |
Set the internal Epetra_Flops() pointer to 0 (no flops counted). More... | |
Epetra_Flops * | GetFlopCounter () const |
Get the pointer to the Epetra_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-vector. More... | |
double | Flops () const |
Returns the number of floating point operations with this multi-vector. More... | |
void | UpdateFlops (int Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (long int Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (long long Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (double Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (float Flops_in) const |
Increment Flop count for this object. More... | |
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 void | SetLabel (const char *const Label) |
Epetra_Object Label definition using char *. More... | |
virtual int | ReportError (const std::string Message, int ErrorCode) const |
Error reporting method. More... | |
Public Member Functions inherited from Epetra_SerialDenseOperator | |
virtual | ~Epetra_SerialDenseOperator () |
Destructor. More... | |
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. More... | |
float | ASUM (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS one norm function (SASUM). More... | |
double | ASUM (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS one norm function (DASUM). More... | |
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). More... | |
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). More... | |
float | NRM2 (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS norm function (SNRM2). More... | |
double | NRM2 (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS norm function (DNRM2). More... | |
void | SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const |
Epetra_BLAS vector scale function (SSCAL) More... | |
void | SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const |
Epetra_BLAS vector scale function (DSCAL) More... | |
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) More... | |
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) More... | |
int | IAMAX (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (ISAMAX) More... | |
int | IAMAX (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (IDAMAX) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
Private Attributes | |
bool | Upper_ |
char | UPLO_ |
Constructor/Destructor Methods | |
Epetra_SerialSymDenseMatrix (void) | |
Default constructor; defines a zero size object. More... | |
Epetra_SerialSymDenseMatrix (Epetra_DataAccess CV, double *A, int LDA, int NumRowsCols) | |
Set object values from two-dimensional array. More... | |
Epetra_SerialSymDenseMatrix (const Epetra_SerialSymDenseMatrix &Source) | |
Epetra_SerialSymDenseMatrix copy constructor. More... | |
virtual | ~Epetra_SerialSymDenseMatrix () |
Epetra_SerialSymDenseMatrix destructor. More... | |
Set Methods | |
int | Shape (int NumRowsCols) |
Set dimensions of a Epetra_SerialSymDenseMatrix object; init values to zero. More... | |
int | Reshape (int NumRowsCols) |
Reshape a Epetra_SerialSymDenseMatrix object. More... | |
void | SetLower () |
Specify that the lower triangle of the this matrix should be used. More... | |
void | SetUpper () |
Specify that the upper triangle of the this matrix should be used. More... | |
Query methods | |
bool | Upper () const |
Returns true if upper triangle of this matrix has and will be used. More... | |
char | UPLO () const |
Returns character value of UPLO used by LAPACK routines. More... | |
Mathematical Methods | |
int | Scale (double ScalarA) |
Inplace scalar-matrix product A = a A. More... | |
double | NormOne () const |
Computes the 1-Norm of the this matrix. More... | |
double | NormInf () const |
Computes the Infinity-Norm of the this matrix. More... | |
Deprecated methods (will be removed in later versions of this class) | |
double | OneNorm () const |
Computes the 1-Norm of the this matrix (identical to NormOne() method). More... | |
double | InfNorm () const |
Computes the Infinity-Norm of the this matrix (identical to NormInf() method). More... | |
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. More... | |
static std::ostream & | GetTracebackStream () |
Get the output stream for error reporting. More... | |
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_Flops * | FlopCounter_ |
Epetra_SerialSymDenseMatrix: A class for constructing and using symmetric positive definite dense matrices.
The Epetra_SerialSymDenseMatrix class enables the construction and use of real-valued, symmetric positive definite, double-precision dense matrices. It is built on the Epetra_SerialDenseMatrix class which in turn is built on the BLAS via the Epetra_BLAS class.
The Epetra_SerialSymDenseMatrix class is intended to provide full-featured support for solving linear and eigen system problems for symmetric positive definite matrices. It is written on top of BLAS and LAPACK and thus has excellent performance and numerical capabilities. Using this class, one can either perform simple factorizations and solves or apply all the tricks available in LAPACK to get the best possible solution for very ill-conditioned problems.
Epetra_SerialSymDenseMatrix vs. Epetra_LAPACK
The Epetra_LAPACK class provides access to most of the same functionality as Epetra_SerialSymDenseMatrix. The primary difference is that Epetra_LAPACK is a "thin" layer on top of LAPACK and Epetra_SerialSymDenseMatrix attempts to provide easy access to the more sophisticated aspects of solving dense linear and eigensystems.
When you should use Epetra_SerialSymDenseMatrix: If you want to (or potentially want to) solve ill-conditioned problems or want to work with a more object-oriented interface, you should probably use Epetra_SerialSymDenseMatrix.
Constructing Epetra_SerialSymDenseMatrix Objects
There are three Epetra_DenseMatrix constructors. The first constructs a zero-sized object which should be made to appropriate length using the Shape() or Reshape() functions and then filled with the [] or () operators. The second is a constructor that accepts user data as a 2D array, the third is a copy constructor. The second constructor has two data access modes (specified by the Epetra_DataAccess argument):
Extracting Data from Epetra_SerialSymDenseMatrix Objects
Once a Epetra_SerialSymDenseMatrix is constructed, it is possible to view the data via access functions.
Vector and Utility Functions
Once a Epetra_SerialSymDenseMatrix is constructed, several mathematical functions can be applied to the object. Specifically:
Counting floating point operations The Epetra_SerialSymDenseMatrix class has Epetra_CompObject as a base class. Thus, floating point operations are counted and accumulated in the Epetra_Flop object (if any) that was set using the SetFlopCounter() method in the Epetra_CompObject base class.
Definition at line 130 of file Epetra_SerialSymDenseMatrix.h.
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix | ( | void | ) |
Default constructor; defines a zero size object.
Epetra_SerialSymDenseMatrix objects defined by the default constructor should be sized with the Shape() or Reshape() functions. Values should be defined by using the [] or ()operators.
Note: By default the active part of the matrix is assumed to be in the lower triangle. To set the upper part as active, call SetUpper(). See Detailed Description section for further discussion.
Definition at line 45 of file Epetra_SerialSymDenseMatrix.cpp.
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix | ( | Epetra_DataAccess | CV, |
double * | A, | ||
int | LDA, | ||
int | NumRowsCols | ||
) |
Set object values from two-dimensional array.
In | Epetra_DataAccess - Enumerated type set to Copy or View. |
In | A - Pointer to an array of double precision numbers. The first vector starts at A. The second vector starts at A+LDA, the third at A+2*LDA, and so on. |
In | LDA - The "Leading Dimension", or stride between vectors in memory. |
In | NumRowsCols - Number of rows and columns in object. |
Note: By default the active part of the matrix is assumed to be in the lower triangle. To set the upper part as active, call SetUpper(). See Detailed Description section for further discussion.
Definition at line 53 of file Epetra_SerialSymDenseMatrix.cpp.
Epetra_SerialSymDenseMatrix::Epetra_SerialSymDenseMatrix | ( | const Epetra_SerialSymDenseMatrix & | Source | ) |
Epetra_SerialSymDenseMatrix copy constructor.
Definition at line 61 of file Epetra_SerialSymDenseMatrix.cpp.
|
virtual |
Epetra_SerialSymDenseMatrix destructor.
Definition at line 68 of file Epetra_SerialSymDenseMatrix.cpp.
|
inline |
Set dimensions of a Epetra_SerialSymDenseMatrix object; init values to zero.
In | NumRowsCols - Number of rows and columns in object. |
Allows user to define the dimensions of a Epetra_DenseMatrix at any point. This function can be called at any point after construction. Any values that were previously in this object are destroyed and the resized matrix starts off with all zero values.
Definition at line 192 of file Epetra_SerialSymDenseMatrix.h.
|
inline |
Reshape a Epetra_SerialSymDenseMatrix object.
In | NumRowsCols - Number of rows and columns in object. |
Allows user to define the dimensions of a Epetra_SerialSymDenseMatrix at any point. This function can be called at any point after construction. Any values that were previously in this object are copied into the new shape. If the new shape is smaller than the original, the upper left portion of the original matrix (the principal submatrix) is copied to the new matrix.
Definition at line 211 of file Epetra_SerialSymDenseMatrix.h.
|
inline |
Specify that the lower triangle of the this matrix should be used.
Definition at line 215 of file Epetra_SerialSymDenseMatrix.h.
|
inline |
Specify that the upper triangle of the this matrix should be used.
Definition at line 218 of file Epetra_SerialSymDenseMatrix.h.
|
inline |
Returns true if upper triangle of this matrix has and will be used.
Definition at line 225 of file Epetra_SerialSymDenseMatrix.h.
|
inline |
Returns character value of UPLO used by LAPACK routines.
Definition at line 228 of file Epetra_SerialSymDenseMatrix.h.
int Epetra_SerialSymDenseMatrix::Scale | ( | double | ScalarA | ) |
Inplace scalar-matrix product A = a A.
Scale a matrix, entry-by-entry using the value ScalarA. This method is sensitive to the UPLO() parameter.
ScalarA | (In) Scalar to multiply with A. |
Definition at line 143 of file Epetra_SerialSymDenseMatrix.cpp.
|
virtual |
Computes the 1-Norm of the this matrix.
Reimplemented from Epetra_SerialDenseMatrix.
Definition at line 100 of file Epetra_SerialSymDenseMatrix.cpp.
|
virtual |
Computes the Infinity-Norm of the this matrix.
Reimplemented from Epetra_SerialDenseMatrix.
Definition at line 105 of file Epetra_SerialSymDenseMatrix.cpp.
void Epetra_SerialSymDenseMatrix::CopyUPLOMat | ( | bool | Upper, |
double * | A, | ||
int | LDA, | ||
int | NumRows | ||
) |
Definition at line 72 of file Epetra_SerialSymDenseMatrix.cpp.
|
inlinevirtual |
Computes the 1-Norm of the this matrix (identical to NormOne() method).
Reimplemented from Epetra_SerialDenseMatrix.
Definition at line 267 of file Epetra_SerialSymDenseMatrix.h.
|
inlinevirtual |
Computes the Infinity-Norm of the this matrix (identical to NormInf() method).
Reimplemented from Epetra_SerialDenseMatrix.
Definition at line 270 of file Epetra_SerialSymDenseMatrix.h.
|
private |
Definition at line 270 of file Epetra_SerialSymDenseMatrix.h.
|
private |
Definition at line 277 of file Epetra_SerialSymDenseMatrix.h.