001:       SUBROUTINE DSYGVX( ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB,
002:      $                   VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
003:      $                   LWORK, IWORK, IFAIL, INFO )
004: *
005: *  -- LAPACK driver routine (version 3.2) --
006: *     Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
007: *     November 2006
008: *
009: *     .. Scalar Arguments ..
010:       CHARACTER          JOBZ, RANGE, UPLO
011:       INTEGER            IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N
012:       DOUBLE PRECISION   ABSTOL, VL, VU
013: *     ..
014: *     .. Array Arguments ..
015:       INTEGER            IFAIL( * ), IWORK( * )
016:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), W( * ), WORK( * ),
017:      $                   Z( LDZ, * )
018: *     ..
019: *
020: *  Purpose
021: *  =======
022: *
023: *  DSYGVX computes selected eigenvalues, and optionally, eigenvectors
024: *  of a real generalized symmetric-definite eigenproblem, of the form
025: *  A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.  Here A
026: *  and B are assumed to be symmetric and B is also positive definite.
027: *  Eigenvalues and eigenvectors can be selected by specifying either a
028: *  range of values or a range of indices for the desired eigenvalues.
029: *
030: *  Arguments
031: *  =========
032: *
033: *  ITYPE   (input) INTEGER
034: *          Specifies the problem type to be solved:
035: *          = 1:  A*x = (lambda)*B*x
036: *          = 2:  A*B*x = (lambda)*x
037: *          = 3:  B*A*x = (lambda)*x
038: *
039: *  JOBZ    (input) CHARACTER*1
040: *          = 'N':  Compute eigenvalues only;
041: *          = 'V':  Compute eigenvalues and eigenvectors.
042: *
043: *  RANGE   (input) CHARACTER*1
044: *          = 'A': all eigenvalues will be found.
045: *          = 'V': all eigenvalues in the half-open interval (VL,VU]
046: *                 will be found.
047: *          = 'I': the IL-th through IU-th eigenvalues will be found.
048: *
049: *  UPLO    (input) CHARACTER*1
050: *          = 'U':  Upper triangle of A and B are stored;
051: *          = 'L':  Lower triangle of A and B are stored.
052: *
053: *  N       (input) INTEGER
054: *          The order of the matrix pencil (A,B).  N >= 0.
055: *
056: *  A       (input/output) DOUBLE PRECISION array, dimension (LDA, N)
057: *          On entry, the symmetric matrix A.  If UPLO = 'U', the
058: *          leading N-by-N upper triangular part of A contains the
059: *          upper triangular part of the matrix A.  If UPLO = 'L',
060: *          the leading N-by-N lower triangular part of A contains
061: *          the lower triangular part of the matrix A.
062: *
063: *          On exit, the lower triangle (if UPLO='L') or the upper
064: *          triangle (if UPLO='U') of A, including the diagonal, is
065: *          destroyed.
066: *
067: *  LDA     (input) INTEGER
068: *          The leading dimension of the array A.  LDA >= max(1,N).
069: *
070: *  B       (input/output) DOUBLE PRECISION array, dimension (LDA, N)
071: *          On entry, the symmetric matrix B.  If UPLO = 'U', the
072: *          leading N-by-N upper triangular part of B contains the
073: *          upper triangular part of the matrix B.  If UPLO = 'L',
074: *          the leading N-by-N lower triangular part of B contains
075: *          the lower triangular part of the matrix B.
076: *
077: *          On exit, if INFO <= N, the part of B containing the matrix is
078: *          overwritten by the triangular factor U or L from the Cholesky
079: *          factorization B = U**T*U or B = L*L**T.
080: *
081: *  LDB     (input) INTEGER
082: *          The leading dimension of the array B.  LDB >= max(1,N).
083: *
084: *  VL      (input) DOUBLE PRECISION
085: *  VU      (input) DOUBLE PRECISION
086: *          If RANGE='V', the lower and upper bounds of the interval to
087: *          be searched for eigenvalues. VL < VU.
088: *          Not referenced if RANGE = 'A' or 'I'.
089: *
090: *  IL      (input) INTEGER
091: *  IU      (input) INTEGER
092: *          If RANGE='I', the indices (in ascending order) of the
093: *          smallest and largest eigenvalues to be returned.
094: *          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
095: *          Not referenced if RANGE = 'A' or 'V'.
096: *
097: *  ABSTOL  (input) DOUBLE PRECISION
098: *          The absolute error tolerance for the eigenvalues.
099: *          An approximate eigenvalue is accepted as converged
100: *          when it is determined to lie in an interval [a,b]
101: *          of width less than or equal to
102: *
103: *                  ABSTOL + EPS *   max( |a|,|b| ) ,
104: *
105: *          where EPS is the machine precision.  If ABSTOL is less than
106: *          or equal to zero, then  EPS*|T|  will be used in its place,
107: *          where |T| is the 1-norm of the tridiagonal matrix obtained
108: *          by reducing A to tridiagonal form.
109: *
110: *          Eigenvalues will be computed most accurately when ABSTOL is
111: *          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
112: *          If this routine returns with INFO>0, indicating that some
113: *          eigenvectors did not converge, try setting ABSTOL to
114: *          2*DLAMCH('S').
115: *
116: *  M       (output) INTEGER
117: *          The total number of eigenvalues found.  0 <= M <= N.
118: *          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
119: *
120: *  W       (output) DOUBLE PRECISION array, dimension (N)
121: *          On normal exit, the first M elements contain the selected
122: *          eigenvalues in ascending order.
123: *
124: *  Z       (output) DOUBLE PRECISION array, dimension (LDZ, max(1,M))
125: *          If JOBZ = 'N', then Z is not referenced.
126: *          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
127: *          contain the orthonormal eigenvectors of the matrix A
128: *          corresponding to the selected eigenvalues, with the i-th
129: *          column of Z holding the eigenvector associated with W(i).
130: *          The eigenvectors are normalized as follows:
131: *          if ITYPE = 1 or 2, Z**T*B*Z = I;
132: *          if ITYPE = 3, Z**T*inv(B)*Z = I.
133: *
134: *          If an eigenvector fails to converge, then that column of Z
135: *          contains the latest approximation to the eigenvector, and the
136: *          index of the eigenvector is returned in IFAIL.
137: *          Note: the user must ensure that at least max(1,M) columns are
138: *          supplied in the array Z; if RANGE = 'V', the exact value of M
139: *          is not known in advance and an upper bound must be used.
140: *
141: *  LDZ     (input) INTEGER
142: *          The leading dimension of the array Z.  LDZ >= 1, and if
143: *          JOBZ = 'V', LDZ >= max(1,N).
144: *
145: *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK))
146: *          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
147: *
148: *  LWORK   (input) INTEGER
149: *          The length of the array WORK.  LWORK >= max(1,8*N).
150: *          For optimal efficiency, LWORK >= (NB+3)*N,
151: *          where NB is the blocksize for DSYTRD returned by ILAENV.
152: *
153: *          If LWORK = -1, then a workspace query is assumed; the routine
154: *          only calculates the optimal size of the WORK array, returns
155: *          this value as the first entry of the WORK array, and no error
156: *          message related to LWORK is issued by XERBLA.
157: *
158: *  IWORK   (workspace) INTEGER array, dimension (5*N)
159: *
160: *  IFAIL   (output) INTEGER array, dimension (N)
161: *          If JOBZ = 'V', then if INFO = 0, the first M elements of
162: *          IFAIL are zero.  If INFO > 0, then IFAIL contains the
163: *          indices of the eigenvectors that failed to converge.
164: *          If JOBZ = 'N', then IFAIL is not referenced.
165: *
166: *  INFO    (output) INTEGER
167: *          = 0:  successful exit
168: *          < 0:  if INFO = -i, the i-th argument had an illegal value
169: *          > 0:  DPOTRF or DSYEVX returned an error code:
170: *             <= N:  if INFO = i, DSYEVX failed to converge;
171: *                    i eigenvectors failed to converge.  Their indices
172: *                    are stored in array IFAIL.
173: *             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
174: *                    minor of order i of B is not positive definite.
175: *                    The factorization of B could not be completed and
176: *                    no eigenvalues or eigenvectors were computed.
177: *
178: *  Further Details
179: *  ===============
180: *
181: *  Based on contributions by
182: *     Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
183: *
184: * =====================================================================
185: *
186: *     .. Parameters ..
187:       DOUBLE PRECISION   ONE
188:       PARAMETER          ( ONE = 1.0D+0 )
189: *     ..
190: *     .. Local Scalars ..
191:       LOGICAL            ALLEIG, INDEIG, LQUERY, UPPER, VALEIG, WANTZ
192:       CHARACTER          TRANS
193:       INTEGER            LWKMIN, LWKOPT, NB
194: *     ..
195: *     .. External Functions ..
196:       LOGICAL            LSAME
197:       INTEGER            ILAENV
198:       EXTERNAL           LSAME, ILAENV
199: *     ..
200: *     .. External Subroutines ..
201:       EXTERNAL           DPOTRF, DSYEVX, DSYGST, DTRMM, DTRSM, XERBLA
202: *     ..
203: *     .. Intrinsic Functions ..
204:       INTRINSIC          MAX, MIN
205: *     ..
206: *     .. Executable Statements ..
207: *
208: *     Test the input parameters.
209: *
210:       UPPER = LSAME( UPLO, 'U' )
211:       WANTZ = LSAME( JOBZ, 'V' )
212:       ALLEIG = LSAME( RANGE, 'A' )
213:       VALEIG = LSAME( RANGE, 'V' )
214:       INDEIG = LSAME( RANGE, 'I' )
215:       LQUERY = ( LWORK.EQ.-1 )
216: *
217:       INFO = 0
218:       IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
219:          INFO = -1
220:       ELSE IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
221:          INFO = -2
222:       ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
223:          INFO = -3
224:       ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
225:          INFO = -4
226:       ELSE IF( N.LT.0 ) THEN
227:          INFO = -5
228:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
229:          INFO = -7
230:       ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
231:          INFO = -9
232:       ELSE
233:          IF( VALEIG ) THEN
234:             IF( N.GT.0 .AND. VU.LE.VL )
235:      $         INFO = -11
236:          ELSE IF( INDEIG ) THEN
237:             IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
238:                INFO = -12
239:             ELSE IF( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
240:                INFO = -13
241:             END IF
242:          END IF
243:       END IF
244:       IF (INFO.EQ.0) THEN
245:          IF (LDZ.LT.1 .OR. (WANTZ .AND. LDZ.LT.N)) THEN
246:             INFO = -18
247:          END IF
248:       END IF
249: *
250:       IF( INFO.EQ.0 ) THEN
251:          LWKMIN = MAX( 1, 8*N )
252:          NB = ILAENV( 1, 'DSYTRD', UPLO, N, -1, -1, -1 )
253:          LWKOPT = MAX( LWKMIN, ( NB + 3 )*N )
254:          WORK( 1 ) = LWKOPT
255: *
256:          IF( LWORK.LT.LWKMIN .AND. .NOT.LQUERY ) THEN
257:             INFO = -20
258:          END IF
259:       END IF
260: *
261:       IF( INFO.NE.0 ) THEN
262:          CALL XERBLA( 'DSYGVX', -INFO )
263:          RETURN
264:       ELSE IF( LQUERY ) THEN
265:          RETURN
266:       END IF
267: *
268: *     Quick return if possible
269: *
270:       M = 0
271:       IF( N.EQ.0 ) THEN
272:          RETURN
273:       END IF
274: *
275: *     Form a Cholesky factorization of B.
276: *
277:       CALL DPOTRF( UPLO, N, B, LDB, INFO )
278:       IF( INFO.NE.0 ) THEN
279:          INFO = N + INFO
280:          RETURN
281:       END IF
282: *
283: *     Transform problem to standard eigenvalue problem and solve.
284: *
285:       CALL DSYGST( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
286:       CALL DSYEVX( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL,
287:      $             M, W, Z, LDZ, WORK, LWORK, IWORK, IFAIL, INFO )
288: *
289:       IF( WANTZ ) THEN
290: *
291: *        Backtransform eigenvectors to the original problem.
292: *
293:          IF( INFO.GT.0 )
294:      $      M = INFO - 1
295:          IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
296: *
297: *           For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
298: *           backtransform eigenvectors: x = inv(L)'*y or inv(U)*y
299: *
300:             IF( UPPER ) THEN
301:                TRANS = 'N'
302:             ELSE
303:                TRANS = 'T'
304:             END IF
305: *
306:             CALL DTRSM( 'Left', UPLO, TRANS, 'Non-unit', N, M, ONE, B,
307:      $                  LDB, Z, LDZ )
308: *
309:          ELSE IF( ITYPE.EQ.3 ) THEN
310: *
311: *           For B*A*x=(lambda)*x;
312: *           backtransform eigenvectors: x = L*y or U'*y
313: *
314:             IF( UPPER ) THEN
315:                TRANS = 'T'
316:             ELSE
317:                TRANS = 'N'
318:             END IF
319: *
320:             CALL DTRMM( 'Left', UPLO, TRANS, 'Non-unit', N, M, ONE, B,
321:      $                  LDB, Z, LDZ )
322:          END IF
323:       END IF
324: *
325: *     Set WORK(1) to optimal workspace size.
326: *
327:       WORK( 1 ) = LWKOPT
328: *
329:       RETURN
330: *
331: *     End of DSYGVX
332: *
333:       END
334: