1: SUBROUTINE ZHBGVD( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
2: $ Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK,
3: $ LIWORK, 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: * @precisions normal z -> c
10: *
11: * .. Scalar Arguments ..
12: CHARACTER JOBZ, UPLO
13: INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LRWORK,
14: $ LWORK, N
15: * ..
16: * .. Array Arguments ..
17: INTEGER IWORK( * )
18: DOUBLE PRECISION RWORK( * ), W( * )
19: COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
20: $ Z( LDZ, * )
21: * ..
22: *
23: * Purpose
24: * =======
25: *
26: * ZHBGVD computes all the eigenvalues, and optionally, the eigenvectors
27: * of a complex generalized Hermitian-definite banded eigenproblem, of
28: * the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
29: * and banded, and B is also positive definite. If eigenvectors are
30: * desired, it uses a divide and conquer algorithm.
31: *
32: * The divide and conquer algorithm makes very mild assumptions about
33: * floating point arithmetic. It will work on machines with a guard
34: * digit in add/subtract, or on those binary machines without guard
35: * digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
36: * Cray-2. It could conceivably fail on hexadecimal or decimal machines
37: * without guard digits, but we know of none.
38: *
39: * Arguments
40: * =========
41: *
42: * JOBZ (input) CHARACTER*1
43: * = 'N': Compute eigenvalues only;
44: * = 'V': Compute eigenvalues and eigenvectors.
45: *
46: * UPLO (input) CHARACTER*1
47: * = 'U': Upper triangles of A and B are stored;
48: * = 'L': Lower triangles of A and B are stored.
49: *
50: * N (input) INTEGER
51: * The order of the matrices A and B. N >= 0.
52: *
53: * KA (input) INTEGER
54: * The number of superdiagonals of the matrix A if UPLO = 'U',
55: * or the number of subdiagonals if UPLO = 'L'. KA >= 0.
56: *
57: * KB (input) INTEGER
58: * The number of superdiagonals of the matrix B if UPLO = 'U',
59: * or the number of subdiagonals if UPLO = 'L'. KB >= 0.
60: *
61: * AB (input/output) COMPLEX*16 array, dimension (LDAB, N)
62: * On entry, the upper or lower triangle of the Hermitian band
63: * matrix A, stored in the first ka+1 rows of the array. The
64: * j-th column of A is stored in the j-th column of the array AB
65: * as follows:
66: * if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
67: * if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
68: *
69: * On exit, the contents of AB are destroyed.
70: *
71: * LDAB (input) INTEGER
72: * The leading dimension of the array AB. LDAB >= KA+1.
73: *
74: * BB (input/output) COMPLEX*16 array, dimension (LDBB, N)
75: * On entry, the upper or lower triangle of the Hermitian band
76: * matrix B, stored in the first kb+1 rows of the array. The
77: * j-th column of B is stored in the j-th column of the array BB
78: * as follows:
79: * if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
80: * if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
81: *
82: * On exit, the factor S from the split Cholesky factorization
83: * B = S**H*S, as returned by ZPBSTF.
84: *
85: * LDBB (input) INTEGER
86: * The leading dimension of the array BB. LDBB >= KB+1.
87: *
88: * W (output) DOUBLE PRECISION array, dimension (N)
89: * If INFO = 0, the eigenvalues in ascending order.
90: *
91: * Z (output) COMPLEX*16 array, dimension (LDZ, N)
92: * If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
93: * eigenvectors, with the i-th column of Z holding the
94: * eigenvector associated with W(i). The eigenvectors are
95: * normalized so that Z**H*B*Z = I.
96: * If JOBZ = 'N', then Z is not referenced.
97: *
98: * LDZ (input) INTEGER
99: * The leading dimension of the array Z. LDZ >= 1, and if
100: * JOBZ = 'V', LDZ >= N.
101: *
102: * WORK (workspace/output) COMPLEX*16 array, dimension (MAX(1,LWORK))
103: * On exit, if INFO=0, WORK(1) returns the optimal LWORK.
104: *
105: * LWORK (input) INTEGER
106: * The dimension of the array WORK.
107: * If N <= 1, LWORK >= 1.
108: * If JOBZ = 'N' and N > 1, LWORK >= N.
109: * If JOBZ = 'V' and N > 1, LWORK >= 2*N**2.
110: *
111: * If LWORK = -1, then a workspace query is assumed; the routine
112: * only calculates the optimal sizes of the WORK, RWORK and
113: * IWORK arrays, returns these values as the first entries of
114: * the WORK, RWORK and IWORK arrays, and no error message
115: * related to LWORK or LRWORK or LIWORK is issued by XERBLA.
116: *
117: * RWORK (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
118: * On exit, if INFO=0, RWORK(1) returns the optimal LRWORK.
119: *
120: * LRWORK (input) INTEGER
121: * The dimension of array RWORK.
122: * If N <= 1, LRWORK >= 1.
123: * If JOBZ = 'N' and N > 1, LRWORK >= N.
124: * If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
125: *
126: * If LRWORK = -1, then a workspace query is assumed; the
127: * routine only calculates the optimal sizes of the WORK, RWORK
128: * and IWORK arrays, returns these values as the first entries
129: * of the WORK, RWORK and IWORK arrays, and no error message
130: * related to LWORK or LRWORK or LIWORK is issued by XERBLA.
131: *
132: * IWORK (workspace/output) INTEGER array, dimension (MAX(1,LIWORK))
133: * On exit, if INFO=0, IWORK(1) returns the optimal LIWORK.
134: *
135: * LIWORK (input) INTEGER
136: * The dimension of array IWORK.
137: * If JOBZ = 'N' or N <= 1, LIWORK >= 1.
138: * If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.
139: *
140: * If LIWORK = -1, then a workspace query is assumed; the
141: * routine only calculates the optimal sizes of the WORK, RWORK
142: * and IWORK arrays, returns these values as the first entries
143: * of the WORK, RWORK and IWORK arrays, and no error message
144: * related to LWORK or LRWORK or LIWORK is issued by XERBLA.
145: *
146: * INFO (output) INTEGER
147: * = 0: successful exit
148: * < 0: if INFO = -i, the i-th argument had an illegal value
149: * > 0: if INFO = i, and i is:
150: * <= N: the algorithm failed to converge:
151: * i off-diagonal elements of an intermediate
152: * tridiagonal form did not converge to zero;
153: * > N: if INFO = N + i, for 1 <= i <= N, then ZPBSTF
154: * returned INFO = i: B is not positive definite.
155: * The factorization of B could not be completed and
156: * no eigenvalues or eigenvectors were computed.
157: *
158: * Further Details
159: * ===============
160: *
161: * Based on contributions by
162: * Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
163: *
164: * =====================================================================
165: *
166: * .. Parameters ..
167: COMPLEX*16 CONE, CZERO
168: PARAMETER ( CONE = ( 1.0D+0, 0.0D+0 ),
169: $ CZERO = ( 0.0D+0, 0.0D+0 ) )
170: * ..
171: * .. Local Scalars ..
172: LOGICAL LQUERY, UPPER, WANTZ
173: CHARACTER VECT
174: INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLRWK,
175: $ LLWK2, LRWMIN, LWMIN
176: * ..
177: * .. External Functions ..
178: LOGICAL LSAME
179: EXTERNAL LSAME
180: * ..
181: * .. External Subroutines ..
182: EXTERNAL DSTERF, XERBLA, ZGEMM, ZHBGST, ZHBTRD, ZLACPY,
183: $ ZPBSTF, ZSTEDC
184: * ..
185: * .. Executable Statements ..
186: *
187: * Test the input parameters.
188: *
189: WANTZ = LSAME( JOBZ, 'V' )
190: UPPER = LSAME( UPLO, 'U' )
191: LQUERY = ( LWORK.EQ.-1 .OR. LRWORK.EQ.-1 .OR. LIWORK.EQ.-1 )
192: *
193: INFO = 0
194: IF( N.LE.1 ) THEN
195: LWMIN = 1+N
196: LRWMIN = 1+N
197: LIWMIN = 1
198: ELSE IF( WANTZ ) THEN
199: LWMIN = 2*N**2
200: LRWMIN = 1 + 5*N + 2*N**2
201: LIWMIN = 3 + 5*N
202: ELSE
203: LWMIN = N
204: LRWMIN = N
205: LIWMIN = 1
206: END IF
207: IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
208: INFO = -1
209: ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
210: INFO = -2
211: ELSE IF( N.LT.0 ) THEN
212: INFO = -3
213: ELSE IF( KA.LT.0 ) THEN
214: INFO = -4
215: ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN
216: INFO = -5
217: ELSE IF( LDAB.LT.KA+1 ) THEN
218: INFO = -7
219: ELSE IF( LDBB.LT.KB+1 ) THEN
220: INFO = -9
221: ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
222: INFO = -12
223: END IF
224: *
225: IF( INFO.EQ.0 ) THEN
226: WORK( 1 ) = LWMIN
227: RWORK( 1 ) = LRWMIN
228: IWORK( 1 ) = LIWMIN
229: *
230: IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
231: INFO = -14
232: ELSE IF( LRWORK.LT.LRWMIN .AND. .NOT.LQUERY ) THEN
233: INFO = -16
234: ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN
235: INFO = -18
236: END IF
237: END IF
238: *
239: IF( INFO.NE.0 ) THEN
240: CALL XERBLA( 'ZHBGVD', -INFO )
241: RETURN
242: ELSE IF( LQUERY ) THEN
243: RETURN
244: END IF
245: *
246: * Quick return if possible
247: *
248: IF( N.EQ.0 )
249: $ RETURN
250: *
251: * Form a split Cholesky factorization of B.
252: *
253: CALL ZPBSTF( UPLO, N, KB, BB, LDBB, INFO )
254: IF( INFO.NE.0 ) THEN
255: INFO = N + INFO
256: RETURN
257: END IF
258: *
259: * Transform problem to standard eigenvalue problem.
260: *
261: INDE = 1
262: INDWRK = INDE + N
263: INDWK2 = 1 + N*N
264: LLWK2 = LWORK - INDWK2 + 2
265: LLRWK = LRWORK - INDWRK + 2
266: CALL ZHBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Z, LDZ,
267: $ WORK, RWORK( INDWRK ), IINFO )
268: *
269: * Reduce Hermitian band matrix to tridiagonal form.
270: *
271: IF( WANTZ ) THEN
272: VECT = 'U'
273: ELSE
274: VECT = 'N'
275: END IF
276: CALL ZHBTRD( VECT, UPLO, N, KA, AB, LDAB, W, RWORK( INDE ), Z,
277: $ LDZ, WORK, IINFO )
278: *
279: * For eigenvalues only, call DSTERF. For eigenvectors, call ZSTEDC.
280: *
281: IF( .NOT.WANTZ ) THEN
282: CALL DSTERF( N, W, RWORK( INDE ), INFO )
283: ELSE
284: CALL ZSTEDC( 'I', N, W, RWORK( INDE ), WORK, N, WORK( INDWK2 ),
285: $ LLWK2, RWORK( INDWRK ), LLRWK, IWORK, LIWORK,
286: $ INFO )
287: CALL ZGEMM( 'N', 'N', N, N, N, CONE, Z, LDZ, WORK, N, CZERO,
288: $ WORK( INDWK2 ), N )
289: CALL ZLACPY( 'A', N, N, WORK( INDWK2 ), N, Z, LDZ )
290: END IF
291: *
292: WORK( 1 ) = LWMIN
293: RWORK( 1 ) = LRWMIN
294: IWORK( 1 ) = LIWMIN
295: RETURN
296: *
297: * End of ZHBGVD
298: *
299: END
CVSweb interface <joel.bertrand@systella.fr>