Contents


NAME

     cgegv - routine is deprecated and has been replaced by  rou-
     tine CGGEV

SYNOPSIS

     SUBROUTINE CGEGV(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHA, BETA, VL,
           LDVL, VR, LDVR, WORK, LDWORK, WORK2, INFO)

     CHARACTER * 1 JOBVL, JOBVR
     COMPLEX A(LDA,*), B(LDB,*), ALPHA(*),  BETA(*),  VL(LDVL,*),
     VR(LDVR,*), WORK(*)
     INTEGER N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
     REAL WORK2(*)

     SUBROUTINE CGEGV_64(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHA, BETA, VL,
           LDVL, VR, LDVR, WORK, LDWORK, WORK2, INFO)

     CHARACTER * 1 JOBVL, JOBVR
     COMPLEX A(LDA,*), B(LDB,*), ALPHA(*),  BETA(*),  VL(LDVL,*),
     VR(LDVR,*), WORK(*)
     INTEGER*8 N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
     REAL WORK2(*)

  F95 INTERFACE
     SUBROUTINE GEGV(JOBVL, JOBVR, [N], A, [LDA], B, [LDB], ALPHA, BETA,
            VL, [LDVL], VR, [LDVR], [WORK], [LDWORK], [WORK2], [INFO])

     CHARACTER(LEN=1) :: JOBVL, JOBVR
     COMPLEX, DIMENSION(:) :: ALPHA, BETA, WORK
     COMPLEX, DIMENSION(:,:) :: A, B, VL, VR
     INTEGER :: N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
     REAL, DIMENSION(:) :: WORK2

     SUBROUTINE GEGV_64(JOBVL, JOBVR, [N], A, [LDA], B, [LDB], ALPHA,
            BETA, VL, [LDVL], VR, [LDVR], [WORK], [LDWORK], [WORK2], [INFO])

     CHARACTER(LEN=1) :: JOBVL, JOBVR
     COMPLEX, DIMENSION(:) :: ALPHA, BETA, WORK
     COMPLEX, DIMENSION(:,:) :: A, B, VL, VR
     INTEGER(8) :: N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
     REAL, DIMENSION(:) :: WORK2

  C INTERFACE
     #include <sunperf.h>

     void cgegv(char jobvl, char jobvr, int n,  complex  *a,  int
               lda,  complex *b, int ldb, complex *alpha, complex
               *beta, complex *vl, int  ldvl,  complex  *vr,  int
               ldvr, int *info);

     void cgegv_64(char jobvl, char jobvr, long  n,  complex  *a,
               long  lda,  complex  *b, long ldb, complex *alpha,
               complex *beta, complex  *vl,  long  ldvl,  complex
               *vr, long ldvr, long *info);

PURPOSE

     cgegv routine is deprecated and has been replaced by routine
     CGGEV.

     CGEGV computes for a pair  of  N-by-N  complex  nonsymmetric
     matrices A and B, the generalized eigenvalues (alpha, beta),
     and optionally, the left and/or right generalized  eigenvec-
     tors (VL and VR).

     A generalized eigenvalue for a pair of  matrices  (A,B)  is,
     roughly  speaking,  a  scalar  w or a ratio  alpha/beta = w,
     such that  A - w*B is singular.  It is  usually  represented
     as the pair (alpha,beta), as there is a reasonable interpre-
     tation for beta=0, and even for both  being  zero.   A  good
     beginning  reference  is the book, "Matrix Computations", by
     G. Golub & C. van Loan (Johns Hopkins U. Press)

     A right generalized eigenvector corresponding to a  general-
     ized eigenvalue  w  for a pair of matrices (A,B) is a vector
     r  such that  (A - w B) r = 0 .  A left  generalized  eigen-
     vector  is  a vector l such that l**H * (A - w B) = 0, where
     l**H is the
     conjugate-transpose of l.

     Note: this routine performs "full balancing"  on  A  and  B.
     See "Further Details", below.

ARGUMENTS

     JOBVL (input)
               = 'N':  do not compute the left generalized eigen-
               vectors;
               = 'V':  compute the left generalized eigenvectors.

     JOBVR (input)
               = 'N':   do  not  compute  the  right  generalized
               eigenvectors;
               = 'V':  compute the  right  generalized  eigenvec-
               tors.

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

     A (input/output)
               On entry, the first of the pair of matrices  whose
               generalized  eigenvalues and (optionally) general-
               ized eigenvectors are to be  computed.   On  exit,
               the  contents  will  have  been destroyed.  (For a
               description of the contents  of  A  on  exit,  see
               "Further Details", below.)

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

     B (input/output)
               On entry, the second of the pair of matrices whose
               generalized  eigenvalues and (optionally) general-
               ized eigenvectors are to be  computed.   On  exit,
               the  contents  will  have  been destroyed.  (For a
               description of the contents  of  B  on  exit,  see
               "Further Details", below.)

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

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

               Note:  the  quotients  ALPHA(j)/VL(j)  may  easily
               over-  or  underflow,  and VL(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  VL  always  less than and usually
               comparable with norm(B).

     VL (output)
               If JOBVL = 'V', the left generalized eigenvectors.
               (See  "Purpose", above.)  Each eigenvector will be
               scaled so the largest component will have abs(real
               part)  +  abs(imag.  part)  = 1, *except* that for
               eigenvalues with alpha=beta=0, a zero vector  will
               be returned as the corresponding eigenvector.  Not
               referenced if JOBVL = 'N'.
     BETA (output)
               If JOBVL = 'V', the left generalized eigenvectors.
               (See  "Purpose", above.)  Each eigenvector will be
               scaled so the largest component will have abs(real
               part)  +  abs(imag.  part)  = 1, *except* that for
               eigenvalues with alpha=beta=0, a zero vector  will
               be returned as the corresponding eigenvector.  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  eigenvec-
               tors.   (See  "Purpose", above.)  Each eigenvector
               will be scaled so the largest component will  have
               abs(real  part)  +  abs(imag.  part) = 1, *except*
               that for eigenvalues  with  alpha=beta=0,  a  zero
               vector  will  be  returned  as  the  corresponding
               eigenvector.  Not referenced if JOBVR = 'N'.

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

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

     LDWORK (input)
               The  dimension  of  the  array  WORK.   LDWORK  >=
               max(1,2*N).   For  good  performance,  LDWORK must
               generally be larger.  To compute the optimal value
               of  LDWORK,  call  ILAENV  to  get blocksizes (for
               CGEQRF, CUNMQR, and CUNGQR.)  Then compute:  NB as
               the  MAX of the blocksizes for CGEQRF, CUNMQR, and
               CUNGQR; The optimal LDWORK is  MAX( 2*N,  N*(NB+1)
               ).

               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.
     WORK2 (workspace)
               dimension(8*N)

     INFO (output)
               = 0:  successful exit
               < 0:  if INFO = -i, the i-th argument had an ille-
               gal value.
               =1,...,N:  The QZ iteration failed.  No  eigenvec-
               tors  have been calculated, but ALPHA(j) and VL(j)
               should  be  correct  for  j=INFO+1,...,N.   >   N:
               errors that usually indicate LAPACK problems:
               =N+1: error return from CGGBAL
               =N+2: error return from CGEQRF
               =N+3: error return from CUNMQR
               =N+4: error return from CUNGQR
               =N+5: error return from CGGHRD
               =N+6: error return from CHGEQZ (other than  failed
               iteration) =N+7: error return from CTGEVC
               =N+8: error return from CGGBAK (computing VL)
               =N+9: error return from CGGBAK (computing VR)
               =N+10: error return from CLASCL (various calls)

FURTHER DETAILS

     Balancing
     ---------

     This driver calls CGGBAL to both permute and scale rows  and
     columns  of  A and B.  The permutations PL and PR are chosen
     so that PL*A*PR and PL*B*R will be upper  triangular  except
     for  the  diagonal  blocks A(i:j,i:j) and B(i:j,i:j), with i
     and j as close together as possible.  The  diagonal  scaling
     matrices   DL   and   DR   are   chosen  so  that  the  pair
     DL*PL*A*PR*DR, DL*PL*B*PR*DR  have  elements  close  to  one
     (except for the elements that start out zero.)

     After the  eigenvalues  and  eigenvectors  of  the  balanced
     matrices have been computed, CGGBAK transforms the eigenvec-
     tors back to what they would have been  (in  perfect  arith-
     metic) if they had not been balanced.

     Contents of A and B on Exit
     -------- -- - --- - -- ----

     If  any  eigenvectors  are  computed  (either  JOBVL='V'  or
     JOBVR='V'  or  both),  then  on exit the arrays A and B will
     contain the complex Schur form[*] of the "balanced" versions
     of  A and B.  If no eigenvectors are computed, then only the
     diagonal blocks will be correct.

     [*] In other words, upper triangular form.