Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
Public Member Functions | Public Attributes | List of all members
Epetra_SerialDenseSVD Class Reference

Epetra_SerialDenseSVD: A class for SVDing dense linear problems. More...

#include <Epetra_SerialDenseSVD.h>

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

Public Member Functions

void AllocateWORK ()
 
void AllocateIWORK ()
 
void InitPointers ()
 
void DeleteArrays ()
 
void ResetMatrix ()
 
void ResetVectors ()
 
Constructor/Destructor Methods
 Epetra_SerialDenseSVD ()
 Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors().
 
virtual ~Epetra_SerialDenseSVD ()
 Epetra_SerialDenseSVD destructor.
 
Set Methods
int SetMatrix (Epetra_SerialDenseMatrix &A)
 Sets the pointers for coefficient matrix.
 
int SetVectors (Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
 Sets the pointers for left and right hand side vector(s). More...
 
Strategy modifying Methods
void SolveWithTranspose (bool Flag)
 Causes equilibration to be called just before the matrix factorization as part of the call to Factor. More...
 
Factor/Solve/Invert Methods

Causes all solves to compute solution to best ability using iterative refinement.

virtual int Factor (void)
 
virtual int Solve (void)
 Computes the solution X to AX = B for the this matrix and the B provided to SetVectors().. More...
 
virtual int Invert (double rthresh=0.0, double athresh=0.0)
 Inverts the this matrix. More...
 
Query methods
bool Transpose ()
 Returns true if transpose of this matrix has and will be used.
 
bool Factored ()
 Returns true if matrix is factored (factor available via AF() and LDAF()).
 
bool Inverted ()
 Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
 
bool Solved ()
 Returns true if the current set of vectors has been solved.
 
Data Accessor methods
Epetra_SerialDenseMatrixMatrix () const
 Returns pointer to current matrix.
 
Epetra_SerialDenseMatrixInvertedMatrix () const
 Returns pointer to inverted matrix (assuming inverse has been performed).
 
Epetra_SerialDenseMatrixLHS () const
 Returns pointer to current LHS.
 
Epetra_SerialDenseMatrixRHS () const
 Returns pointer to current RHS.
 
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.
 
int LDA () const
 Returns the leading dimension of the this matrix.
 
double * B () const
 Returns pointer to current RHS.
 
int LDB () const
 Returns the leading dimension of the RHS.
 
int NRHS () const
 Returns the number of current right hand sides and solution vectors.
 
double * X () const
 Returns pointer to current solution.
 
int LDX () const
 Returns the leading dimension of the solution.
 
double * S () const
 
double * AI () const
 Returns pointer to the inverted matrix (may be the same as A() if factorization done in place).
 
int LDAI () const
 Returns the leading dimension of the inverted matrix.
 
double ANORM () const
 Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
 
I/O methods
virtual void Print (std::ostream &os) const
 Print service methods; defines behavior of ostream << operator.
 
Additional methods for support of Epetra_SerialDenseOperator interface
virtual int SetUseTranspose (bool use_transpose)
 If set true, transpose of this operator will be applied. More...
 
virtual int Apply (const Epetra_SerialDenseMatrix &Xmat, Epetra_SerialDenseMatrix &Ymat)
 Returns the result of a Epetra_SerialDenseOperator applied to a Epetra_SerialDenseMatrix X in Y. More...
 
virtual int ApplyInverse (const Epetra_SerialDenseMatrix &Xmat, Epetra_SerialDenseMatrix &Ymat)
 Returns the result of a Epetra_SerialDenseOperator inverse applied to an Epetra_SerialDenseMatrix X in Y. More...
 
virtual double NormInf () const
 Returns the infinity norm of the global matrix.
 
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_SerialDenseOperator
virtual ~Epetra_SerialDenseOperator ()
 Destructor.
 
- 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_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)
 
- Public Member Functions inherited from Epetra_LAPACK
 Epetra_LAPACK (void)
 Epetra_LAPACK Constructor. More...
 
 Epetra_LAPACK (const Epetra_LAPACK &LAPACK)
 Epetra_LAPACK Copy Constructor. More...
 
virtual ~Epetra_LAPACK (void)
 Epetra_LAPACK Destructor.
 
void POTRF (const char UPLO, const int N, float *A, const int LDA, int *INFO) const
 Epetra_LAPACK factorization for positive definite matrix (SPOTRF)
 
void POTRF (const char UPLO, const int N, double *A, const int LDA, int *INFO) const
 Epetra_LAPACK factorization for positive definite matrix (DPOTRF)
 
void POTRS (const char UPLO, const int N, const int NRHS, const float *A, const int LDA, float *X, const int LDX, int *INFO) const
 Epetra_LAPACK solve (after factorization) for positive definite matrix (SPOTRS)
 
void POTRS (const char UPLO, const int N, const int NRHS, const double *A, const int LDA, double *X, const int LDX, int *INFO) const
 Epetra_LAPACK solve (after factorization) for positive definite matrix (DPOTRS)
 
void POTRI (const char UPLO, const int N, float *A, const int LDA, int *INFO) const
 Epetra_LAPACK inversion for positive definite matrix (SPOTRI)
 
void POTRI (const char UPLO, const int N, double *A, const int LDA, int *INFO) const
 Epetra_LAPACK inversion for positive definite matrix (DPOTRI)
 
void POCON (const char UPLO, const int N, const float *A, const int LDA, const float ANORM, float *RCOND, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK condition number estimator for positive definite matrix (SPOCON)
 
void POCON (const char UPLO, const int N, const double *A, const int LDA, const double ANORM, double *RCOND, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK condition number estimator for positive definite matrix (DPOCON)
 
void POSV (const char UPLO, const int N, const int NRHS, float *A, const int LDA, float *X, const int LDX, int *INFO) const
 Epetra_LAPACK factor and solve for positive definite matrix (SPOSV)
 
void POSV (const char UPLO, const int N, const int NRHS, double *A, const int LDA, double *X, const int LDX, int *INFO) const
 Epetra_LAPACK factor and solve for positive definite matrix (DPOSV)
 
void POEQU (const int N, const float *A, const int LDA, float *S, float *SCOND, float *AMAX, int *INFO) const
 Epetra_LAPACK equilibration for positive definite matrix (SPOEQU)
 
void POEQU (const int N, const double *A, const int LDA, double *S, double *SCOND, double *AMAX, int *INFO) const
 Epetra_LAPACK equilibration for positive definite matrix (DPOEQU)
 
void PORFS (const char UPLO, const int N, const int NRHS, const float *A, const int LDA, const float *AF, const int LDAF, const float *B, const int LDB, float *X, const int LDX, float *FERR, float *BERR, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK solve driver for positive definite matrix (SPOSVX)
 
void PORFS (const char UPLO, const int N, const int NRHS, const double *A, const int LDA, const double *AF, const int LDAF, const double *B, const int LDB, double *X, const int LDX, double *FERR, double *BERR, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK solve driver for positive definite matrix (DPOSVX)
 
void POSVX (const char FACT, const char UPLO, const int N, const int NRHS, float *A, const int LDA, float *AF, const int LDAF, const char EQUED, float *S, float *B, const int LDB, float *X, const int LDX, float *RCOND, float *FERR, float *BERR, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK solve driver for positive definite matrix (SPOSVX)
 
void POSVX (const char FACT, const char UPLO, const int N, const int NRHS, double *A, const int LDA, double *AF, const int LDAF, const char EQUED, double *S, double *B, const int LDB, double *X, const int LDX, double *RCOND, double *FERR, double *BERR, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK solve driver for positive definite matrix (DPOSVX)
 
void GELS (const char TRANS, const int M, const int N, const int NRHS, double *A, const int LDA, double *B, const int LDB, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK simple driver to solve least-squares systems.
 
void GETRF (const int M, const int N, float *A, const int LDA, int *IPIV, int *INFO) const
 Epetra_LAPACK factorization for general matrix (SGETRF)
 
void GETRF (const int M, const int N, double *A, const int LDA, int *IPIV, int *INFO) const
 Epetra_LAPACK factorization for general matrix (DGETRF)
 
void GEQRF (const int M, const int N, float *A, const int LDA, float *TAU, float *WORK, const int lwork, int *INFO) const
 Epetra_LAPACK QR factorization for general matrix (SGEQRF)
 
void GEQRF (const int M, const int N, double *A, const int LDA, double *TAU, double *WORK, const int lwork, int *INFO) const
 Epetra_LAPACK factorization for general matrix (DGEQRF)
 
void GETRS (const char TRANS, const int N, const int NRHS, const float *A, const int LDA, const int *IPIV, float *X, const int LDX, int *INFO) const
 Epetra_LAPACK solve (after factorization) for general matrix (SGETRS)
 
void GETRS (const char TRANS, const int N, const int NRHS, const double *A, const int LDA, const int *IPIV, double *X, const int LDX, int *INFO) const
 Epetra_LAPACK solve (after factorization) for general matrix (DGETRS)
 
void GETRI (const int N, float *A, const int LDA, int *IPIV, float *WORK, const int *LWORK, int *INFO) const
 Epetra_LAPACK inversion for general matrix (SGETRI)
 
void GETRI (const int N, double *A, const int LDA, int *IPIV, double *WORK, const int *LWORK, int *INFO) const
 Epetra_LAPACK inversion for general matrix (DGETRI)
 
void GECON (const char NORM, const int N, const float *A, const int LDA, const float ANORM, float *RCOND, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK condition number estimator for general matrix (SGECON)
 
void GECON (const char NORM, const int N, const double *A, const int LDA, const double ANORM, double *RCOND, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK condition number estimator for general matrix (DGECON)
 
void GESV (const int N, const int NRHS, float *A, const int LDA, int *IPIV, float *X, const int LDX, int *INFO) const
 Epetra_LAPACK factor and solve for general matrix (SGESV)
 
void GESV (const int N, const int NRHS, double *A, const int LDA, int *IPIV, double *X, const int LDX, int *INFO) const
 Epetra_LAPACK factor and solve for general matrix (DGESV)
 
void GEEQU (const int M, const int N, const float *A, const int LDA, float *R, float *C, float *ROWCND, float *COLCND, float *AMAX, int *INFO) const
 Epetra_LAPACK equilibration for general matrix (SGEEQU)
 
void GEEQU (const int M, const int N, const double *A, const int LDA, double *R, double *C, double *ROWCND, double *COLCND, double *AMAX, int *INFO) const
 Epetra_LAPACK equilibration for general matrix (DGEEQU)
 
void GERFS (const char TRANS, const int N, const int NRHS, const float *A, const int LDA, const float *AF, const int LDAF, const int *IPIV, const float *B, const int LDB, float *X, const int LDX, float *FERR, float *BERR, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK Refine solution (GERFS)
 
void GERFS (const char TRANS, const int N, const int NRHS, const double *A, const int LDA, const double *AF, const int LDAF, const int *IPIV, const double *B, const int LDB, double *X, const int LDX, double *FERR, double *BERR, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK Refine solution (GERFS)
 
void GESVX (const char FACT, const char TRANS, const int N, const int NRHS, float *A, const int LDA, float *AF, const int LDAF, int *IPIV, const char EQUED, float *R, float *C, float *B, const int LDB, float *X, const int LDX, float *RCOND, float *FERR, float *BERR, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK solve driver for general matrix (SGESVX)
 
void GESVX (const char FACT, const char TRANS, const int N, const int NRHS, double *A, const int LDA, double *AF, const int LDAF, int *IPIV, const char EQUED, double *R, double *C, double *B, const int LDB, double *X, const int LDX, double *RCOND, double *FERR, double *BERR, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK solve driver for general matrix (DGESVX)
 
void GEHRD (const int N, const int ILO, const int IHI, float *A, const int LDA, float *TAU, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for reduction to Hessenberg form (SGEHRD)
 
void GEHRD (const int N, const int ILO, const int IHI, double *A, const int LDA, double *TAU, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for reduction to Hessenberg form (DGEHRD)
 
void HSEQR (const char JOB, const char COMPZ, const int N, const int ILO, const int IHI, float *H, const int LDH, float *WR, float *WI, float *Z, const int LDZ, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for computing the eigenvalues of a real upper Hessenberg matrix (SHSEQR)
 
void HSEQR (const char JOB, const char COMPZ, const int N, const int ILO, const int IHI, double *H, const int LDH, double *WR, double *WI, double *Z, const int LDZ, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for computing the eigenvalues of a real upper Hessenberg matrix (DHSEQR)
 
void ORGQR (const int M, const int N, const int K, float *A, const int LDA, float *TAU, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for generating a m x n real matrix Q with orthonormal columns, defined as the product of k elementary reflectors. (SORGQR)
 
void ORGQR (const int M, const int N, const int K, double *A, const int LDA, double *TAU, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for generating a m x n real matrix Q with orthonormal columns, defined as the product of k elementary reflectors. (DORGQR)
 
void ORGHR (const int N, const int ILO, const int IHI, float *A, const int LDA, float *TAU, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for generating a real orthogonal matrix Q defined by elementary reflectors. (SORGHR)
 
void ORGHR (const int N, const int ILO, const int IHI, double *A, const int LDA, double *TAU, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for generating a real orthogonal matrix Q defined by elementary reflectors. (DORGHR)
 
void ORMHR (const char SIDE, const char TRANS, const int M, const int N, const int ILO, const int IHI, const float *A, const int LDA, const float *TAU, float *C, const int LDC, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for applying an orthogonal matrix in-place (SORMHR)
 
void ORMHR (const char SIDE, const char TRANS, const int M, const int N, const int ILO, const int IHI, const double *A, const int LDA, const double *TAU, double *C, const int LDC, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper for applying an orthogonal matrix in-place (DORMHR)
 
void LARFT (const char DIRECT, const char STOREV, const int N, const int K, double *V, const int LDV, double *TAU, double *T, const int LDT) const
 Epetra_LAPACK for forming the triangular factor of a product of elementary Householder reflectors (SLARFT).
 
void LARFT (const char DIRECT, const char STOREV, const int N, const int K, float *V, const int LDV, float *TAU, float *T, const int LDT) const
 Epetra_LAPACK for forming the triangular factor of a product of elementary Householder reflectors (DLARFT).
 
void TREVC (const char SIDE, const char HOWMNY, int *SELECT, const int N, const float *T, const int LDT, float *VL, const int LDVL, float *VR, const int LDVR, const int MM, int *M, float *WORK, int *INFO) const
 Epetra_LAPACK wrapper for computing eigenvectors of a quasi-triangular/triagnular matrix (STREVC) More...
 
void TREVC (const char SIDE, const char HOWMNY, int *SELECT, const int N, const double *T, const int LDT, double *VL, const int LDVL, double *VR, const int LDVR, const int MM, int *M, double *WORK, int *INFO) const
 Epetra_LAPACK wrapper for computing eigenvectors of a quasi-triangular/triagnular matrix (DTREVC) More...
 
void TREXC (const char COMPQ, const int N, float *T, const int LDT, float *Q, const int LDQ, int IFST, int ILST, float *WORK, int *INFO) const
 Epetra_LAPACK wrapper for reordering the real-Schur/Schur factorization of a matrix (STREXC)
 
void TREXC (const char COMPQ, const int N, double *T, const int LDT, double *Q, const int LDQ, int IFST, int ILST, double *WORK, int *INFO) const
 Epetra_LAPACK wrapper for reordering the real-Schur/Schur factorization of a matrix (DTREXC)
 
void GESVD (const char JOBU, const char JOBVT, const int M, const int N, float *A, const int LDA, float *S, float *U, const int LDU, float *VT, const int LDVT, float *WORK, const int *LWORK, int *INFO) const
 Epetra_LAPACK wrapper for computing the singular value decomposition (SGESVD)
 
void GESVD (const char JOBU, const char JOBVT, const int M, const int N, double *A, const int LDA, double *S, double *U, const int LDU, double *VT, const int LDVT, double *WORK, const int *LWORK, int *INFO) const
 Epetra_LAPACK wrapper for computing the singular value decomposition (DGESVD)
 
void GGSVD (const char JOBU, const char JOBV, const char JOBQ, const int M, const int N, const int P, int *K, int *L, double *A, const int LDA, double *B, const int LDB, double *ALPHA, double *BETA, double *U, const int LDU, double *V, const int LDV, double *Q, const int LDQ, double *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute the generalized singular value decomposition (GSVD) of an M-by-N real matrix A and P-by-N real matrix B.
 
void GGSVD (const char JOBU, const char JOBV, const char JOBQ, const int M, const int N, const int P, int *K, int *L, float *A, const int LDA, float *B, const int LDB, float *ALPHA, float *BETA, float *U, const int LDU, float *V, const int LDV, float *Q, const int LDQ, float *WORK, int *IWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute the generalized singular value decomposition (GSVD) of an M-by-N real matrix A and P-by-N real matrix B.
 
void GEEV (const char JOBVL, const char JOBVR, const int N, double *A, const int LDA, double *WR, double *WI, double *VL, const int LDVL, double *VR, const int LDVR, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
 
void GEEV (const char JOBVL, const char JOBVR, const int N, float *A, const int LDA, float *WR, float *WI, float *VL, const int LDVL, float *VR, const int LDVR, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
 
void SPEV (const char JOBZ, const char UPLO, const int N, double *AP, double *W, double *Z, int LDZ, double *WORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all the eigenvalues and, optionally, eigenvectors of a real symmetric matrix A in packed storage.
 
void SPEV (const char JOBZ, const char UPLO, const int N, float *AP, float *W, float *Z, int LDZ, float *WORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all the eigenvalues and, optionally, eigenvectors of a real symmetric matrix A in packed storage.
 
void SPGV (const int ITYPE, const char JOBZ, const char UPLO, const int N, double *AP, double *BP, double *W, double *Z, const int LDZ, double *WORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all the eigenvalues and, optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
 
void SPGV (const int ITYPE, const char JOBZ, const char UPLO, const int N, float *AP, float *BP, float *W, float *Z, const int LDZ, float *WORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all the eigenvalues and, optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
 
void SYEV (const char JOBZ, const char UPLO, const int N, double *A, const int LDA, double *W, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
 
void SYEV (const char JOBZ, const char UPLO, const int N, float *A, const int LDA, float *W, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
 
void SYEVD (const char JOBZ, const char UPLO, const int N, double *A, const int LDA, double *W, double *WORK, const int LWORK, int *IWORK, const int LIWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
 
void SYEVD (const char JOBZ, const char UPLO, const int N, float *A, const int LDA, float *W, float *WORK, const int LWORK, int *IWORK, const int LIWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
 
void SYEVX (const char JOBZ, const char RANGE, const char UPLO, const int N, double *A, const int LDA, const double *VL, const double *VU, const int *IL, const int *IU, const double ABSTOL, int *M, double *W, double *Z, const int LDZ, double *WORK, const int LWORK, int *IWORK, int *IFAIL, int *INFO) const
 Epetra_LAPACK wrapper to compute selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
 
void SYEVX (const char JOBZ, const char RANGE, const char UPLO, const int N, float *A, const int LDA, const float *VL, const float *VU, const int *IL, const int *IU, const float ABSTOL, int *M, float *W, float *Z, const int LDZ, float *WORK, const int LWORK, int *IWORK, int *IFAIL, int *INFO) const
 Epetra_LAPACK wrapper to compute selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
 
void SYGV (const int ITYPE, const char JOBZ, const char UPLO, const int N, double *A, const int LDA, double *B, const int LDB, double *W, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
 
void SYGV (const int ITYPE, const char JOBZ, const char UPLO, const int N, float *A, const int LDA, float *B, const int LDB, float *W, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
 
void SYGVX (const int ITYPE, const char JOBZ, const char RANGE, const char UPLO, const int N, double *A, const int LDA, double *B, const int LDB, const double *VL, const double *VU, const int *IL, const int *IU, const double ABSTOL, int *M, double *W, double *Z, const int LDZ, double *WORK, const int LWORK, int *IWORK, int *IFAIL, int *INFO) const
 Epetra_LAPACK wrapper to compute selected eigenvalues, and optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
 
void SYGVX (const int ITYPE, const char JOBZ, const char RANGE, const char UPLO, const int N, float *A, const int LDA, float *B, const int LDB, const float *VL, const float *VU, const int *IL, const int *IU, const float ABSTOL, int *M, float *W, float *Z, const int LDZ, float *WORK, const int LWORK, int *IWORK, int *IFAIL, int *INFO) const
 Epetra_LAPACK wrapper to compute selected eigenvalues, and optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
 
void SYEVR (const char JOBZ, const char RANGE, const char UPLO, const int N, double *A, const int LDA, const double *VL, const double *VU, const int *IL, const int *IU, const double ABSTOL, int *M, double *W, double *Z, const int LDZ, int *ISUPPZ, double *WORK, const int LWORK, int *IWORK, const int LIWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix T.
 
void SYEVR (const char JOBZ, const char RANGE, const char UPLO, const int N, float *A, const int LDA, const float *VL, const float *VU, const int *IL, const int *IU, const float ABSTOL, int *M, float *W, float *Z, const int LDZ, int *ISUPPZ, float *WORK, const int LWORK, int *IWORK, const int LIWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix T.
 
void GEEVX (const char BALANC, const char JOBVL, const char JOBVR, const char SENSE, const int N, double *A, const int LDA, double *WR, double *WI, double *VL, const int LDVL, double *VR, const int LDVR, int *ILO, int *IHI, double *SCALE, double *ABNRM, double *RCONDE, double *RCONDV, double *WORK, const int LWORK, int *IWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
 
void GEEVX (const char BALANC, const char JOBVL, const char JOBVR, const char SENSE, const int N, float *A, const int LDA, float *WR, float *WI, float *VL, const int LDVL, float *VR, const int LDVR, int *ILO, int *IHI, float *SCALE, float *ABNRM, float *RCONDE, float *RCONDV, float *WORK, const int LWORK, int *IWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
 
void GESDD (const char JOBZ, const int M, const int N, double *A, const int LDA, double *S, double *U, const int LDU, double *VT, const int LDVT, double *WORK, const int LWORK, int *IWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute the singular value decomposition (SVD) of a real M-by-N matrix A, optionally computing the left and right singular vectors.
 
void GESDD (const char JOBZ, const int M, const int N, float *A, const int LDA, float *S, float *U, const int LDU, float *VT, const int LDVT, float *WORK, const int LWORK, int *IWORK, int *INFO) const
 Epetra_LAPACK wrapper to.
 
void GGEV (const char JOBVL, const char JOBVR, const int N, double *A, const int LDA, double *B, const int LDB, double *ALPHAR, double *ALPHAI, double *BETA, double *VL, const int LDVL, double *VR, const int LDVR, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute for a pair of N-by-N real nonsymmetric matrices (A,B) the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors.
 
void GGEV (const char JOBVL, const char JOBVR, const int N, float *A, const int LDA, float *B, const int LDB, float *ALPHAR, float *ALPHAI, float *BETA, float *VL, const int LDVL, float *VR, const int LDVR, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to compute for a pair of N-by-N real nonsymmetric matrices (A,B) the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors.
 
void GGLSE (const int M, const int N, const int P, double *A, const int LDA, double *B, const int LDB, double *C, double *D, double *X, double *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to solve the linear equality-constrained least squares (LSE) problem.
 
void GGLSE (const int M, const int N, const int P, float *A, const int LDA, float *B, const int LDB, float *C, float *D, float *X, float *WORK, const int LWORK, int *INFO) const
 Epetra_LAPACK wrapper to solve the linear equality-constrained least squares (LSE) problem.
 
void LAMCH (const char CMACH, float &T) const
 Epetra_LAPACK wrapper for DLAMCH routine. On out, T holds machine double precision floating point characteristics. This information is returned by the Lapack routine.
 
void LAMCH (const char CMACH, double &T) const
 Epetra_LAPACK wrapper for SLAMCH routine. On out, T holds machine single precision floating point characteristics. This information is returned by the Lapack routine.
 
void TRTRS (const char UPLO, const char TRANS, const char DIAG, const int N, const int NRHS, const float *A, const int LDA, float *B, const int LDB, int *INFO) const
 Epetra_LAPACK wrapper for TRTRS routine.
 
void TRTRS (const char UPLO, const char TRANS, const char DIAG, const int N, const int NRHS, const double *A, const int LDA, double *B, const int LDB, int *INFO) const
 Epetra_LAPACK wrapper for TRTRS routine.
 

Public Attributes

bool Transpose_
 
bool Factored_
 
bool Solved_
 
bool Inverted_
 
char TRANS_
 
int M_
 
int N_
 
int Min_MN_
 
int NRHS_
 
int LDA_
 
int LDAI_
 
int LDB_
 
int LDX_
 
int INFO_
 
int LWORK_
 
int * IWORK_
 
double ANORM_
 
Epetra_SerialDenseMatrixMatrix_
 
Epetra_SerialDenseMatrixLHS_
 
Epetra_SerialDenseMatrixRHS_
 
Epetra_SerialDenseMatrixInverse_
 
double * A_
 
double * AI_
 
double * WORK_
 
double * U_
 
double * S_
 
double * Vt_
 
double * B_
 
double * X_
 
bool UseTranspose_
 

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_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_CompObject
Epetra_FlopsFlopCounter_
 

Detailed Description

Epetra_SerialDenseSVD: A class for SVDing dense linear problems.

The Epetra_SerialDenseSVD class enables the definition, in terms of Epetra_SerialDenseMatrix
and Epetra_SerialDenseVector objects, of a dense linear problem, followed by the solution of that problem via the
most sophisticated techniques available in LAPACK.

The Epetra_SerialDenseSVD class is intended to provide full-featured support for solving linear problems for general dense rectangular (or square) 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_SerialDenseSVD vs. Epetra_LAPACK

The Epetra_LAPACK class provides access to most of the same functionality as Epetra_SerialDenseSolver. The primary difference is that Epetra_LAPACK is a "thin" layer on top of LAPACK and Epetra_SerialDenseSolver attempts to provide easy access to the more sophisticated aspects of solving dense linear and eigensystems.

Constructing Epetra_SerialDenseSVD Objects

There is a single Epetra_SerialDenseSVD constructor. However, the matrix, right hand side and solution vectors must be set prior to executing most methods in this class.

Setting vectors used for linear solves

The matrix A, the left hand side X and the right hand side B (when solving AX = B, for X), can be set by appropriate set methods. Each of these three objects must be an Epetra_SerialDenseMatrix or and Epetra_SerialDenseVector object. The set methods are as follows:

Vector and Utility Functions

Once a Epetra_SerialDenseSVD is constructed, several mathematical functions can be applied to the object. Specifically:

Counting floating point operations The Epetra_SerialDenseSVD 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.

Examples using Epetra_SerialDenseSVD can be found in the Epetra test directories.

Member Function Documentation

virtual int Epetra_SerialDenseSVD::Apply ( const Epetra_SerialDenseMatrix Xmat,
Epetra_SerialDenseMatrix Ymat 
)
inlinevirtual

Returns the result of a Epetra_SerialDenseOperator applied to a Epetra_SerialDenseMatrix X in Y.

Parameters
InX - A Epetra_SerialDenseMatrix to multiply with operator.
OutY -A Epetra_SerialDenseMatrix containing result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_SerialDenseOperator.

References Epetra_SerialDenseMatrix::Multiply().

virtual int Epetra_SerialDenseSVD::ApplyInverse ( const Epetra_SerialDenseMatrix Xmat,
Epetra_SerialDenseMatrix Ymat 
)
inlinevirtual

Returns the result of a Epetra_SerialDenseOperator inverse applied to an Epetra_SerialDenseMatrix X in Y.

Parameters
InX - A Epetra_SerialDenseMatrix to solve for.
OutY -A Epetra_SerialDenseMatrix containing result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_SerialDenseOperator.

virtual int Epetra_SerialDenseSVD::Invert ( double  rthresh = 0.0,
double  athresh = 0.0 
)
virtual

Inverts the this matrix.

Returns
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.
virtual int Epetra_SerialDenseSVD::SetUseTranspose ( bool  use_transpose)
inlinevirtual

If set true, transpose of this operator will be applied.

This flag allows the transpose of the given operator to be used implicitly.  Setting this flag
affects only the Apply() and ApplyInverse() methods.  If the implementation of this interface

does not support transpose use, this method should return a value of -1.

Parameters
Inuse_transpose -If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Integer error code, set to 0 if successful. Set to -1 if this implementation does not support transpose.

Implements Epetra_SerialDenseOperator.

int Epetra_SerialDenseSVD::SetVectors ( Epetra_SerialDenseMatrix X,
Epetra_SerialDenseMatrix B 
)

Sets the pointers for left and right hand side vector(s).

Row dimension of X must match column dimension of matrix A, row dimension of B must match row dimension of A. X and B must have the same dimensions.

virtual int Epetra_SerialDenseSVD::Solve ( void  )
virtual

Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..

Inverse of Matrix must be formed

Returns
Integer error code, set to 0 if successful.
void Epetra_SerialDenseSVD::SolveWithTranspose ( bool  Flag)
inline

Causes equilibration to be called just before the matrix factorization as part of the call to Factor.

This function must be called before the factorization is performed.If Flag is true, causes all subsequent function calls to work with the transpose of this matrix, otherwise not.


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