zggevx


NAME

zggevx - compute for a pair of N-by-N complex nonsymmetric matrices (A,B) the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors


SYNOPSIS

  SUBROUTINE ZGGEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB, 
 *      ALPHA, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE, RSCALE, ABNRM, 
 *      BBNRM, RCONDE, RCONDV, WORK, LWORK, RWORK, IWORK, BWORK, INFO)
  CHARACTER * 1 BALANC, JOBVL, JOBVR, SENSE
  DOUBLE COMPLEX A(LDA,*), B(LDB,*), ALPHA(*), BETA(*), VL(LDVL,*), VR(LDVR,*), WORK(*)
  INTEGER N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
  INTEGER IWORK(*)
  LOGICAL BWORK(*)
  DOUBLE PRECISION ABNRM, BBNRM
  DOUBLE PRECISION LSCALE(*), RSCALE(*), RCONDE(*), RCONDV(*), RWORK(*)
 
  SUBROUTINE ZGGEVX_64( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, 
 *      LDB, ALPHA, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE, RSCALE, 
 *      ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, RWORK, IWORK, BWORK, 
 *      INFO)
  CHARACTER * 1 BALANC, JOBVL, JOBVR, SENSE
  DOUBLE COMPLEX A(LDA,*), B(LDB,*), ALPHA(*), BETA(*), VL(LDVL,*), VR(LDVR,*), WORK(*)
  INTEGER*8 N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
  INTEGER*8 IWORK(*)
  LOGICAL*8 BWORK(*)
  DOUBLE PRECISION ABNRM, BBNRM
  DOUBLE PRECISION LSCALE(*), RSCALE(*), RCONDE(*), RCONDV(*), RWORK(*)
 

F95 INTERFACE

  SUBROUTINE GGEVX( BALANC, JOBVL, JOBVR, SENSE, [N], A, [LDA], B, 
 *       [LDB], ALPHA, BETA, VL, [LDVL], VR, [LDVR], ILO, IHI, LSCALE, 
 *       RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, [WORK], [LWORK], [RWORK], 
 *       [IWORK], [BWORK], [INFO])
  CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE
  COMPLEX(8), DIMENSION(:) :: ALPHA, BETA, WORK
  COMPLEX(8), DIMENSION(:,:) :: A, B, VL, VR
  INTEGER :: N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
  INTEGER, DIMENSION(:) :: IWORK
  LOGICAL, DIMENSION(:) :: BWORK
  REAL(8) :: ABNRM, BBNRM
  REAL(8), DIMENSION(:) :: LSCALE, RSCALE, RCONDE, RCONDV, RWORK
 
  SUBROUTINE GGEVX_64( BALANC, JOBVL, JOBVR, SENSE, [N], A, [LDA], B, 
 *       [LDB], ALPHA, BETA, VL, [LDVL], VR, [LDVR], ILO, IHI, LSCALE, 
 *       RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, [WORK], [LWORK], [RWORK], 
 *       [IWORK], [BWORK], [INFO])
  CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE
  COMPLEX(8), DIMENSION(:) :: ALPHA, BETA, WORK
  COMPLEX(8), DIMENSION(:,:) :: A, B, VL, VR
  INTEGER(8) :: N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
  INTEGER(8), DIMENSION(:) :: IWORK
  LOGICAL(8), DIMENSION(:) :: BWORK
  REAL(8) :: ABNRM, BBNRM
  REAL(8), DIMENSION(:) :: LSCALE, RSCALE, RCONDE, RCONDV, RWORK
 

C INTERFACE

#include <sunperf.h>

void zggevx(char balanc, char jobvl, char jobvr, char sense, int n, doublecomplex *a, int lda, doublecomplex *b, int ldb, doublecomplex *alpha, doublecomplex *beta, doublecomplex *vl, int ldvl, doublecomplex *vr, int ldvr, int *ilo, int *ihi, double *lscale, double *rscale, double *abnrm, double *bbnrm, double *rconde, double *rcondv, int *info);

void zggevx_64(char balanc, char jobvl, char jobvr, char sense, long n, doublecomplex *a, long lda, doublecomplex *b, long ldb, doublecomplex *alpha, doublecomplex *beta, doublecomplex *vl, long ldvl, doublecomplex *vr, long ldvr, long *ilo, long *ihi, double *lscale, double *rscale, double *abnrm, double *bbnrm, double *rconde, double *rcondv, long *info);


PURPOSE

zggevx computes for a pair of N-by-N complex nonsymmetric matrices (A,B) the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors.

Optionally, it also computes a balancing transformation to improve the conditioning of the eigenvalues and eigenvectors (ILO, IHI, LSCALE, RSCALE, ABNRM, and BBNRM), reciprocal condition numbers for the eigenvalues (RCONDE), and reciprocal condition numbers for the right eigenvectors (RCONDV).

A generalized eigenvalue for a pair of matrices (A,B) is a scalar lambda or a ratio alpha/beta = lambda, such that A - lambda*B is singular. It is usually represented as the pair (alpha,beta), as there is a reasonable interpretation for beta=0, and even for both being zero.

