next up previous contents index
Next: ERRHANDL(3LAS) Up: Manual Pages Previous: INTRO(3LAS)

  EIGENVAL(3LAS)

 

NAME

SetEigenvalAccuracy, GetMaxEigenval, GetMinEigenval -- estimation of extremal eigenvalues

SYNOPSIS

#include <laspack/eigenval.h>

void SetEigenvalAccuracy(double Eps);
double GetMinEigenval(QMatrix *Q, PrecondProcType PrecondProc,
                      double OmegaPrecond);
double GetMaxEigenval(QMatrix *Q, PrecondProcType PrecondProc,
                      double OmegaPrecond);

DESCRIPTION

  
The procedures GetMinEigenval and GetMaxEigenval estimate extremal eigenvalues of the matrix Q. In LASPack both eigenvalues are required internally for the Chebyshev method. Otherwise they could be applied in an application, e.g. for the estimation of the condition number of a matrix (see example below).

The algorithm used is based on the Lanczos method for symmetric matrices. Eigenvalues of non-symmetric matrices are estimated by means of Q^T Q.

For singular matrices, if the null space has been specified, the zero eigenvalue is not considered.

If any procedure for preconditioning is given by the parameter PrecondProc (i.e. if PrecondProc != NULL), eigenvalues are computed for the preconditioned matrix. The procedure used should have the prototype

  Vector *(*PrecondProcType)(QMatrix *, Vector *, Vector *, double)
The parameter OmegaPrecond will be passed as a relaxation parameter to the preconditioner.

 
The accuracy for stopping criterion can be set by SetEigenvalAccuracy. The default is 1e-4.

Each of the procedures GetMinEigenval and GetMaxEigenval causes the estimation of both, minimum and maximum eigenvalues. They are stored together with the parameters PrecondProc and OmegaPrecond in connection with the matrix Q, and are returned at the next call of above routines if the same parameters are used. Therefore, the extremal eigenvalues need not to be stored additionally, but could be queried on demand without extra computational costs.

REFERENCES

The implementation of the Lanczos method follows algorithms described in:

G. Golub, C. van Loan: Matrix Computations, second edition, The Johns Hopkins University Press, Baltimore, 1989.

FILES

eigenval.h ... header file
eigenval.c ... source file

EXAMPLES

In the following code fragment, condition numbers of a matrix and their preconditioned variants are estimated.

QMatrix A;
double LambdaMax, LambdaMin;
size_t Dim;


Dim = ... ;
Q_Constr(&A, "A", Dim, True, Rowws, Normal, True);

/* generation of the matrix A */

...
    
printf("Condition numbers of the matrix %s:\n\n", Q_GetName(&A));

LambdaMax = GetMaxEigenval(&A, NULL, 1.0);
LambdaMin = GetMinEigenval(&A, NULL, 1.0);
printf("without preconditioning:      k = %12.5e\n", LambdaMax / LambdaMin); 

LambdaMax = GetMaxEigenval(&A, JacobiPrecond, 1.0);
LambdaMin = GetMinEigenval(&A, JacobiPrecond, 1.0);
printf("with Jacobi preconditioning:  k = %12.5e\n", LambdaMax / LambdaMin); 

LambdaMax = GetMaxEigenval(&A, SSORPrecond, 1.0);
LambdaMin = GetMinEigenval(&A, SSORPrecond, 1.0);
printf("with SSOR preconditioning:    k = %12.5e\n", LambdaMax / LambdaMin); 

LambdaMax = GetMaxEigenval(&A, ILUPrecond, 1.0);
LambdaMin = GetMinEigenval(&A, ILUPrecond, 1.0);
printf("with ILU preconditioning:     k = %12.5e\n", LambdaMax / LambdaMin); 
    
Q_Destr(&A);

SEE ALSO

qmatrix(3LAS), precond(3LAS)



next up previous contents index
Next: ERRHANDL(3LAS) Up: Manual Pages Previous: INTRO(3LAS)



Tomas Skalicky (skalicky@msmfs1.mw.tu-dresden.de)