Next: , Previous: Complex Hermitian Matrices, Up: Eigensystems


15.3 Real Nonsymmetric Matrices

The solution of the real nonsymmetric eigensystem problem for a matrix A involves computing the Schur decomposition

     A = Z T Z^T

where Z is an orthogonal matrix of Schur vectors and T, the Schur form, is quasi upper triangular with diagonal 1-by-1 blocks which are real eigenvalues of A, and diagonal 2-by-2 blocks whose eigenvalues are complex conjugate eigenvalues of A. The algorithm used is the double-shift Francis method.

— Function: gsl_eigen_nonsymm_workspace * gsl_eigen_nonsymm_alloc (const size_t n)

This function allocates a workspace for computing eigenvalues of n-by-n real nonsymmetric matrices. The size of the workspace is O(2n).

— Function: void gsl_eigen_nonsymm_free (gsl_eigen_nonsymm_workspace * w)

This function frees the memory associated with the workspace w.

— Function: void gsl_eigen_nonsymm_params (const int compute_t, const int balance, gsl_eigen_nonsymm_workspace * w)

This function sets some parameters which determine how the eigenvalue problem is solved in subsequent calls to gsl_eigen_nonsymm.

If compute_t is set to 1, the full Schur form T will be computed by gsl_eigen_nonsymm. If it is set to 0, T will not be computed (this is the default setting). Computing the full Schur form T requires approximately 1.5–2 times the number of flops.

If balance is set to 1, a balancing transformation is applied to the matrix prior to computing eigenvalues. This transformation is designed to make the rows and columns of the matrix have comparable norms, and can result in more accurate eigenvalues for matrices whose entries vary widely in magnitude. See Balancing for more information. Note that the balancing transformation does not preserve the orthogonality of the Schur vectors, so if you wish to compute the Schur vectors with gsl_eigen_nonsymm_Z you will obtain the Schur vectors of the balanced matrix instead of the original matrix. The relationship will be

          T = Q^t D^(-1) A D Q

where Q is the matrix of Schur vectors for the balanced matrix, and D is the balancing transformation. Then gsl_eigen_nonsymm_Z will compute a matrix Z which satisfies

          T = Z^(-1) A Z

with Z = D Q. Note that Z will not be orthogonal. For this reason, balancing is not performed by default.

— Function: int gsl_eigen_nonsymm (gsl_matrix * A, gsl_vector_complex * eval, gsl_eigen_nonsymm_workspace * w)

This function computes the eigenvalues of the real nonsymmetric matrix A and stores them in the vector eval. If T is desired, it is stored in the upper portion of A on output. Otherwise, on output, the diagonal of A will contain the 1-by-1 real eigenvalues and 2-by-2 complex conjugate eigenvalue systems, and the rest of A is destroyed. In rare cases, this function may fail to find all eigenvalues. If this happens, an error code is returned and the number of converged eigenvalues is stored in w->n_evals. The converged eigenvalues are stored in the beginning of eval.

— Function: int gsl_eigen_nonsymm_Z (gsl_matrix * A, gsl_vector_complex * eval, gsl_matrix * Z, gsl_eigen_nonsymm_workspace * w)

This function is identical to gsl_eigen_nonsymm except that it also computes the Schur vectors and stores them into Z.

— Function: gsl_eigen_nonsymmv_workspace * gsl_eigen_nonsymmv_alloc (const size_t n)

This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n real nonsymmetric matrices. The size of the workspace is O(5n).

— Function: void gsl_eigen_nonsymmv_free (gsl_eigen_nonsymmv_workspace * w)

This function frees the memory associated with the workspace w.

— Function: void gsl_eigen_nonsymmv_params (const int balance, gsl_eigen_nonsymm_workspace * w)

This function sets parameters which determine how the eigenvalue problem is solved in subsequent calls to gsl_eigen_nonsymmv. If balance is set to 1, a balancing transformation is applied to the matrix. See gsl_eigen_nonsymm_params for more information. Balancing is turned off by default since it does not preserve the orthogonality of the Schur vectors.

— Function: int gsl_eigen_nonsymmv (gsl_matrix * A, gsl_vector_complex * eval, gsl_matrix_complex * evec, gsl_eigen_nonsymmv_workspace * w)

This function computes eigenvalues and right eigenvectors of the n-by-n real nonsymmetric matrix A. It first calls gsl_eigen_nonsymm to compute the eigenvalues, Schur form T, and Schur vectors. Then it finds eigenvectors of T and backtransforms them using the Schur vectors. The Schur vectors are destroyed in the process, but can be saved by using gsl_eigen_nonsymmv_Z. The computed eigenvectors are normalized to have unit magnitude. On output, the upper portion of A contains the Schur form T. If gsl_eigen_nonsymm fails, no eigenvectors are computed, and an error code is returned.

— Function: int gsl_eigen_nonsymmv_Z (gsl_matrix * A, gsl_vector_complex * eval, gsl_matrix_complex * evec, gsl_matrix * Z, gsl_eigen_nonsymmv_workspace * w)

This function is identical to gsl_eigen_nonsymmv except that it also saves the Schur vectors into Z.