Epetra_SerialDenseSolver Class Reference

Epetra_SerialDenseSolver: A class for solving dense linear problems. More...

#include <Epetra_SerialDenseSolver.h>

Inheritance diagram for Epetra_SerialDenseSolver:

Inheritance graph
[legend]
Collaboration diagram for Epetra_SerialDenseSolver:

Collaboration graph
[legend]
List of all members.

Public Member Functions

Constructor/Destructor Methods
 Epetra_SerialDenseSolver ()
 Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors().
virtual ~Epetra_SerialDenseSolver ()
 Epetra_SerialDenseSolver 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).
Strategy modifying Methods
void FactorWithEquilibration (bool Flag)
 Causes equilibration to be called just before the matrix factorization as part of the call to Factor.
void SolveWithTranspose (bool Flag)
 If Flag is true, causes all subsequent function calls to work with the transpose of this matrix, otherwise not.
void SolveToRefinedSolution (bool Flag)
 Causes all solves to compute solution to best ability using iterative refinement.
void EstimateSolutionErrors (bool Flag)
 Causes all solves to estimate the forward and backward solution error.
Factor/Solve/Invert Methods
virtual int Factor (void)
 Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF.
virtual int Solve (void)
 Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
virtual int Invert (void)
 Inverts the this matrix.
virtual int ComputeEquilibrateScaling (void)
 Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix.
virtual int EquilibrateMatrix (void)
 Equilibrates the this matrix.
int EquilibrateRHS (void)
 Equilibrates the current RHS.
virtual int ApplyRefinement (void)
 Apply Iterative Refinement.
int UnequilibrateLHS (void)
 Unscales the solution vectors if equilibration was used to solve the system.
virtual int ReciprocalConditionEstimate (double &Value)
 Returns the reciprocal of the 1-norm condition number of the this matrix.
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 A_Equilibrated ()
 Returns true if factor is equilibrated (factor available via AF() and LDAF()).
bool B_Equilibrated ()
 Returns true if RHS is equilibrated (RHS available via B() and LDB()).
virtual bool ShouldEquilibrate ()
 Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
bool SolutionErrorsEstimated ()
 Returns true if forward and backward error estimated have been computed (available via FERR() and BERR()).
bool Inverted ()
 Returns true if matrix inverse has been computed (inverse available via AF() and LDAF()).
bool ReciprocalConditionEstimated ()
 Returns true if the condition number of the this matrix has been computed (value available via ReciprocalConditionEstimate()).
bool Solved ()
 Returns true if the current set of vectors has been solved.
bool SolutionRefined ()
 Returns true if the current set of vectors has been refined.
Data Accessor methods
Epetra_SerialDenseMatrixMatrix () const
 Returns pointer to current matrix.
Epetra_SerialDenseMatrixFactoredMatrix () const
 Returns pointer to factored matrix (assuming factorization 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 * AF () const
 Returns pointer to the factored matrix (may be the same as A() if factorization done in place).
int LDAF () const
 Returns the leading dimension of the factored matrix.
int * IPIV () const
 Returns pointer to pivot vector (if factorization has been computed), zero otherwise.
double ANORM () const
 Returns the 1-Norm of the this matrix (returns -1 if not yet computed).
double RCOND () const
 Returns the reciprocal of the condition number of the this matrix (returns -1 if not yet computed).
double ROWCND () const
 Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).
double COLCND () const
 Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed).
double AMAX () const
 Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).
double * FERR () const
 Returns a pointer to the forward error estimates computed by LAPACK.
double * BERR () const
 Returns a pointer to the backward error estimates computed by LAPACK.
double * R () const
 Returns a pointer to the row scaling vector used for equilibration.
double * C () const
 Returns a pointer to the column scale vector used for equilibration.
I/O methods
virtual void Print (ostream &os) const
 Print service methods; defines behavior of ostream << operator.

Protected Member Functions

void AllocateWORK ()
void AllocateIWORK ()
void InitPointers ()
void DeleteArrays ()
void ResetMatrix ()
void ResetVectors ()

Protected Attributes

bool Equilibrate_
bool ShouldEquilibrate_
bool A_Equilibrated_
bool B_Equilibrated_
bool Transpose_
bool Factored_
bool EstimateSolutionErrors_
bool SolutionErrorsEstimated_
bool Solved_
bool Inverted_
bool ReciprocalConditionEstimated_
bool RefineSolution_
bool SolutionRefined_
char TRANS_
int M_
int N_
int Min_MN_
int NRHS_
int LDA_
int LDAF_
int LDB_
int LDX_
int INFO_
int LWORK_
int * IPIV_
int * IWORK_
double ANORM_
double RCOND_
double ROWCND_
double COLCND_
double AMAX_
Epetra_SerialDenseMatrixMatrix_
Epetra_SerialDenseMatrixLHS_
Epetra_SerialDenseMatrixRHS_
Epetra_SerialDenseMatrixFactor_
double * A_
double * FERR_
double * BERR_
double * AF_
double * WORK_
double * R_
double * C_
double * B_
double * X_

