Google

Go to the first, previous, next, last section, table of contents.


Eigensystems

This chapter describes functions for computing eigenvalues and eigenvectors of matrices. There are routines for real symmetric and complex hermitian matrices, and eigenvalues can be computed with or without eigenvectors. The algorithms used are symmetric bidiagonalization followed by QR reduction.

These routines are intended for "small" systems where simple algorithms are acceptable. Anyone interested finding eigenvalues and eigenvectors of large matrices will want to use the sophisticated routines found in LAPACK. The Fortran version of LAPACK is recommended as the standard package for linear algebra.

The functions described in this chapter are declared in the header file `gsl_eigen.h'.

Real Symmetric Matrices

Function: gsl_eigen_symm_workspace * gsl_eigen_symm_alloc (const size_t n)
This function allocates a workspace for computing eigenvalues of n-by-n real symmetric matrices. The size of the workspace is O(2n).

Function: void gsl_eigen_symm_free (gsl_eigen_symm_workspace * w)
This function frees the memory associated with the workspace w.

Function: int gsl_eigen_symm (gsl_matrix * A, gsl_vector * eval, gsl_eigen_symm_workspace * w)
This function computes the eigenvalues of the real symmetric matrix A. Additional workspace of the appropriate size must be provided in w. The diagonal and lower triangular part of A are destroyed during the computation, but the strict upper triangular part is not referenced. The eigenvalues are stored in the vector eval and are unordered.

Function: gsl_eigen_symmv_workspace * gsl_eigen_symmv_alloc (const size_t n)
This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n real symmetric matrices. The size of the workspace is O(4n).

Function: void gsl_eigen_symmv_free (gsl_eigen_symmv_workspace * w)
This function frees the memory associated with the workspace w.

Function: int gsl_eigen_symmv (gsl_matrix * A, gsl_vector * eval, gsl_matrix * evec, gsl_eigen_symmv_workspace * w)
This function computes the eigenvalues and eigenvectors of the real symmetric matrix A. Additional workspace of the appropriate size must be provided in w. The diagonal and lower triangular part of A are destroyed during the computation, but the strict upper triangular part is not referenced. The eigenvalues are stored in the vector eval and are unordered. The corresponding eigenvectors are stored in the columns of the matrix evec. For example, the eigenvector in the first column corresponds to the first eigenvalue. The eigenvectors are guaranteed to be mutually orthogonal and normalised to unit magnitude.

Complex Hermitian Matrices

Function: gsl_eigen_herm_workspace * gsl_eigen_herm_alloc (const size_t n)
This function allocates a workspace for computing eigenvalues of n-by-n complex hermitian matrices. The size of the workspace is O(3n).

Function: void gsl_eigen_herm_free (gsl_eigen_herm_workspace * w)
This function frees the memory associated with the workspace w.

Function: int gsl_eigen_herm (gsl_matrix_complex * A, gsl_vector * eval, gsl_eigen_herm_workspace * w)
This function computes the eigenvalues of the complex hermitian matrix A. Additional workspace of the appropriate size must be provided in w. The diagonal and lower triangular part of A are destroyed during the computation, but the strict upper triangular part is not referenced. The imaginary parts of the diagonal are assumed to be zero and are not referenced. The eigenvalues are stored in the vector eval and are unordered.

Function: gsl_eigen_hermv_workspace * gsl_eigen_hermv_alloc (const size_t n)
This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n complex hermitian matrices. The size of the workspace is O(5n).

Function: void gsl_eigen_hermv_free (gsl_eigen_hermv_workspace * w)
This function frees the memory associated with the workspace w.

Function: int gsl_eigen_hermv (gsl_matrix_complex * A, gsl_vector * eval, gsl_matrix_complex * evec, gsl_eigen_hermv_workspace * w)
This function computes the eigenvalues and eigenvectors of the complex hermitian matrix A. Additional workspace of the appropriate size must be provided in w. The diagonal and lower triangular part of A are destroyed during the computation, but the strict upper triangular part is not referenced. The imaginary parts of the diagonal are assumed to be zero and are not referenced. The eigenvalues are stored in the vector eval and are unordered. The corresponding complex eigenvectors are stored in the columns of the matrix evec. For example, the eigenvector in the first column corresponds to the first eigenvalue. The eigenvectors are guaranteed to be mutually orthogonal and normalised to unit magnitude.

Sorting Eigenvalues and Eigenvectors

Function: int gsl_eigen_symmv_sort (gsl_vector * eval, gsl_matrix * evec, gsl_eigen_sort_t sort_type)
This function simultaneously sorts the eigenvalues stored in the vector eval and the corresponding real eigenvectors stored in the columns of the matrix evec into ascending or descending order according to the value of the parameter sort_type,

GSL_EIGEN_SORT_VAL_ASC
ascending order in numerical value
GSL_EIGEN_SORT_VAL_DESC
descending order in numerical value
GSL_EIGEN_SORT_ABS_ASC
ascending order in magnitude
GSL_EIGEN_SORT_ABS_DESC
descending order in magnitude

Function: int gsl_eigen_hermv_sort (gsl_vector * eval, gsl_matrix_complex * evec, gsl_eigen_sort_t sort_type)
This function simultaneously sorts the eigenvalues stored in the vector eval and the corresponding complex eigenvectors stored in the columns of the matrix evec into ascending or descending order according to the value of the parameter sort_type as shown above.

Examples

The following program computes the eigenvalues and eigenvectors of the 4-th order Hilbert matrix, H(i,j) = 1/(i + j + 1).

#include <stdio.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_eigen.h>

int
main (void)
{
  double data[] = { 1.0  , 1/2.0, 1/3.0, 1/4.0,
                    1/2.0, 1/3.0, 1/4.0, 1/5.0,
                    1/3.0, 1/4.0, 1/5.0, 1/6.0,
                    1/4.0, 1/5.0, 1/6.0, 1/7.0 };

  gsl_matrix_view m 
    = gsl_matrix_view_array(data, 4, 4);

  gsl_vector *eval = gsl_vector_alloc (4);
  gsl_matrix *evec = gsl_matrix_alloc (4, 4);

  gsl_eigen_symmv_workspace * w = 
    gsl_eigen_symmv_alloc (4);
  
  gsl_eigen_symmv (&m.matrix, eval, evec, w);

  gsl_eigen_symmv_free(w);

  gsl_eigen_symmv_sort (eval, evec, 
                        GSL_EIGEN_SORT_ABS_ASC);
  
  {
    int i;

    for (i = 0; i < 4; i++)
      {
        double eval_i 
           = gsl_vector_get(eval, i);
        gsl_vector_view evec_i 
           = gsl_matrix_column(evec, i);

        printf("eigenvalue = %g\n", eval_i);
        printf("eigenvector = \n");
        gsl_vector_fprintf(stdout, 
                           &evec_i.vector, "%g");
      }
  }

  return 0;
}

Here is the beginning of the output from the program,

$ ./a.out 
eigenvalue = 9.67023e-05
eigenvector = 
-0.0291933
0.328712
-0.791411
0.514553
...

This can be compared with the corresponding output from GNU OCTAVE,

octave> [v,d] = eig(hilb(4));
octave> diag(d)  
ans =

   9.6702e-05
   6.7383e-03
   1.6914e-01
   1.5002e+00

octave> v 
v =

   0.029193   0.179186  -0.582076   0.792608
  -0.328712  -0.741918   0.370502   0.451923
   0.791411   0.100228   0.509579   0.322416
  -0.514553   0.638283   0.514048   0.252161

Note that the eigenvectors can differ by a change of sign, since the sign of an eigenvector is arbitrary.

References and Further Reading

Further information on the algorithms described in this section can be found in the following book,

  • G. H. Golub, C. F. Van Loan, Matrix Computations (3rd Ed, 1996), Johns Hopkins University Press, ISBN 0-8018-5414-8.

The LAPACK library is described in,

The LAPACK source code can be found at the website above along with an online copy of the users guide.


Go to the first, previous, next, last section, table of contents.