dgeevx - compute for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors
SUBROUTINE DGEEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, WR, WI, * VL, LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONE, RCONV, WORK, * LDWORK, IWORK2, INFO) CHARACTER * 1 BALANC, JOBVL, JOBVR, SENSE INTEGER N, LDA, LDVL, LDVR, ILO, IHI, LDWORK, INFO INTEGER IWORK2(*) DOUBLE PRECISION ABNRM DOUBLE PRECISION A(LDA,*), WR(*), WI(*), VL(LDVL,*), VR(LDVR,*), SCALE(*), RCONE(*), RCONV(*), WORK(*)
SUBROUTINE DGEEVX_64( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, WR, * WI, VL, LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONE, RCONV, * WORK, LDWORK, IWORK2, INFO) CHARACTER * 1 BALANC, JOBVL, JOBVR, SENSE INTEGER*8 N, LDA, LDVL, LDVR, ILO, IHI, LDWORK, INFO INTEGER*8 IWORK2(*) DOUBLE PRECISION ABNRM DOUBLE PRECISION A(LDA,*), WR(*), WI(*), VL(LDVL,*), VR(LDVR,*), SCALE(*), RCONE(*), RCONV(*), WORK(*)
SUBROUTINE GEEVX( BALANC, JOBVL, JOBVR, SENSE, [N], A, [LDA], WR, * WI, VL, [LDVL], VR, [LDVR], ILO, IHI, SCALE, ABNRM, RCONE, RCONV, * WORK, [LDWORK], [IWORK2], [INFO]) CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE INTEGER :: N, LDA, LDVL, LDVR, ILO, IHI, LDWORK, INFO INTEGER, DIMENSION(:) :: IWORK2 REAL(8) :: ABNRM REAL(8), DIMENSION(:) :: WR, WI, SCALE, RCONE, RCONV, WORK REAL(8), DIMENSION(:,:) :: A, VL, VR
SUBROUTINE GEEVX_64( BALANC, JOBVL, JOBVR, SENSE, [N], A, [LDA], WR, * WI, VL, [LDVL], VR, [LDVR], ILO, IHI, SCALE, ABNRM, RCONE, RCONV, * WORK, [LDWORK], [IWORK2], [INFO]) CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE INTEGER(8) :: N, LDA, LDVL, LDVR, ILO, IHI, LDWORK, INFO INTEGER(8), DIMENSION(:) :: IWORK2 REAL(8) :: ABNRM REAL(8), DIMENSION(:) :: WR, WI, SCALE, RCONE, RCONV, WORK REAL(8), DIMENSION(:,:) :: A, VL, VR
#include <sunperf.h>
void dgeevx(char balanc, char jobvl, char jobvr, char sense, int n, double *a, int lda, double *wr, double *wi, double *vl, int ldvl, double *vr, int ldvr, int *ilo, int *ihi, double *scale, double *abnrm, double *rcone, double *rconv, double *work, int ldwork, int *info);
void dgeevx_64(char balanc, char jobvl, char jobvr, char sense, long n, double *a, long lda, double *wr, double *wi, double *vl, long ldvl, double *vr, long ldvr, long *ilo, long *ihi, double *scale, double *abnrm, double *rcone, double *rconv, double *work, long ldwork, long *info);
dgeevx computes for an N-by-N real nonsymmetric matrix A, the eigenvalues and, optionally, the left and/or right eigenvectors.
Optionally also, it computes a balancing transformation to improve the conditioning of the eigenvalues and eigenvectors (ILO, IHI, SCALE, and ABNRM), reciprocal condition numbers for the eigenvalues (RCONDE), and reciprocal condition numbers for the right
eigenvectors (RCONDV).
The right eigenvector v(j)
of A satisfies
A * v(j) = lambda(j) * v(j)
where lambda(j)
is its eigenvalue.
The left eigenvector u(j)
of A satisfies
u(j)**H * A = lambda(j) * u(j)**H
where u(j)**H denotes the conjugate transpose of u(j).
The computed eigenvectors are normalized to have Euclidean norm equal to 1 and largest component real.
Balancing a matrix means permuting the rows and columns to make it more nearly upper triangular, and applying a diagonal similarity transformation D * A * D**(-1), where D is a diagonal matrix, to make its rows and columns closer in norm and the condition numbers of its eigenvalues and eigenvectors smaller. The computed reciprocal condition numbers correspond to the balanced matrix. Permuting rows and columns will not change the condition numbers (in exact arithmetic) but diagonal scaling will. For further explanation of balancing, see section 4.10.2 of the LAPACK Users' Guide.
= 'P': Perform permutations to make the matrix more nearly upper triangular. Do not diagonally scale; = 'S': Diagonally scale the matrix, i.e. replace A by D*A*D**(-1), where D is a diagonal matrix chosen to make the rows and columns of A more equal in norm. Do not permute; = 'B': Both diagonally scale and permute A.
Computed reciprocal condition numbers will be for the matrix after balancing and/or permuting. Permuting does not change condition numbers (in exact arithmetic), but balancing does.
= 'N': left eigenvectors of A are not computed;
= 'V': left eigenvectors of A are computed. If SENSE = 'E' or 'B', JOBVL must = 'V'.
= 'N': right eigenvectors of A are not computed;
= 'V': right eigenvectors of A are computed. If SENSE = 'E' or 'B', JOBVR must = 'V'.
= 'E': Computed for eigenvalues only;
= 'V': Computed for right eigenvectors only;
= 'B': Computed for eigenvalues and right eigenvectors.
If SENSE = 'E' or 'B', both left and right eigenvectors must also be computed (JOBVL = 'V' and JOBVR = 'V').
u(j)
are stored one
after another in the columns of VL, in the same order
as their eigenvalues.
If JOBVL = 'N', VL is not referenced.
If the j-th eigenvalue is real, then u(j)
= VL(:,j),
the j-th column of VL.
If the j-th and (j+1)-st eigenvalues form a complex
conjugate pair, then u(j)
= VL(:,j)
+ i*VL(:,j+1) and
u(j+1)
= VL(:,j)
- i*VL(:,j+1).
v(j)
are stored one
after another in the columns of VR, in the same order
as their eigenvalues.
If JOBVR = 'N', VR is not referenced.
If the j-th eigenvalue is real, then v(j)
= VR(:,j),
the j-th column of VR.
If the j-th and (j+1)-st eigenvalues form a complex
conjugate pair, then v(j)
= VR(:,j)
+ i*VR(:,j+1) and
v(j+1)
= VR(:,j)
- i*VR(:,j+1).
A(i,j)
= 0 if I > J and
J = 1,...,ILO-1 or I = IHI+1,...,N.
P(j)
is the index of the row and column
interchanged with row and column j, and D(j)
is the scaling
factor applied to row and column j, then
SCALE(J)
= P(J), for J = 1,...,ILO-1
= D(J), for J = ILO,...,IHI
= P(J)
for J = IHI+1,...,N.
The order in which the interchanges are made is N to IHI+1,
then 1 to ILO-1.
RCONE(j)
is the reciprocal condition number of the j-th
eigenvalue.
RCONV(j)
is the reciprocal condition number of the j-th
right eigenvector.
WORK(1)
returns the optimal LDWORK.
If LDWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LDWORK is issued by XERBLA.
dimension(2*N-2)
If SENSE = 'N' or 'E', not referenced.
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value.
> 0: if INFO = i, the QR algorithm failed to compute all the eigenvalues, and no eigenvectors or condition numbers have been computed; elements 1:ILO-1 and i+1:N of WR and WI contain eigenvalues which have converged.