1: SUBROUTINE ZHEGVX( ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB,
2: $ VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
3: $ LWORK, RWORK, IWORK, IFAIL, INFO )
4: *
5: * -- LAPACK driver routine (version 3.3.1) --
6: * -- LAPACK is a software package provided by Univ. of Tennessee, --
7: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
8: * -- April 2011 --
9: *
10: * .. Scalar Arguments ..
11: CHARACTER JOBZ, RANGE, UPLO
12: INTEGER IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N
13: DOUBLE PRECISION ABSTOL, VL, VU
14: * ..
15: * .. Array Arguments ..
16: INTEGER IFAIL( * ), IWORK( * )
17: DOUBLE PRECISION RWORK( * ), W( * )
18: COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * ),
19: $ Z( LDZ, * )
20: * ..
21: *
22: * Purpose
23: * =======
24: *
25: * ZHEGVX computes selected eigenvalues, and optionally, eigenvectors
26: * of a complex generalized Hermitian-definite eigenproblem, of the form
27: * A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. Here A and
28: * B are assumed to be Hermitian and B is also positive definite.
29: * Eigenvalues and eigenvectors can be selected by specifying either a
30: * range of values or a range of indices for the desired eigenvalues.
31: *
32: * Arguments
33: * =========
34: *
35: * ITYPE (input) INTEGER
36: * Specifies the problem type to be solved:
37: * = 1: A*x = (lambda)*B*x
38: * = 2: A*B*x = (lambda)*x
39: * = 3: B*A*x = (lambda)*x
40: *
41: * JOBZ (input) CHARACTER*1
42: * = 'N': Compute eigenvalues only;
43: * = 'V': Compute eigenvalues and eigenvectors.
44: *
45: * RANGE (input) CHARACTER*1
46: * = 'A': all eigenvalues will be found.
47: * = 'V': all eigenvalues in the half-open interval (VL,VU]
48: * will be found.
49: * = 'I': the IL-th through IU-th eigenvalues will be found.
50: **
51: * UPLO (input) CHARACTER*1
52: * = 'U': Upper triangles of A and B are stored;
53: * = 'L': Lower triangles of A and B are stored.
54: *
55: * N (input) INTEGER
56: * The order of the matrices A and B. N >= 0.
57: *
58: * A (input/output) COMPLEX*16 array, dimension (LDA, N)
59: * On entry, the Hermitian matrix A. If UPLO = 'U', the
60: * leading N-by-N upper triangular part of A contains the
61: * upper triangular part of the matrix A. If UPLO = 'L',
62: * the leading N-by-N lower triangular part of A contains
63: * the lower triangular part of the matrix A.
64: *
65: * On exit, the lower triangle (if UPLO='L') or the upper
66: * triangle (if UPLO='U') of A, including the diagonal, is
67: * destroyed.
68: *
69: * LDA (input) INTEGER
70: * The leading dimension of the array A. LDA >= max(1,N).
71: *
72: * B (input/output) COMPLEX*16 array, dimension (LDB, N)
73: * On entry, the Hermitian matrix B. If UPLO = 'U', the
74: * leading N-by-N upper triangular part of B contains the
75: * upper triangular part of the matrix B. If UPLO = 'L',
76: * the leading N-by-N lower triangular part of B contains
77: * the lower triangular part of the matrix B.
78: *
79: * On exit, if INFO <= N, the part of B containing the matrix is
80: * overwritten by the triangular factor U or L from the Cholesky
81: * factorization B = U**H*U or B = L*L**H.
82: *
83: * LDB (input) INTEGER
84: * The leading dimension of the array B. LDB >= max(1,N).
85: *
86: * VL (input) DOUBLE PRECISION
87: * VU (input) DOUBLE PRECISION
88: * If RANGE='V', the lower and upper bounds of the interval to
89: * be searched for eigenvalues. VL < VU.
90: * Not referenced if RANGE = 'A' or 'I'.
91: *
92: * IL (input) INTEGER
93: * IU (input) INTEGER
94: * If RANGE='I', the indices (in ascending order) of the
95: * smallest and largest eigenvalues to be returned.
96: * 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
97: * Not referenced if RANGE = 'A' or 'V'.
98: *
99: * ABSTOL (input) DOUBLE PRECISION
100: * The absolute error tolerance for the eigenvalues.
101: * An approximate eigenvalue is accepted as converged
102: * when it is determined to lie in an interval [a,b]
103: * of width less than or equal to
104: *
105: * ABSTOL + EPS * max( |a|,|b| ) ,
106: *
107: * where EPS is the machine precision. If ABSTOL is less than
108: * or equal to zero, then EPS*|T| will be used in its place,
109: * where |T| is the 1-norm of the tridiagonal matrix obtained
110: * by reducing A to tridiagonal form.
111: *
112: * Eigenvalues will be computed most accurately when ABSTOL is
113: * set to twice the underflow threshold 2*DLAMCH('S'), not zero.
114: * If this routine returns with INFO>0, indicating that some
115: * eigenvectors did not converge, try setting ABSTOL to
116: * 2*DLAMCH('S').
117: *
118: * M (output) INTEGER
119: * The total number of eigenvalues found. 0 <= M <= N.
120: * If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
121: *
122: * W (output) DOUBLE PRECISION array, dimension (N)
123: * The first M elements contain the selected
124: * eigenvalues in ascending order.
125: *
126: * Z (output) COMPLEX*16 array, dimension (LDZ, max(1,M))
127: * If JOBZ = 'N', then Z is not referenced.
128: * If JOBZ = 'V', then if INFO = 0, the first M columns of Z
129: * contain the orthonormal eigenvectors of the matrix A
130: * corresponding to the selected eigenvalues, with the i-th
131: * column of Z holding the eigenvector associated with W(i).
132: * The eigenvectors are normalized as follows:
133: * if ITYPE = 1 or 2, Z**T*B*Z = I;
134: * if ITYPE = 3, Z**T*inv(B)*Z = I.
135: *
136: * If an eigenvector fails to converge, then that column of Z
137: * contains the latest approximation to the eigenvector, and the
138: * index of the eigenvector is returned in IFAIL.
139: * Note: the user must ensure that at least max(1,M) columns are
140: * supplied in the array Z; if RANGE = 'V', the exact value of M
141: * is not known in advance and an upper bound must be used.
142: *
143: * LDZ (input) INTEGER
144: * The leading dimension of the array Z. LDZ >= 1, and if
145: * JOBZ = 'V', LDZ >= max(1,N).
146: *
147: * WORK (workspace/output) COMPLEX*16 array, dimension (MAX(1,LWORK))
148: * On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
149: *
150: * LWORK (input) INTEGER
151: * The length of the array WORK. LWORK >= max(1,2*N).
152: * For optimal efficiency, LWORK >= (NB+1)*N,
153: * where NB is the blocksize for ZHETRD returned by ILAENV.
154: *
155: * If LWORK = -1, then a workspace query is assumed; the routine
156: * only calculates the optimal size of the WORK array, returns
157: * this value as the first entry of the WORK array, and no error
158: * message related to LWORK is issued by XERBLA.
159: *
160: * RWORK (workspace) DOUBLE PRECISION array, dimension (7*N)
161: *
162: * IWORK (workspace) INTEGER array, dimension (5*N)
163: *
164: * IFAIL (output) INTEGER array, dimension (N)
165: * If JOBZ = 'V', then if INFO = 0, the first M elements of
166: * IFAIL are zero. If INFO > 0, then IFAIL contains the
167: * indices of the eigenvectors that failed to converge.
168: * If JOBZ = 'N', then IFAIL is not referenced.
169: *
170: * INFO (output) INTEGER
171: * = 0: successful exit
172: * < 0: if INFO = -i, the i-th argument had an illegal value
173: * > 0: ZPOTRF or ZHEEVX returned an error code:
174: * <= N: if INFO = i, ZHEEVX failed to converge;
175: * i eigenvectors failed to converge. Their indices
176: * are stored in array IFAIL.
177: * > N: if INFO = N + i, for 1 <= i <= N, then the leading
178: * minor of order i of B is not positive definite.
179: * The factorization of B could not be completed and
180: * no eigenvalues or eigenvectors were computed.
181: *
182: * Further Details
183: * ===============
184: *
185: * Based on contributions by
186: * Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
187: *
188: * =====================================================================
189: *
190: * .. Parameters ..
191: COMPLEX*16 CONE
192: PARAMETER ( CONE = ( 1.0D+0, 0.0D+0 ) )
193: * ..
194: * .. Local Scalars ..
195: LOGICAL ALLEIG, INDEIG, LQUERY, UPPER, VALEIG, WANTZ
196: CHARACTER TRANS
197: INTEGER LWKOPT, NB
198: * ..
199: * .. External Functions ..
200: LOGICAL LSAME
201: INTEGER ILAENV
202: EXTERNAL LSAME, ILAENV
203: * ..
204: * .. External Subroutines ..
205: EXTERNAL XERBLA, ZHEEVX, ZHEGST, ZPOTRF, ZTRMM, ZTRSM
206: * ..
207: * .. Intrinsic Functions ..
208: INTRINSIC MAX, MIN
209: * ..
210: * .. Executable Statements ..
211: *
212: * Test the input parameters.
213: *
214: WANTZ = LSAME( JOBZ, 'V' )
215: UPPER = LSAME( UPLO, 'U' )
216: ALLEIG = LSAME( RANGE, 'A' )
217: VALEIG = LSAME( RANGE, 'V' )
218: INDEIG = LSAME( RANGE, 'I' )
219: LQUERY = ( LWORK.EQ.-1 )
220: *
221: INFO = 0
222: IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
223: INFO = -1
224: ELSE IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
225: INFO = -2
226: ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
227: INFO = -3
228: ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
229: INFO = -4
230: ELSE IF( N.LT.0 ) THEN
231: INFO = -5
232: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
233: INFO = -7
234: ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
235: INFO = -9
236: ELSE
237: IF( VALEIG ) THEN
238: IF( N.GT.0 .AND. VU.LE.VL )
239: $ INFO = -11
240: ELSE IF( INDEIG ) THEN
241: IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
242: INFO = -12
243: ELSE IF( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
244: INFO = -13
245: END IF
246: END IF
247: END IF
248: IF (INFO.EQ.0) THEN
249: IF (LDZ.LT.1 .OR. (WANTZ .AND. LDZ.LT.N)) THEN
250: INFO = -18
251: END IF
252: END IF
253: *
254: IF( INFO.EQ.0 ) THEN
255: NB = ILAENV( 1, 'ZHETRD', UPLO, N, -1, -1, -1 )
256: LWKOPT = MAX( 1, ( NB + 1 )*N )
257: WORK( 1 ) = LWKOPT
258: *
259: IF( LWORK.LT.MAX( 1, 2*N ) .AND. .NOT.LQUERY ) THEN
260: INFO = -20
261: END IF
262: END IF
263: *
264: IF( INFO.NE.0 ) THEN
265: CALL XERBLA( 'ZHEGVX', -INFO )
266: RETURN
267: ELSE IF( LQUERY ) THEN
268: RETURN
269: END IF
270: *
271: * Quick return if possible
272: *
273: M = 0
274: IF( N.EQ.0 ) THEN
275: RETURN
276: END IF
277: *
278: * Form a Cholesky factorization of B.
279: *
280: CALL ZPOTRF( UPLO, N, B, LDB, INFO )
281: IF( INFO.NE.0 ) THEN
282: INFO = N + INFO
283: RETURN
284: END IF
285: *
286: * Transform problem to standard eigenvalue problem and solve.
287: *
288: CALL ZHEGST( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
289: CALL ZHEEVX( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL,
290: $ M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL,
291: $ INFO )
292: *
293: IF( WANTZ ) THEN
294: *
295: * Backtransform eigenvectors to the original problem.
296: *
297: IF( INFO.GT.0 )
298: $ M = INFO - 1
299: IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
300: *
301: * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
302: * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
303: *
304: IF( UPPER ) THEN
305: TRANS = 'N'
306: ELSE
307: TRANS = 'C'
308: END IF
309: *
310: CALL ZTRSM( 'Left', UPLO, TRANS, 'Non-unit', N, M, CONE, B,
311: $ LDB, Z, LDZ )
312: *
313: ELSE IF( ITYPE.EQ.3 ) THEN
314: *
315: * For B*A*x=(lambda)*x;
316: * backtransform eigenvectors: x = L*y or U**H *y
317: *
318: IF( UPPER ) THEN
319: TRANS = 'C'
320: ELSE
321: TRANS = 'N'
322: END IF
323: *
324: CALL ZTRMM( 'Left', UPLO, TRANS, 'Non-unit', N, M, CONE, B,
325: $ LDB, Z, LDZ )
326: END IF
327: END IF
328: *
329: * Set WORK(1) to optimal complex workspace size.
330: *
331: WORK( 1 ) = LWKOPT
332: *
333: RETURN
334: *
335: * End of ZHEGVX
336: *
337: END
CVSweb interface <joel.bertrand@systella.fr>