The right eigenvector v(j) corresponding to the eigenvalue lambda(j) of (A,B) satisfies

                 A * v(j) = lambda(j) * B * v(j) .

The left eigenvector u(j) corresponding to the eigenvalue lambda(j) of (A,B) satisfies

                 u(j)**H * A  = lambda(j) * u(j)**H * B.

where u(j)**H is the conjugate-transpose of u(j).


ARGUMENTS

* BALANC (input)
Specifies the balance option to be performed:

* JOBVL (input)
* JOBVR (input)

* SENSE (input)
Determines which reciprocal condition numbers are computed.

* N (input)
The order of the matrices A, B, VL, and VR. N >= 0.

* A (input/output)
On entry, the matrix A in the pair (A,B). On exit, A has been overwritten. If JOBVL='V' or JOBVR='V' or both, then A contains the first part of the complex Schur form of the ``balanced'' versions of the input A and B.

* LDA (input)
The leading dimension of A. LDA >= max(1,N).

* B (input/output)
On entry, the matrix B in the pair (A,B). On exit, B has been overwritten. If JOBVL='V' or JOBVR='V' or both, then B contains the second part of the complex Schur form of the ``balanced'' versions of the input A and B.

* LDB (input)
The leading dimension of B. LDB >= max(1,N).

* ALPHA (output)
On exit, ALPHA(j)/BETA(j), j=1,...,N, will be the generalized eigenvalues.

Note: the quotient ALPHA(j)/BETA(j) ) may easily over- or underflow, and BETA(j) may even be zero. Thus, the user should avoid naively computing the ratio ALPHA/BETA. However, ALPHA will be always less than and usually comparable with norm(A) in magnitude, and BETA always less than and usually comparable with norm(B).

* BETA (output)
See description of ALPHA.

* VL (output)
If JOBVL = 'V', the left generalized eigenvectors u(j) are stored one after another in the columns of VL, in the same order as their eigenvalues. Each eigenvector will be scaled so the largest component will have abs(real part) + abs(imag. part) = 1. Not referenced if JOBVL = 'N'.

* LDVL (input)
The leading dimension of the matrix VL. LDVL >= 1, and if JOBVL = 'V', LDVL >= N.

* VR (output)
If JOBVR = 'V', the right generalized eigenvectors v(j) are stored one after another in the columns of VR, in the same order as their eigenvalues. Each eigenvector will be scaled so the largest component will have abs(real part) + abs(imag. part) = 1. Not referenced if JOBVR = 'N'.

* LDVR (input)
The leading dimension of the matrix VR. LDVR >= 1, and if JOBVR = 'V', LDVR >= N.

* ILO (output)
ILO is an integer value such that on exit A(i,j) = 0 and B(i,j) = 0 if i > j and j = 1,...,ILO-1 or i = IHI+1,...,N. If BALANC = 'N' or 'S', ILO = 1 and IHI = N.

* IHI (output)
IHI is an integer value such that on exit A(i,j) = 0 and B(i,j) = 0 if i > j and j = 1,...,ILO-1 or i = IHI+1,...,N. If BALANC = 'N' or 'S', ILO = 1 and IHI = N.

* LSCALE (output)
Details of the permutations and scaling factors applied to the left side of A and B. If PL(j) is the index of the row interchanged with row j, and DL(j) is the scaling factor applied to row j, then LSCALE(j) = PL(j) for j = 1,...,ILO-1

* RSCALE (output)
Details of the permutations and scaling factors applied to the right side of A and B. If PR(j) is the index of the column interchanged with column j, and DR(j) is the scaling factor applied to column j, then RSCALE(j) = PR(j) for j = 1,...,ILO-1

* ABNRM (output)
The one-norm of the balanced matrix A.

* BBNRM (output)
The one-norm of the balanced matrix B.

* RCONDE (output)
If SENSE = 'E' or 'B', the reciprocal condition numbers of the selected eigenvalues, stored in consecutive elements of the array. If SENSE = 'V', RCONDE is not referenced.

* RCONDV (output)
If JOB = 'V' or 'B', the estimated reciprocal condition numbers of the selected eigenvectors, stored in consecutive elements of the array. If the eigenvalues cannot be reordered to compute RCONDV(j), RCONDV(j) is set to 0; this can only occur when the true value would be very small anyway. If SENSE = 'E', RCONDV is not referenced. Not referenced if JOB = 'E'.

* WORK (workspace)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

* LWORK (input)
The dimension of the array WORK. LWORK >= max(1,2*N). If SENSE = 'N' or 'E', LWORK >= 2*N. If SENSE = 'V' or 'B', LWORK >= 2*N*N+2*N.

If LWORK = -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 LWORK is issued by XERBLA.

* RWORK (workspace)
Real workspace.

* IWORK (workspace)
If SENSE = 'E', IWORK is not referenced.

* BWORK (workspace)
If SENSE = 'N', BWORK is not referenced.

* INFO (output)