Epetra_SerialSpdDenseSolver Class Reference

Epetra_SerialSpdDenseSolver: A class for constructing and using symmetric positive definite dense matrices. More...

#include <Epetra_SerialSpdDenseSolver.h>

Inheritance diagram for Epetra_SerialSpdDenseSolver:

[legend]
Collaboration diagram for Epetra_SerialSpdDenseSolver:
[legend]
List of all members.

Public Member Functions

Constructor/Destructor Methods
 Epetra_SerialSpdDenseSolver ()
 Default constructor; matrix should be set using SetMatrix(), LHS and RHS set with SetVectors().
virtual ~Epetra_SerialSpdDenseSolver ()
 Epetra_SerialDenseSolver destructor.
Set Methods
int SetMatrix (Epetra_SerialSymDenseMatrix &A)
 Sets the pointers for coefficient matrix; special version for symmetric matrices.
Factor/Solve/Invert Methods
int Factor (void)
 Computes the in-place Cholesky factorization of the matrix using the LAPACK routine DPOTRF.
int Solve (void)
 Computes the solution X to AX = B for the this matrix and the B provided to SetVectors()..
int Invert (void)
 Inverts the this matrix.
int ComputeEquilibrateScaling (void)
 Computes the scaling vector S(i) = 1/sqrt(A(i,i) of the this matrix.
int EquilibrateMatrix (void)
 Equilibrates the this matrix.
int EquilibrateRHS (void)
 Equilibrates the current RHS.
int ApplyRefinement (void)
 Apply Iterative Refinement.
int UnequilibrateLHS (void)
 Unscales the solution vectors if equilibration was used to solve the system.
int ReciprocalConditionEstimate (double &Value)
 Returns the reciprocal of the 1-norm condition number of the this matrix.
Query methods
bool ShouldEquilibrate ()
 Returns true if the LAPACK general rules for equilibration suggest you should equilibrate the system.
Data Accessor methods
Epetra_SerialSymDenseMatrixSymMatrix () const
 Returns pointer to current matrix.
Epetra_SerialSymDenseMatrixSymFactoredMatrix () const
 Returns pointer to factored matrix (assuming factorization has been performed).
double SCOND ()
 Ratio of smallest to largest equilibration scale factors for the this matrix (returns -1 if not yet computed).
double AMAX ()
 Returns the absolute value of the largest entry of the this matrix (returns -1 if not yet computed).

Detailed Description

Epetra_SerialSpdDenseSolver: A class for constructing and using symmetric positive definite dense matrices.

The Epetra_SerialSpdDenseSolver class enables the construction and use of real-valued, symmetric positive definite, double-precision dense matrices. It is built on the Epetra_DenseMatrix class which in turn is built on the BLAS and LAPACK via the Epetra_BLAS and Epetra_LAPACK classes.

The Epetra_SerialSpdDenseSolver 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_SerialSpdDenseSolver vs. Epetra_LAPACK

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

Constructing Epetra_SerialSpdDenseSolver 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):

  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.
Setting vectors used for linear solves

Setting the X and B vectors (which are Epetra_DenseMatrix objects) used for solving linear systems is done separately from the constructor. This allows a single matrix factor to be used for multiple solves. Similar to the constructor, the vectors X and B can be copied or viewed using the Epetra_DataAccess argument.

Extracting Data from Epetra_SerialSpdDenseSolver Objects

Once a Epetra_SerialSpdDenseSolver 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.
Vector and Utility Functions

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

The final useful function is Flops(). Each Epetra_SerialSpdDenseSolver 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.

Strategies for Solving Linear Systems In many cases, linear systems can be accurately solved by simply computing the Cholesky 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 the simple approach may not work. In these situations, equilibration and iterative refinement may improve the accuracy, or prevent a breakdown in the factorization.

Epetra_SerialSpdDenseSolver 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_SerialSpdDenseSolver 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_SerialSpdDenseSolver can be found in the Epetra test directories.


Member Function Documentation

int Epetra_SerialSpdDenseSolver::ApplyRefinement void   )  [virtual]
 

Apply Iterative Refinement.

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

Reimplemented from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::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 from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::EquilibrateMatrix void   )  [virtual]
 

Equilibrates the this matrix.

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

Reimplemented from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::EquilibrateRHS void   ) 
 

Equilibrates the current RHS.

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

Reimplemented from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::Factor void   )  [virtual]
 

Computes the in-place Cholesky factorization of the matrix using the LAPACK routine DPOTRF.

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

Reimplemented from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::Invert void   )  [virtual]
 

Inverts the this matrix.

Note: This function works a little differently that DPOTRI in that it fills the entire matrix with the inverse, independent of the UPLO specification.

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

Reimplemented from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::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 from Epetra_SerialDenseSolver.

double Epetra_SerialSpdDenseSolver::SCOND  )  [inline]
 

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

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

int Epetra_SerialSpdDenseSolver::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 from Epetra_SerialDenseSolver.

int Epetra_SerialSpdDenseSolver::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 from Epetra_SerialDenseSolver.


The documentation for this class was generated from the following file:
Generated on Thu Sep 18 12:42:07 2008 for Epetra by doxygen 1.3.9.1