Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

dgegv (3p)

Name

dgegv - routine is deprecated and has been replaced by routine DGGEV

Synopsis

SUBROUTINE DGEGV(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

CHARACTER*1 JOBVL, JOBVR
INTEGER N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
DOUBLE  PRECISION  A(LDA,*),  B(LDB,*),  ALPHAR(*), ALPHAI(*), BETA(*),
VL(LDVL,*), VR(LDVR,*), WORK(*)

SUBROUTINE DGEGV_64(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

CHARACTER*1 JOBVL, JOBVR
INTEGER*8 N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
DOUBLE PRECISION A(LDA,*),  B(LDB,*),  ALPHAR(*),  ALPHAI(*),  BETA(*),
VL(LDVL,*), VR(LDVR,*), WORK(*)




F95 INTERFACE
SUBROUTINE GEGV(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR,
ALPHAI, BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

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

SUBROUTINE GEGV_64(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR,
ALPHAI, BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

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




C INTERFACE
#include <sunperf.h>

void  dgegv(char  jobvl,  char jobvr, int n, double *a, int lda, double
*b, int ldb, double *alphar, double  *alphai,  double  *beta,
double *vl, int ldvl, double *vr, int ldvr, int *info);

void dgegv_64(char jobvl, char jobvr, long n, double *a, long lda, dou-
ble *b, long ldb,  double  *alphar,  double  *alphai,  double
*beta,  double  *vl,  long  ldvl, double *vr, long ldvr, long
*info);

Description

Oracle Solaris Studio Performance Library                            dgegv(3P)



NAME
       dgegv - routine is deprecated and has been replaced by routine DGGEV


SYNOPSIS
       SUBROUTINE DGEGV(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
             BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

       CHARACTER*1 JOBVL, JOBVR
       INTEGER N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
       DOUBLE  PRECISION  A(LDA,*),  B(LDB,*),  ALPHAR(*), ALPHAI(*), BETA(*),
       VL(LDVL,*), VR(LDVR,*), WORK(*)

       SUBROUTINE DGEGV_64(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
             BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

       CHARACTER*1 JOBVL, JOBVR
       INTEGER*8 N, LDA, LDB, LDVL, LDVR, LDWORK, INFO
       DOUBLE PRECISION A(LDA,*),  B(LDB,*),  ALPHAR(*),  ALPHAI(*),  BETA(*),
       VL(LDVL,*), VR(LDVR,*), WORK(*)




   F95 INTERFACE
       SUBROUTINE GEGV(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR,
              ALPHAI, BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

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

       SUBROUTINE GEGV_64(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR,
              ALPHAI, BETA, VL, LDVL, VR, LDVR, WORK, LDWORK, INFO)

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




   C INTERFACE
       #include <sunperf.h>

       void  dgegv(char  jobvl,  char jobvr, int n, double *a, int lda, double
                 *b, int ldb, double *alphar, double  *alphai,  double  *beta,
                 double *vl, int ldvl, double *vr, int ldvr, int *info);

       void dgegv_64(char jobvl, char jobvr, long n, double *a, long lda, dou-
                 ble *b, long ldb,  double  *alphar,  double  *alphai,  double
                 *beta,  double  *vl,  long  ldvl, double *vr, long ldvr, long
                 *info);



PURPOSE
       dgegv routine is deprecated and has been replaced by routine DGGEV.

       DGEGV computes for a pair of n-by-n real nonsymmetric matrices A and B,
       the  generalized  eigenvalues  (alphar +/- alphai*i, beta), and option-
       ally, the left and/or right generalized eigenvectors (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  interpretation  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 generalized eigen-
       value  w  for a pair of matrices (A,B) is a vector  r  such that  (A  -
       w  B)  r  = 0 .  A left generalized eigenvector 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 eigenvectors;
                 = 'V':  compute the left generalized eigenvectors.


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


       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) generalized 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 general-
                 ized eigenvalues and  (optionally)  generalized  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).


       ALPHAR (output)
                 On  exit,  (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
                 be the generalized eigenvalues.  If ALPHAI(j) is  zero,  then
                 the  j-th  eigenvalue is real; if positive, then the j-th and
                 (j+1)-st eigenvalues  are  a  complex  conjugate  pair,  with
                 ALPHAI(j+1) negative.

                 Note:  the  quotients ALPHAR(j)/BETA(j) and ALPHAI(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, ALPHAR and ALPHAI will be  always  less
                 than  and  usually  comparable with norm(A) in magnitude, and
                 BETA always less than and usually comparable with norm(B).


       ALPHAI (output)
                 See the description of ALPHAR.


       BETA (output)
                 See the description of ALPHAR.


       VL (output)
                 If JOBVL = 'V',  the  left  generalized  eigenvectors.   (See
                 "Purpose",  above.)   Real eigenvectors take one column, com-
                 plex take two columns, the first for the real  part  and  the
                 second  for  the imaginary part.  Complex eigenvectors corre-
                 spond to an eigenvalue with positive  imaginary  part.   Each
                 eigenvector will be scaled so the largest component will have
                 abs(real part) + abs(imag. part) = 1, *except* that  for  ei-
                 genvalues  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  eigenvectors.   (See
                 "Purpose",  above.)   Real eigenvectors take one column, com-
                 plex take two columns, the first for the real  part  and  the
                 second  for  the imaginary part.  Complex eigenvectors corre-
                 spond to an eigenvalue with positive  imaginary  part.   Each
                 eigenvector will be scaled so the largest component will have
                 abs(real part) + abs(imag. part) = 1, *except* that  for  ei-
                 genvalues  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,8*N).   For
                 good  performance,  LDWORK must generally be larger.  To com-
                 pute the optimal value of LDWORK, call ILAENV to  get  block-
                 sizes  (for DGEQRF, DORMQR, and DORGQR.)  Then compute: NB --
                 MAX of the blocksizes for DGEQRF,  DORMQR,  and  DORGQR;  The
                 optimal LDWORK is: 2*N + MAX( 6*N, N*(NB+1) ).

                 If  LDWORK  = -1, then a workspace query is assumed; the rou-
                 tine 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.


       INFO (output)
                 = 0:  successful exit
                 < 0:  if INFO = -i, the i-th argument had an illegal value.
                 = 1,...,N: The QZ iteration  failed.   No  eigenvectors  have
                 been calculated, but ALPHAR(j), ALPHAI(j), and BETA(j) should
                 be correct for j=INFO+1,...,N.  >  N:   errors  that  usually
                 indicate LAPACK problems:
                 =N+1: error return from DGGBAL
                 =N+2: error return from DGEQRF
                 =N+3: error return from DORMQR
                 =N+4: error return from DORGQR
                 =N+5: error return from DGGHRD
                 =N+6:  error return from DHGEQZ (other than failed iteration)
                 =N+7: error return from DTGEVC
                 =N+8: error return from DGGBAK (computing VL)
                 =N+9: error return from DGGBAK (computing VR)
                 =N+10: error return from DLASCL (various calls)

FURTHER DETAILS
       Balancing
       ---------

       This driver calls DGGBAL 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,  DGGBAK  transforms  the eigenvectors back to what they
       would have been (in perfect arithmetic) 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  real  Schur
       form[*]  of the "balanced" versions of A and B.  If no eigenvectors are
       computed, then only the diagonal blocks will be correct.

       [*] See DHGEQZ, DGEGS, or read the book "Matrix Computations",
           by Golub & van Loan, pub. by Johns Hopkins U. Press.




                                  7 Nov 2015                         dgegv(3P)