LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhpgvd ( 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,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

ZHPGVD

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

Purpose:
 ZHPGVD 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.
 If eigenvectors are desired, it uses a divide and conquer algorithm.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
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,LWORK))
          On exit, if INFO = 0, WORK(1) returns the required LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= N.
          If JOBZ = 'V' and N > 1, LWORK >= 2*N.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the required sizes of the WORK, RWORK and
          IWORK arrays, returns these values as the first entries of
          the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
          On exit, if INFO = 0, RWORK(1) returns the required LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The dimension of array RWORK.
          If N <= 1,               LRWORK >= 1.
          If JOBZ = 'N' and N > 1, LRWORK >= N.
          If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.

          If LRWORK = -1, then a workspace query is assumed; the
          routine only calculates the required sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the required LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of array IWORK.
          If JOBZ  = 'N' or N <= 1, LIWORK >= 1.
          If JOBZ  = 'V' and N > 1, LIWORK >= 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the required sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  ZPPTRF or ZHPEVD returned an error code:
             <= N:  if INFO = i, ZHPEVD 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
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 233 of file zhpgvd.f.

233 *
234 * -- LAPACK driver routine (version 3.6.0) --
235 * -- LAPACK is a software package provided by Univ. of Tennessee, --
236 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
237 * November 2015
238 *
239 * .. Scalar Arguments ..
240  CHARACTER jobz, uplo
241  INTEGER info, itype, ldz, liwork, lrwork, lwork, n
242 * ..
243 * .. Array Arguments ..
244  INTEGER iwork( * )
245  DOUBLE PRECISION rwork( * ), w( * )
246  COMPLEX*16 ap( * ), bp( * ), work( * ), z( ldz, * )
247 * ..
248 *
249 * =====================================================================
250 *
251 * .. Local Scalars ..
252  LOGICAL lquery, upper, wantz
253  CHARACTER trans
254  INTEGER j, liwmin, lrwmin, lwmin, neig
255 * ..
256 * .. External Functions ..
257  LOGICAL lsame
258  EXTERNAL lsame
259 * ..
260 * .. External Subroutines ..
261  EXTERNAL xerbla, zhpevd, zhpgst, zpptrf, ztpmv, ztpsv
262 * ..
263 * .. Intrinsic Functions ..
264  INTRINSIC dble, max
265 * ..
266 * .. Executable Statements ..
267 *
268 * Test the input parameters.
269 *
270  wantz = lsame( jobz, 'V' )
271  upper = lsame( uplo, 'U' )
272  lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
273 *
274  info = 0
275  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
276  info = -1
277  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
278  info = -2
279  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
280  info = -3
281  ELSE IF( n.LT.0 ) THEN
282  info = -4
283  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
284  info = -9
285  END IF
286 *
287  IF( info.EQ.0 ) THEN
288  IF( n.LE.1 ) THEN
289  lwmin = 1
290  liwmin = 1
291  lrwmin = 1
292  ELSE
293  IF( wantz ) THEN
294  lwmin = 2*n
295  lrwmin = 1 + 5*n + 2*n**2
296  liwmin = 3 + 5*n
297  ELSE
298  lwmin = n
299  lrwmin = n
300  liwmin = 1
301  END IF
302  END IF
303 *
304  work( 1 ) = lwmin
305  rwork( 1 ) = lrwmin
306  iwork( 1 ) = liwmin
307  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
308  info = -11
309  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
310  info = -13
311  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
312  info = -15
313  END IF
314  END IF
315 *
316  IF( info.NE.0 ) THEN
317  CALL xerbla( 'ZHPGVD', -info )
318  RETURN
319  ELSE IF( lquery ) THEN
320  RETURN
321  END IF
322 *
323 * Quick return if possible
324 *
325  IF( n.EQ.0 )
326  $ RETURN
327 *
328 * Form a Cholesky factorization of B.
329 *
330  CALL zpptrf( uplo, n, bp, info )
331  IF( info.NE.0 ) THEN
332  info = n + info
333  RETURN
334  END IF
335 *
336 * Transform problem to standard eigenvalue problem and solve.
337 *
338  CALL zhpgst( itype, uplo, n, ap, bp, info )
339  CALL zhpevd( jobz, uplo, n, ap, w, z, ldz, work, lwork, rwork,
340  $ lrwork, iwork, liwork, info )
341  lwmin = max( dble( lwmin ), dble( work( 1 ) ) )
342  lrwmin = max( dble( lrwmin ), dble( rwork( 1 ) ) )
343  liwmin = max( dble( liwmin ), dble( iwork( 1 ) ) )
344 *
345  IF( wantz ) THEN
346 *
347 * Backtransform eigenvectors to the original problem.
348 *
349  neig = n
350  IF( info.GT.0 )
351  $ neig = info - 1
352  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
353 *
354 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
355 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
356 *
357  IF( upper ) THEN
358  trans = 'N'
359  ELSE
360  trans = 'C'
361  END IF
362 *
363  DO 10 j = 1, neig
364  CALL ztpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
365  $ 1 )
366  10 CONTINUE
367 *
368  ELSE IF( itype.EQ.3 ) THEN
369 *
370 * For B*A*x=(lambda)*x;
371 * backtransform eigenvectors: x = L*y or U**H *y
372 *
373  IF( upper ) THEN
374  trans = 'C'
375  ELSE
376  trans = 'N'
377  END IF
378 *
379  DO 20 j = 1, neig
380  CALL ztpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
381  $ 1 )
382  20 CONTINUE
383  END IF
384  END IF
385 *
386  work( 1 ) = lwmin
387  rwork( 1 ) = lrwmin
388  iwork( 1 ) = liwmin
389  RETURN
390 *
391 * End of ZHPGVD
392 *
subroutine zhpevd(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZHPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matric...
Definition: zhpevd.f:203
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
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: