LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhpgv ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
complex*16, dimension( * )  AP,
complex*16, dimension( * )  BP,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer  INFO 
)

ZHPGV

Download ZHPGV + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZHPGV computes all the eigenvalues and, optionally, the eigenvectors
 of a complex generalized Hermitian-definite eigenproblem, of the form
 A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
 Here A and B are assumed to be Hermitian, stored in packed format,
 and B is also positive definite.
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the problem type to be solved:
          = 1:  A*x = (lambda)*B*x
          = 2:  A*B*x = (lambda)*x
          = 3:  B*A*x = (lambda)*x
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in,out]AP
          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, the contents of AP are destroyed.
[in,out]BP
          BP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          B, packed columnwise in a linear array.  The j-th column of B
          is stored in the array BP as follows:
          if UPLO = 'U', BP(i + (j-1)*j/2) = B(i,j) for 1<=i<=j;
          if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.

          On exit, the triangular factor U or L from the Cholesky
          factorization B = U**H*U or B = L*L**H, in the same storage
          format as B.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
          eigenvectors.  The eigenvectors are normalized as follows:
          if ITYPE = 1 or 2, Z**H*B*Z = I;
          if ITYPE = 3, Z**H*inv(B)*Z = I.
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is COMPLEX*16 array, dimension (max(1, 2*N-1))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(1, 3*N-2))
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  ZPPTRF or ZHPEV returned an error code:
             <= N:  if INFO = i, ZHPEV failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not convergeto zero;
             > N:   if INFO = N + i, for 1 <= i <= n, then the leading
                    minor of order i of B is not positive definite.
                    The factorization of B could not be completed and
                    no eigenvalues or eigenvectors were computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 167 of file zhpgv.f.

167 *
168 * -- LAPACK driver routine (version 3.6.0) --
169 * -- LAPACK is a software package provided by Univ. of Tennessee, --
170 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
171 * November 2015
172 *
173 * .. Scalar Arguments ..
174  CHARACTER jobz, uplo
175  INTEGER info, itype, ldz, n
176 * ..
177 * .. Array Arguments ..
178  DOUBLE PRECISION rwork( * ), w( * )
179  COMPLEX*16 ap( * ), bp( * ), work( * ), z( ldz, * )
180 * ..
181 *
182 * =====================================================================
183 *
184 * .. Local Scalars ..
185  LOGICAL upper, wantz
186  CHARACTER trans
187  INTEGER j, neig
188 * ..
189 * .. External Functions ..
190  LOGICAL lsame
191  EXTERNAL lsame
192 * ..
193 * .. External Subroutines ..
194  EXTERNAL xerbla, zhpev, zhpgst, zpptrf, ztpmv, ztpsv
195 * ..
196 * .. Executable Statements ..
197 *
198 * Test the input parameters.
199 *
200  wantz = lsame( jobz, 'V' )
201  upper = lsame( uplo, 'U' )
202 *
203  info = 0
204  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
205  info = -1
206  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
207  info = -2
208  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
209  info = -3
210  ELSE IF( n.LT.0 ) THEN
211  info = -4
212  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
213  info = -9
214  END IF
215  IF( info.NE.0 ) THEN
216  CALL xerbla( 'ZHPGV ', -info )
217  RETURN
218  END IF
219 *
220 * Quick return if possible
221 *
222  IF( n.EQ.0 )
223  $ RETURN
224 *
225 * Form a Cholesky factorization of B.
226 *
227  CALL zpptrf( uplo, n, bp, info )
228  IF( info.NE.0 ) THEN
229  info = n + info
230  RETURN
231  END IF
232 *
233 * Transform problem to standard eigenvalue problem and solve.
234 *
235  CALL zhpgst( itype, uplo, n, ap, bp, info )
236  CALL zhpev( jobz, uplo, n, ap, w, z, ldz, work, rwork, info )
237 *
238  IF( wantz ) THEN
239 *
240 * Backtransform eigenvectors to the original problem.
241 *
242  neig = n
243  IF( info.GT.0 )
244  $ neig = info - 1
245  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
246 *
247 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
248 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
249 *
250  IF( upper ) THEN
251  trans = 'N'
252  ELSE
253  trans = 'C'
254  END IF
255 *
256  DO 10 j = 1, neig
257  CALL ztpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
258  $ 1 )
259  10 CONTINUE
260 *
261  ELSE IF( itype.EQ.3 ) THEN
262 *
263 * For B*A*x=(lambda)*x;
264 * backtransform eigenvectors: x = L*y or U**H *y
265 *
266  IF( upper ) THEN
267  trans = 'C'
268  ELSE
269  trans = 'N'
270  END IF
271 *
272  DO 20 j = 1, neig
273  CALL ztpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
274  $ 1 )
275  20 CONTINUE
276  END IF
277  END IF
278  RETURN
279 *
280 * End of ZHPGV
281 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ztpsv(UPLO, TRANS, DIAG, N, AP, X, INCX)
ZTPSV
Definition: ztpsv.f:146
subroutine ztpmv(UPLO, TRANS, DIAG, N, AP, X, INCX)
ZTPMV
Definition: ztpmv.f:144
subroutine zhpev(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, RWORK, INFO)
ZHPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: zhpev.f:140
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zpptrf(UPLO, N, AP, INFO)
ZPPTRF
Definition: zpptrf.f:121
subroutine zhpgst(ITYPE, UPLO, N, AP, BP, INFO)
ZHPGST
Definition: zhpgst.f:115

Here is the call graph for this function:

Here is the caller graph for this function: