Contents


NAME

     dsygv - compute all the  eigenvalues,  and  optionally,  the
     eigenvectors of a real generalized symmetric-definite eigen-
     problem, of the form A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or
     B*A*x=(lambda)*x

SYNOPSIS

     SUBROUTINE DSYGV(ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
           LDWORK, INFO)

     CHARACTER * 1 JOBZ, UPLO
     INTEGER ITYPE, N, LDA, LDB, LDWORK, INFO
     DOUBLE PRECISION A(LDA,*), B(LDB,*), W(*), WORK(*)

     SUBROUTINE DSYGV_64(ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
           LDWORK, INFO)

     CHARACTER * 1 JOBZ, UPLO
     INTEGER*8 ITYPE, N, LDA, LDB, LDWORK, INFO
     DOUBLE PRECISION A(LDA,*), B(LDB,*), W(*), WORK(*)

  F95 INTERFACE
     SUBROUTINE SYGV(ITYPE, JOBZ, UPLO, N, A, [LDA], B, [LDB], W, [WORK],
            [LDWORK], [INFO])

     CHARACTER(LEN=1) :: JOBZ, UPLO
     INTEGER :: ITYPE, N, LDA, LDB, LDWORK, INFO
     REAL(8), DIMENSION(:) :: W, WORK
     REAL(8), DIMENSION(:,:) :: A, B

     SUBROUTINE SYGV_64(ITYPE, JOBZ, UPLO, N, A, [LDA], B, [LDB], W, [WORK],
            [LDWORK], [INFO])

     CHARACTER(LEN=1) :: JOBZ, UPLO
     INTEGER(8) :: ITYPE, N, LDA, LDB, LDWORK, INFO
     REAL(8), DIMENSION(:) :: W, WORK
     REAL(8), DIMENSION(:,:) :: A, B

  C INTERFACE
     #include <sunperf.h>

     void dsygv(int itype, char jobz, char uplo,  int  n,  double
               *a,  int  lda,  double *b, int ldb, double *w, int
               *info);

     void dsygv_64(long itype, char jobz, char uplo, long n, dou-
               ble  *a, long lda, double *b, long ldb, double *w,
               long *info);

PURPOSE

     dsygv computes all  the  eigenvalues,  and  optionally,  the
     eigenvectors of a real generalized symmetric-definite eigen-
     problem, of the form A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or
     B*A*x=(lambda)*x.   Here A and B are assumed to be symmetric
     and B is also
     positive definite.

ARGUMENTS

     ITYPE (input)
               Specifies the problem type to be solved:
               = 1:  A*x = (lambda)*B*x
               = 2:  A*B*x = (lambda)*x
               = 3:  B*A*x = (lambda)*x

     JOBZ (input)
               = 'N':  Compute eigenvalues only;
               = 'V':  Compute eigenvalues and eigenvectors.

     UPLO (input)
               = 'U':  Upper triangles of A and B are stored;
               = 'L':  Lower triangles of A and B are stored.

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

     A (input/output)
               On entry, the symmetric matrix A.  If UPLO =  'U',
               the leading N-by-N upper triangular part of A con-
               tains the upper triangular part of the  matrix  A.
               If UPLO = 'L', the leading N-by-N lower triangular
               part of A contains the lower  triangular  part  of
               the matrix A.

               On exit, if JOBZ = 'V', then if INFO = 0,  A  con-
               tains the matrix Z of eigenvectors.  The eigenvec-
               tors are normalized as follows:  if ITYPE =  1  or
               2,  Z**T*B*Z = I; if ITYPE = 3, Z**T*inv(B)*Z = I.
               If JOBZ = 'N', then on exit the upper triangle (if
               UPLO='U')  or  the lower triangle (if UPLO='L') of
               A, including the diagonal, is destroyed.

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

     B (input/output)
               On entry, the symmetric positive  definite  matrix
               B.   If  UPLO = 'U', the leading N-by-N upper tri-
               angular part of B contains  the  upper  triangular
               part  of the matrix B.  If UPLO = 'L', the leading
               N-by-N lower triangular part  of  B  contains  the
               lower triangular part of the matrix B.

               On exit, if INFO <= N, the part  of  B  containing
               the matrix is overwritten by the triangular factor
               U or L from the Cholesky factorization B =  U**T*U
               or B = L*L**T.

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

     W (output)
               If INFO = 0, the eigenvalues in ascending order.

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

     LDWORK (input)
               The  length  of  the  array   WORK.    LDWORK   >=
               max(1,3*N-1).   For  optimal efficiency, LDWORK >=
               (NB+2)*N, where NB is  the  blocksize  for  SSYTRD
               returned by ILAENV.

               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.

     INFO (output)
               = 0:  successful exit
               < 0:  if INFO = -i, the i-th argument had an ille-
               gal value
               > 0:  SPOTRF or SSYEV returned an error code:
               <= N:  if INFO = i, SSYEV failed  to  converge;  i
               off-diagonal    elements    of   an   intermediate
               tridiagonal form did not converge to  zero;  >  N:
               if INFO = N + i, for 1 <= i <= N, then the leading
               minor of order i of B is  not  positive  definite.
               The  factorization of B could not be completed and
               no eigenvalues or eigenvectors were computed.