Detailed Description

Epetra_SerialDenseSolver: A class for solving dense linear problems.

The Epetra_SerialDenseSolver 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_SerialDenseSolver 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_SerialDenseSolver 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_SerialDenseSolver Objects

There is a single Epetra_SerialDenseSolver 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_SerialDenseSolver is constructed, several mathematical functions can be applied to the object. Specifically:

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

Strategies for Solving Linear Systems In many cases, linear systems can be accurately solved by simply computing the LU factorization of the matrix and then performing a forward back solve with a given set of right hand side vectors. However, in some instances, the factorization may be very poorly conditioned and this simple approach may not work. In these situations, equilibration and iterative refinement may improve the accuracy, or prevent a breakdown in the factorization.

Epetra_SerialDenseSolver will use equilibration with the factorization if, once the object is constructed and before it is factored, you call the function FactorWithEquilibration(true) to force equilibration to be used. If you are uncertain if equilibration should be used, you may call the function ShouldEquilibrate() which will return true if equilibration could possibly help. ShouldEquilibrate() uses guidelines specified in the LAPACK User Guide, namely if SCOND < 0.1 and AMAX < Underflow or AMAX > Overflow, to determine if equilibration might be useful.

Epetra_SerialDenseSolver will use iterative refinement after a forward/back solve if you call SolveToRefinedSolution(true). It will also compute forward and backward error estimates if you call EstimateSolutionErrors(true). Access to the forward (back) error estimates is available via FERR() (BERR()).

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


Member Function Documentation

virtual int Epetra_SerialDenseSolver::ApplyRefinement ( void   )  [virtual]

Apply Iterative Refinement.

Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.

double Epetra_SerialDenseSolver::COLCND (  )  const [inline]

Ratio of smallest to largest column scale factors for the this matrix (returns -1 if not yet computed).

If COLCND() is >= 0.1 then equilibration is not needed.

virtual int Epetra_SerialDenseSolver::ComputeEquilibrateScaling ( void   )  [virtual]

Computes the scaling vector S(i) = 1/sqrt(A(i,i)) of the this matrix.

Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.

virtual int Epetra_SerialDenseSolver::EquilibrateMatrix ( void   )  [virtual]

Equilibrates the this matrix.

Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.

int Epetra_SerialDenseSolver::EquilibrateRHS ( void   ) 

Equilibrates the current RHS.

Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.

void Epetra_SerialDenseSolver::EstimateSolutionErrors ( bool  Flag  ) 

Causes all solves to estimate the forward and backward solution error.

Error estimates will be in the arrays FERR and BERR, resp, after the solve step is complete. These arrays are accessible via the FERR() and BERR() access functions.

virtual int Epetra_SerialDenseSolver::Factor ( void   )  [virtual]

Computes the in-place LU factorization of the matrix using the LAPACK routine DGETRF.

Returns:
Integer error code, set to 0 if successful.

Reimplemented in Epetra_SerialSpdDenseSolver.

void Epetra_SerialDenseSolver::FactorWithEquilibration ( 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.

virtual int Epetra_SerialDenseSolver::Invert ( void   )  [virtual]

Inverts the this matrix.

Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.

virtual int Epetra_SerialDenseSolver::ReciprocalConditionEstimate ( double &  Value  )  [virtual]

Returns the reciprocal of the 1-norm condition number of the this matrix.

Parameters:
Value Out On return contains the reciprocal of the 1-norm condition number of the this matrix.
Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.

double Epetra_SerialDenseSolver::ROWCND (  )  const [inline]

Ratio of smallest to largest row scale factors for the this matrix (returns -1 if not yet computed).

If ROWCND() is >= 0.1 and AMAX() is not close to overflow or underflow, then equilibration is not needed.

int Epetra_SerialDenseSolver::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_SerialDenseSolver::Solve ( void   )  [virtual]

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

Returns:
Integer error code, set to 0 if successful.

Reimplemented in Epetra_SerialSpdDenseSolver.

int Epetra_SerialDenseSolver::UnequilibrateLHS ( void   ) 

Unscales the solution vectors if equilibration was used to solve the system.

Returns:
Integer error code, set to 0 if successful. Otherwise returns the LAPACK error code INFO.

Reimplemented in Epetra_SerialSpdDenseSolver.


The documentation for this class was generated from the following file:
Generated on Tue Jul 13 09:23:31 2010 for Epetra by  doxygen 1.4.7