1: SUBROUTINE ZHPGV( ITYPE, JOBZ, UPLO, N, AP, BP, W, Z, LDZ, WORK,
2: $ RWORK, INFO )
3: *
4: * -- LAPACK driver routine (version 3.2) --
5: * -- LAPACK is a software package provided by Univ. of Tennessee, --
6: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
7: * November 2006
8: *
9: * .. Scalar Arguments ..
10: CHARACTER JOBZ, UPLO
11: INTEGER INFO, ITYPE, LDZ, N
12: * ..
13: * .. Array Arguments ..
14: DOUBLE PRECISION RWORK( * ), W( * )
15: COMPLEX*16 AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
16: * ..
17: *
18: * Purpose
19: * =======
20: *
21: * ZHPGV computes all the eigenvalues and, optionally, the eigenvectors
22: * of a complex generalized Hermitian-definite eigenproblem, of the form
23: * A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
24: * Here A and B are assumed to be Hermitian, stored in packed format,
25: * and B is also positive definite.
26: *
27: * Arguments
28: * =========
29: *
30: * ITYPE (input) INTEGER
31: * Specifies the problem type to be solved:
32: * = 1: A*x = (lambda)*B*x
33: * = 2: A*B*x = (lambda)*x
34: * = 3: B*A*x = (lambda)*x
35: *
36: * JOBZ (input) CHARACTER*1
37: * = 'N': Compute eigenvalues only;
38: * = 'V': Compute eigenvalues and eigenvectors.
39: *
40: * UPLO (input) CHARACTER*1
41: * = 'U': Upper triangles of A and B are stored;
42: * = 'L': Lower triangles of A and B are stored.
43: *
44: * N (input) INTEGER
45: * The order of the matrices A and B. N >= 0.
46: *
47: * AP (input/output) COMPLEX*16 array, dimension (N*(N+1)/2)
48: * On entry, the upper or lower triangle of the Hermitian matrix
49: * A, packed columnwise in a linear array. The j-th column of A
50: * is stored in the array AP as follows:
51: * if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
52: * if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
53: *
54: * On exit, the contents of AP are destroyed.
55: *
56: * BP (input/output) COMPLEX*16 array, dimension (N*(N+1)/2)
57: * On entry, the upper or lower triangle of the Hermitian matrix
58: * B, packed columnwise in a linear array. The j-th column of B
59: * is stored in the array BP as follows:
60: * if UPLO = 'U', BP(i + (j-1)*j/2) = B(i,j) for 1<=i<=j;
61: * if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.
62: *
63: * On exit, the triangular factor U or L from the Cholesky
64: * factorization B = U**H*U or B = L*L**H, in the same storage
65: * format as B.
66: *
67: * W (output) DOUBLE PRECISION array, dimension (N)
68: * If INFO = 0, the eigenvalues in ascending order.
69: *
70: * Z (output) COMPLEX*16 array, dimension (LDZ, N)
71: * If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
72: * eigenvectors. The eigenvectors are normalized as follows:
73: * if ITYPE = 1 or 2, Z**H*B*Z = I;
74: * if ITYPE = 3, Z**H*inv(B)*Z = I.
75: * If JOBZ = 'N', then Z is not referenced.
76: *
77: * LDZ (input) INTEGER
78: * The leading dimension of the array Z. LDZ >= 1, and if
79: * JOBZ = 'V', LDZ >= max(1,N).
80: *
81: * WORK (workspace) COMPLEX*16 array, dimension (max(1, 2*N-1))
82: *
83: * RWORK (workspace) DOUBLE PRECISION array, dimension (max(1, 3*N-2))
84: *
85: * INFO (output) INTEGER
86: * = 0: successful exit
87: * < 0: if INFO = -i, the i-th argument had an illegal value
88: * > 0: ZPPTRF or ZHPEV returned an error code:
89: * <= N: if INFO = i, ZHPEV failed to converge;
90: * i off-diagonal elements of an intermediate
91: * tridiagonal form did not convergeto zero;
92: * > N: if INFO = N + i, for 1 <= i <= n, then the leading
93: * minor of order i of B is not positive definite.
94: * The factorization of B could not be completed and
95: * no eigenvalues or eigenvectors were computed.
96: *
97: * =====================================================================
98: *
99: * .. Local Scalars ..
100: LOGICAL UPPER, WANTZ
101: CHARACTER TRANS
102: INTEGER J, NEIG
103: * ..
104: * .. External Functions ..
105: LOGICAL LSAME
106: EXTERNAL LSAME
107: * ..
108: * .. External Subroutines ..
109: EXTERNAL XERBLA, ZHPEV, ZHPGST, ZPPTRF, ZTPMV, ZTPSV
110: * ..
111: * .. Executable Statements ..
112: *
113: * Test the input parameters.
114: *
115: WANTZ = LSAME( JOBZ, 'V' )
116: UPPER = LSAME( UPLO, 'U' )
117: *
118: INFO = 0
119: IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
120: INFO = -1
121: ELSE IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
122: INFO = -2
123: ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
124: INFO = -3
125: ELSE IF( N.LT.0 ) THEN
126: INFO = -4
127: ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
128: INFO = -9
129: END IF
130: IF( INFO.NE.0 ) THEN
131: CALL XERBLA( 'ZHPGV ', -INFO )
132: RETURN
133: END IF
134: *
135: * Quick return if possible
136: *
137: IF( N.EQ.0 )
138: $ RETURN
139: *
140: * Form a Cholesky factorization of B.
141: *
142: CALL ZPPTRF( UPLO, N, BP, INFO )
143: IF( INFO.NE.0 ) THEN
144: INFO = N + INFO
145: RETURN
146: END IF
147: *
148: * Transform problem to standard eigenvalue problem and solve.
149: *
150: CALL ZHPGST( ITYPE, UPLO, N, AP, BP, INFO )
151: CALL ZHPEV( JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, RWORK, INFO )
152: *
153: IF( WANTZ ) THEN
154: *
155: * Backtransform eigenvectors to the original problem.
156: *
157: NEIG = N
158: IF( INFO.GT.0 )
159: $ NEIG = INFO - 1
160: IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
161: *
162: * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
163: * backtransform eigenvectors: x = inv(L)'*y or inv(U)*y
164: *
165: IF( UPPER ) THEN
166: TRANS = 'N'
167: ELSE
168: TRANS = 'C'
169: END IF
170: *
171: DO 10 J = 1, NEIG
172: CALL ZTPSV( UPLO, TRANS, 'Non-unit', N, BP, Z( 1, J ),
173: $ 1 )
174: 10 CONTINUE
175: *
176: ELSE IF( ITYPE.EQ.3 ) THEN
177: *
178: * For B*A*x=(lambda)*x;
179: * backtransform eigenvectors: x = L*y or U'*y
180: *
181: IF( UPPER ) THEN
182: TRANS = 'C'
183: ELSE
184: TRANS = 'N'
185: END IF
186: *
187: DO 20 J = 1, NEIG
188: CALL ZTPMV( UPLO, TRANS, 'Non-unit', N, BP, Z( 1, J ),
189: $ 1 )
190: 20 CONTINUE
191: END IF
192: END IF
193: RETURN
194: *
195: * End of ZHPGV
196: *
197: END
CVSweb interface <joel.bertrand@systella.fr>