1*> \brief \b ZHPGV
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZHPGV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhpgv.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhpgv.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhpgv.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18*  Definition:
19*  ===========
20*
21*       SUBROUTINE ZHPGV( ITYPE, JOBZ, UPLO, N, AP, BP, W, Z, LDZ, WORK,
22*                         RWORK, INFO )
23*
24*       .. Scalar Arguments ..
25*       CHARACTER          JOBZ, UPLO
26*       INTEGER            INFO, ITYPE, LDZ, N
27*       ..
28*       .. Array Arguments ..
29*       DOUBLE PRECISION   RWORK( * ), W( * )
30*       COMPLEX*16         AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
31*       ..
32*
33*
34*> \par Purpose:
35*  =============
36*>
37*> \verbatim
38*>
39*> ZHPGV computes all the eigenvalues and, optionally, the eigenvectors
40*> of a complex generalized Hermitian-definite eigenproblem, of the form
41*> A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
42*> Here A and B are assumed to be Hermitian, stored in packed format,
43*> and B is also positive definite.
44*> \endverbatim
45*
46*  Arguments:
47*  ==========
48*
49*> \param[in] ITYPE
50*> \verbatim
51*>          ITYPE is INTEGER
52*>          Specifies the problem type to be solved:
53*>          = 1:  A*x = (lambda)*B*x
54*>          = 2:  A*B*x = (lambda)*x
55*>          = 3:  B*A*x = (lambda)*x
56*> \endverbatim
57*>
58*> \param[in] JOBZ
59*> \verbatim
60*>          JOBZ is CHARACTER*1
61*>          = 'N':  Compute eigenvalues only;
62*>          = 'V':  Compute eigenvalues and eigenvectors.
63*> \endverbatim
64*>
65*> \param[in] UPLO
66*> \verbatim
67*>          UPLO is CHARACTER*1
68*>          = 'U':  Upper triangles of A and B are stored;
69*>          = 'L':  Lower triangles of A and B are stored.
70*> \endverbatim
71*>
72*> \param[in] N
73*> \verbatim
74*>          N is INTEGER
75*>          The order of the matrices A and B.  N >= 0.
76*> \endverbatim
77*>
78*> \param[in,out] AP
79*> \verbatim
80*>          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
81*>          On entry, the upper or lower triangle of the Hermitian matrix
82*>          A, packed columnwise in a linear array.  The j-th column of A
83*>          is stored in the array AP as follows:
84*>          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
85*>          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
86*>
87*>          On exit, the contents of AP are destroyed.
88*> \endverbatim
89*>
90*> \param[in,out] BP
91*> \verbatim
92*>          BP is COMPLEX*16 array, dimension (N*(N+1)/2)
93*>          On entry, the upper or lower triangle of the Hermitian matrix
94*>          B, packed columnwise in a linear array.  The j-th column of B
95*>          is stored in the array BP as follows:
96*>          if UPLO = 'U', BP(i + (j-1)*j/2) = B(i,j) for 1<=i<=j;
97*>          if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.
98*>
99*>          On exit, the triangular factor U or L from the Cholesky
100*>          factorization B = U**H*U or B = L*L**H, in the same storage
101*>          format as B.
102*> \endverbatim
103*>
104*> \param[out] W
105*> \verbatim
106*>          W is DOUBLE PRECISION array, dimension (N)
107*>          If INFO = 0, the eigenvalues in ascending order.
108*> \endverbatim
109*>
110*> \param[out] Z
111*> \verbatim
112*>          Z is COMPLEX*16 array, dimension (LDZ, N)
113*>          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
114*>          eigenvectors.  The eigenvectors are normalized as follows:
115*>          if ITYPE = 1 or 2, Z**H*B*Z = I;
116*>          if ITYPE = 3, Z**H*inv(B)*Z = I.
117*>          If JOBZ = 'N', then Z is not referenced.
118*> \endverbatim
119*>
120*> \param[in] LDZ
121*> \verbatim
122*>          LDZ is INTEGER
123*>          The leading dimension of the array Z.  LDZ >= 1, and if
124*>          JOBZ = 'V', LDZ >= max(1,N).
125*> \endverbatim
126*>
127*> \param[out] WORK
128*> \verbatim
129*>          WORK is COMPLEX*16 array, dimension (max(1, 2*N-1))
130*> \endverbatim
131*>
132*> \param[out] RWORK
133*> \verbatim
134*>          RWORK is DOUBLE PRECISION array, dimension (max(1, 3*N-2))
135*> \endverbatim
136*>
137*> \param[out] INFO
138*> \verbatim
139*>          INFO is INTEGER
140*>          = 0:  successful exit
141*>          < 0:  if INFO = -i, the i-th argument had an illegal value
142*>          > 0:  ZPPTRF or ZHPEV returned an error code:
143*>             <= N:  if INFO = i, ZHPEV failed to converge;
144*>                    i off-diagonal elements of an intermediate
145*>                    tridiagonal form did not convergeto zero;
146*>             > N:   if INFO = N + i, for 1 <= i <= n, then the leading
147*>                    minor of order i of B is not positive definite.
148*>                    The factorization of B could not be completed and
149*>                    no eigenvalues or eigenvectors were computed.
150*> \endverbatim
151*
152*  Authors:
153*  ========
154*
155*> \author Univ. of Tennessee
156*> \author Univ. of California Berkeley
157*> \author Univ. of Colorado Denver
158*> \author NAG Ltd.
159*
160*> \ingroup complex16OTHEReigen
161*
162*  =====================================================================
163      SUBROUTINE ZHPGV( ITYPE, JOBZ, UPLO, N, AP, BP, W, Z, LDZ, WORK,
164     $                  RWORK, INFO )
165*
166*  -- LAPACK driver routine --
167*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
168*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
169*
170*     .. Scalar Arguments ..
171      CHARACTER          JOBZ, UPLO
172      INTEGER            INFO, ITYPE, LDZ, N
173*     ..
174*     .. Array Arguments ..
175      DOUBLE PRECISION   RWORK( * ), W( * )
176      COMPLEX*16         AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
177*     ..
178*
179*  =====================================================================
180*
181*     .. Local Scalars ..
182      LOGICAL            UPPER, WANTZ
183      CHARACTER          TRANS
184      INTEGER            J, NEIG
185*     ..
186*     .. External Functions ..
187      LOGICAL            LSAME
188      EXTERNAL           LSAME
189*     ..
190*     .. External Subroutines ..
191      EXTERNAL           XERBLA, ZHPEV, ZHPGST, ZPPTRF, ZTPMV, ZTPSV
192*     ..
193*     .. Executable Statements ..
194*
195*     Test the input parameters.
196*
197      WANTZ = LSAME( JOBZ, 'V' )
198      UPPER = LSAME( UPLO, 'U' )
199*
200      INFO = 0
201      IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
202         INFO = -1
203      ELSE IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
204         INFO = -2
205      ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
206         INFO = -3
207      ELSE IF( N.LT.0 ) THEN
208         INFO = -4
209      ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
210         INFO = -9
211      END IF
212      IF( INFO.NE.0 ) THEN
213         CALL XERBLA( 'ZHPGV ', -INFO )
214         RETURN
215      END IF
216*
217*     Quick return if possible
218*
219      IF( N.EQ.0 )
220     $   RETURN
221*
222*     Form a Cholesky factorization of B.
223*
224      CALL ZPPTRF( UPLO, N, BP, INFO )
225      IF( INFO.NE.0 ) THEN
226         INFO = N + INFO
227         RETURN
228      END IF
229*
230*     Transform problem to standard eigenvalue problem and solve.
231*
232      CALL ZHPGST( ITYPE, UPLO, N, AP, BP, INFO )
233      CALL ZHPEV( JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, RWORK, INFO )
234*
235      IF( WANTZ ) THEN
236*
237*        Backtransform eigenvectors to the original problem.
238*
239         NEIG = N
240         IF( INFO.GT.0 )
241     $      NEIG = INFO - 1
242         IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
243*
244*           For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
245*           backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
246*
247            IF( UPPER ) THEN
248               TRANS = 'N'
249            ELSE
250               TRANS = 'C'
251            END IF
252*
253            DO 10 J = 1, NEIG
254               CALL ZTPSV( UPLO, TRANS, 'Non-unit', N, BP, Z( 1, J ),
255     $                     1 )
256   10       CONTINUE
257*
258         ELSE IF( ITYPE.EQ.3 ) THEN
259*
260*           For B*A*x=(lambda)*x;
261*           backtransform eigenvectors: x = L*y or U**H *y
262*
263            IF( UPPER ) THEN
264               TRANS = 'C'
265            ELSE
266               TRANS = 'N'
267            END IF
268*
269            DO 20 J = 1, NEIG
270               CALL ZTPMV( UPLO, TRANS, 'Non-unit', N, BP, Z( 1, J ),
271     $                     1 )
272   20       CONTINUE
273         END IF
274      END IF
275      RETURN
276*
277*     End of ZHPGV
278*
279      END
280