001:       SUBROUTINE DSYGV( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W, WORK,
002:      $                  LWORK, INFO )
003: *
004: *  -- LAPACK driver routine (version 3.2) --
005: *     Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
006: *     November 2006
007: *
008: *     .. Scalar Arguments ..
009:       CHARACTER          JOBZ, UPLO
010:       INTEGER            INFO, ITYPE, LDA, LDB, LWORK, N
011: *     ..
012: *     .. Array Arguments ..
013:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
014: *     ..
015: *
016: *  Purpose
017: *  =======
018: *
019: *  DSYGV computes all the eigenvalues, and optionally, the eigenvectors
020: *  of a real generalized symmetric-definite eigenproblem, of the form
021: *  A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
022: *  Here A and B are assumed to be symmetric and B is also
023: *  positive definite.
024: *
025: *  Arguments
026: *  =========
027: *
028: *  ITYPE   (input) INTEGER
029: *          Specifies the problem type to be solved:
030: *          = 1:  A*x = (lambda)*B*x
031: *          = 2:  A*B*x = (lambda)*x
032: *          = 3:  B*A*x = (lambda)*x
033: *
034: *  JOBZ    (input) CHARACTER*1
035: *          = 'N':  Compute eigenvalues only;
036: *          = 'V':  Compute eigenvalues and eigenvectors.
037: *
038: *  UPLO    (input) CHARACTER*1
039: *          = 'U':  Upper triangles of A and B are stored;
040: *          = 'L':  Lower triangles of A and B are stored.
041: *
042: *  N       (input) INTEGER
043: *          The order of the matrices A and B.  N >= 0.
044: *
045: *  A       (input/output) DOUBLE PRECISION array, dimension (LDA, N)
046: *          On entry, the symmetric matrix A.  If UPLO = 'U', the
047: *          leading N-by-N upper triangular part of A contains the
048: *          upper triangular part of the matrix A.  If UPLO = 'L',
049: *          the leading N-by-N lower triangular part of A contains
050: *          the lower triangular part of the matrix A.
051: *
052: *          On exit, if JOBZ = 'V', then if INFO = 0, A contains the
053: *          matrix Z of eigenvectors.  The eigenvectors are normalized
054: *          as follows:
055: *          if ITYPE = 1 or 2, Z**T*B*Z = I;
056: *          if ITYPE = 3, Z**T*inv(B)*Z = I.
057: *          If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
058: *          or the lower triangle (if UPLO='L') of A, including the
059: *          diagonal, is destroyed.
060: *
061: *  LDA     (input) INTEGER
062: *          The leading dimension of the array A.  LDA >= max(1,N).
063: *
064: *  B       (input/output) DOUBLE PRECISION array, dimension (LDB, N)
065: *          On entry, the symmetric positive definite matrix B.
066: *          If UPLO = 'U', the leading N-by-N upper triangular part of B
067: *          contains the upper triangular part of the matrix B.
068: *          If UPLO = 'L', the leading N-by-N lower triangular part of B
069: *          contains the lower triangular part of the matrix B.
070: *
071: *          On exit, if INFO <= N, the part of B containing the matrix is
072: *          overwritten by the triangular factor U or L from the Cholesky
073: *          factorization B = U**T*U or B = L*L**T.
074: *
075: *  LDB     (input) INTEGER
076: *          The leading dimension of the array B.  LDB >= max(1,N).
077: *
078: *  W       (output) DOUBLE PRECISION array, dimension (N)
079: *          If INFO = 0, the eigenvalues in ascending order.
080: *
081: *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK))
082: *          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
083: *
084: *  LWORK   (input) INTEGER
085: *          The length of the array WORK.  LWORK >= max(1,3*N-1).
086: *          For optimal efficiency, LWORK >= (NB+2)*N,
087: *          where NB is the blocksize for DSYTRD returned by ILAENV.
088: *
089: *          If LWORK = -1, then a workspace query is assumed; the routine
090: *          only calculates the optimal size of the WORK array, returns
091: *          this value as the first entry of the WORK array, and no error
092: *          message related to LWORK is issued by XERBLA.
093: *
094: *  INFO    (output) INTEGER
095: *          = 0:  successful exit
096: *          < 0:  if INFO = -i, the i-th argument had an illegal value
097: *          > 0:  DPOTRF or DSYEV returned an error code:
098: *             <= N:  if INFO = i, DSYEV failed to converge;
099: *                    i off-diagonal elements of an intermediate
100: *                    tridiagonal form did not converge to zero;
101: *             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
102: *                    minor of order i of B is not positive definite.
103: *                    The factorization of B could not be completed and
104: *                    no eigenvalues or eigenvectors were computed.
105: *
106: *  =====================================================================
107: *
108: *     .. Parameters ..
109:       DOUBLE PRECISION   ONE
110:       PARAMETER          ( ONE = 1.0D+0 )
111: *     ..
112: *     .. Local Scalars ..
113:       LOGICAL            LQUERY, UPPER, WANTZ
114:       CHARACTER          TRANS
115:       INTEGER            LWKMIN, LWKOPT, NB, NEIG
116: *     ..
117: *     .. External Functions ..
118:       LOGICAL            LSAME
119:       INTEGER            ILAENV
120:       EXTERNAL           LSAME, ILAENV
121: *     ..
122: *     .. External Subroutines ..
123:       EXTERNAL           DPOTRF, DSYEV, DSYGST, DTRMM, DTRSM, XERBLA
124: *     ..
125: *     .. Intrinsic Functions ..
126:       INTRINSIC          MAX
127: *     ..
128: *     .. Executable Statements ..
129: *
130: *     Test the input parameters.
131: *
132:       WANTZ = LSAME( JOBZ, 'V' )
133:       UPPER = LSAME( UPLO, 'U' )
134:       LQUERY = ( LWORK.EQ.-1 )
135: *
136:       INFO = 0
137:       IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
138:          INFO = -1
139:       ELSE IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
140:          INFO = -2
141:       ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
142:          INFO = -3
143:       ELSE IF( N.LT.0 ) THEN
144:          INFO = -4
145:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
146:          INFO = -6
147:       ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
148:          INFO = -8
149:       END IF
150: *
151:       IF( INFO.EQ.0 ) THEN
152:          LWKMIN = MAX( 1, 3*N - 1 )
153:          NB = ILAENV( 1, 'DSYTRD', UPLO, N, -1, -1, -1 )
154:          LWKOPT = MAX( LWKMIN, ( NB + 2 )*N )
155:          WORK( 1 ) = LWKOPT
156: *
157:          IF( LWORK.LT.LWKMIN .AND. .NOT.LQUERY ) THEN
158:             INFO = -11
159:          END IF
160:       END IF
161: *
162:       IF( INFO.NE.0 ) THEN
163:          CALL XERBLA( 'DSYGV ', -INFO )
164:          RETURN
165:       ELSE IF( LQUERY ) THEN
166:          RETURN
167:       END IF
168: *
169: *     Quick return if possible
170: *
171:       IF( N.EQ.0 )
172:      $   RETURN
173: *
174: *     Form a Cholesky factorization of B.
175: *
176:       CALL DPOTRF( UPLO, N, B, LDB, INFO )
177:       IF( INFO.NE.0 ) THEN
178:          INFO = N + INFO
179:          RETURN
180:       END IF
181: *
182: *     Transform problem to standard eigenvalue problem and solve.
183: *
184:       CALL DSYGST( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
185:       CALL DSYEV( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO )
186: *
187:       IF( WANTZ ) THEN
188: *
189: *        Backtransform eigenvectors to the original problem.
190: *
191:          NEIG = N
192:          IF( INFO.GT.0 )
193:      $      NEIG = INFO - 1
194:          IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
195: *
196: *           For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
197: *           backtransform eigenvectors: x = inv(L)'*y or inv(U)*y
198: *
199:             IF( UPPER ) THEN
200:                TRANS = 'N'
201:             ELSE
202:                TRANS = 'T'
203:             END IF
204: *
205:             CALL DTRSM( 'Left', UPLO, TRANS, 'Non-unit', N, NEIG, ONE,
206:      $                  B, LDB, A, LDA )
207: *
208:          ELSE IF( ITYPE.EQ.3 ) THEN
209: *
210: *           For B*A*x=(lambda)*x;
211: *           backtransform eigenvectors: x = L*y or U'*y
212: *
213:             IF( UPPER ) THEN
214:                TRANS = 'T'
215:             ELSE
216:                TRANS = 'N'
217:             END IF
218: *
219:             CALL DTRMM( 'Left', UPLO, TRANS, 'Non-unit', N, NEIG, ONE,
220:      $                  B, LDB, A, LDA )
221:          END IF
222:       END IF
223: *
224:       WORK( 1 ) = LWKOPT
225:       RETURN
226: *
227: *     End of DSYGV
228: *
229:       END
230: