LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhegvx ( integer  ITYPE,
character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldb, * )  B,
integer  LDB,
double precision  VL,
double precision  VU,
integer  IL,
integer  IU,
double precision  ABSTOL,
integer  M,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer, dimension( * )  IFAIL,
integer  INFO 
)

ZHEGVX

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

Purpose:
 ZHEGVX computes selected eigenvalues, and optionally, 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 and B is also positive definite.
 Eigenvalues and eigenvectors can be selected by specifying either a
 range of values or a range of indices for the desired eigenvalues.
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]RANGE
          RANGE is CHARACTER*1
          = 'A': all eigenvalues will be found.
          = 'V': all eigenvalues in the half-open interval (VL,VU]
                 will be found.
          = 'I': the IL-th through IU-th eigenvalues will be found.
[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]A
          A is COMPLEX*16 array, dimension (LDA, N)
          On entry, the Hermitian matrix A.  If UPLO = 'U', the
          leading N-by-N upper triangular part of A contains the
          upper triangular part of the matrix A.  If UPLO = 'L',
          the leading N-by-N lower triangular part of A contains
          the lower triangular part of the matrix A.

          On exit,  the lower triangle (if UPLO='L') or the upper
          triangle (if UPLO='U') of A, including the diagonal, is
          destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX*16 array, dimension (LDB, N)
          On entry, the Hermitian matrix B.  If UPLO = 'U', the
          leading N-by-N upper triangular part of B contains the
          upper triangular part of the matrix B.  If UPLO = 'L',
          the leading N-by-N lower triangular part of B contains
          the lower triangular part of the matrix B.

          On exit, if INFO <= N, the part of B containing the matrix is
          overwritten by the triangular factor U or L from the Cholesky
          factorization B = U**H*U or B = L*L**H.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in]VL
          VL is DOUBLE PRECISION

          If RANGE='V', the lower bound of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]VU
          VU is DOUBLE PRECISION

          If RANGE='V', the upper bound of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]IL
          IL is INTEGER

          If RANGE='I', the index of the
          smallest eigenvalue to be returned.
          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
          Not referenced if RANGE = 'A' or 'V'.
[in]IU
          IU is INTEGER

          If RANGE='I', the index of the
          largest eigenvalue to be returned.
          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
          Not referenced if RANGE = 'A' or 'V'.
[in]ABSTOL
          ABSTOL is DOUBLE PRECISION
          The absolute error tolerance for the eigenvalues.
          An approximate eigenvalue is accepted as converged
          when it is determined to lie in an interval [a,b]
          of width less than or equal to

                  ABSTOL + EPS *   max( |a|,|b| ) ,

          where EPS is the machine precision.  If ABSTOL is less than
          or equal to zero, then  EPS*|T|  will be used in its place,
          where |T| is the 1-norm of the tridiagonal matrix obtained
          by reducing C to tridiagonal form, where C is the symmetric
          matrix of the standard symmetric problem to which the
          generalized problem is transformed.

          Eigenvalues will be computed most accurately when ABSTOL is
          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
          If this routine returns with INFO>0, indicating that some
          eigenvectors did not converge, try setting ABSTOL to
          2*DLAMCH('S').
[out]M
          M is INTEGER
          The total number of eigenvalues found.  0 <= M <= N.
          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          The first M elements contain the selected
          eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, max(1,M))
          If JOBZ = 'N', then Z is not referenced.
          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
          contain the orthonormal eigenvectors of the matrix A
          corresponding to the selected eigenvalues, with the i-th
          column of Z holding the eigenvector associated with W(i).
          The eigenvectors are normalized as follows:
          if ITYPE = 1 or 2, Z**T*B*Z = I;
          if ITYPE = 3, Z**T*inv(B)*Z = I.

          If an eigenvector fails to converge, then that column of Z
          contains the latest approximation to the eigenvector, and the
          index of the eigenvector is returned in IFAIL.
          Note: the user must ensure that at least max(1,M) columns are
          supplied in the array Z; if RANGE = 'V', the exact value of M
          is not known in advance and an upper bound must be used.
[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 optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK.  LWORK >= max(1,2*N).
          For optimal efficiency, LWORK >= (NB+1)*N,
          where NB is the blocksize for ZHETRD returned by ILAENV.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (7*N)
[out]IWORK
          IWORK is INTEGER array, dimension (5*N)
[out]IFAIL
          IFAIL is INTEGER array, dimension (N)
          If JOBZ = 'V', then if INFO = 0, the first M elements of
          IFAIL are zero.  If INFO > 0, then IFAIL contains the
          indices of the eigenvectors that failed to converge.
          If JOBZ = 'N', then IFAIL is not referenced.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  ZPOTRF or ZHEEVX returned an error code:
             <= N:  if INFO = i, ZHEEVX failed to converge;
                    i eigenvectors failed to converge.  Their indices
                    are stored in array IFAIL.
             > 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
June 2016
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 309 of file zhegvx.f.

309 *
310 * -- LAPACK driver routine (version 3.6.1) --
311 * -- LAPACK is a software package provided by Univ. of Tennessee, --
312 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
313 * June 2016
314 *
315 * .. Scalar Arguments ..
316  CHARACTER jobz, range, uplo
317  INTEGER il, info, itype, iu, lda, ldb, ldz, lwork, m, n
318  DOUBLE PRECISION abstol, vl, vu
319 * ..
320 * .. Array Arguments ..
321  INTEGER ifail( * ), iwork( * )
322  DOUBLE PRECISION rwork( * ), w( * )
323  COMPLEX*16 a( lda, * ), b( ldb, * ), work( * ),
324  $ z( ldz, * )
325 * ..
326 *
327 * =====================================================================
328 *
329 * .. Parameters ..
330  COMPLEX*16 cone
331  parameter ( cone = ( 1.0d+0, 0.0d+0 ) )
332 * ..
333 * .. Local Scalars ..
334  LOGICAL alleig, indeig, lquery, upper, valeig, wantz
335  CHARACTER trans
336  INTEGER lwkopt, nb
337 * ..
338 * .. External Functions ..
339  LOGICAL lsame
340  INTEGER ilaenv
341  EXTERNAL lsame, ilaenv
342 * ..
343 * .. External Subroutines ..
344  EXTERNAL xerbla, zheevx, zhegst, zpotrf, ztrmm, ztrsm
345 * ..
346 * .. Intrinsic Functions ..
347  INTRINSIC max, min
348 * ..
349 * .. Executable Statements ..
350 *
351 * Test the input parameters.
352 *
353  wantz = lsame( jobz, 'V' )
354  upper = lsame( uplo, 'U' )
355  alleig = lsame( range, 'A' )
356  valeig = lsame( range, 'V' )
357  indeig = lsame( range, 'I' )
358  lquery = ( lwork.EQ.-1 )
359 *
360  info = 0
361  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
362  info = -1
363  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
364  info = -2
365  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
366  info = -3
367  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
368  info = -4
369  ELSE IF( n.LT.0 ) THEN
370  info = -5
371  ELSE IF( lda.LT.max( 1, n ) ) THEN
372  info = -7
373  ELSE IF( ldb.LT.max( 1, n ) ) THEN
374  info = -9
375  ELSE
376  IF( valeig ) THEN
377  IF( n.GT.0 .AND. vu.LE.vl )
378  $ info = -11
379  ELSE IF( indeig ) THEN
380  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
381  info = -12
382  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
383  info = -13
384  END IF
385  END IF
386  END IF
387  IF (info.EQ.0) THEN
388  IF (ldz.LT.1 .OR. (wantz .AND. ldz.LT.n)) THEN
389  info = -18
390  END IF
391  END IF
392 *
393  IF( info.EQ.0 ) THEN
394  nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
395  lwkopt = max( 1, ( nb + 1 )*n )
396  work( 1 ) = lwkopt
397 *
398  IF( lwork.LT.max( 1, 2*n ) .AND. .NOT.lquery ) THEN
399  info = -20
400  END IF
401  END IF
402 *
403  IF( info.NE.0 ) THEN
404  CALL xerbla( 'ZHEGVX', -info )
405  RETURN
406  ELSE IF( lquery ) THEN
407  RETURN
408  END IF
409 *
410 * Quick return if possible
411 *
412  m = 0
413  IF( n.EQ.0 ) THEN
414  RETURN
415  END IF
416 *
417 * Form a Cholesky factorization of B.
418 *
419  CALL zpotrf( uplo, n, b, ldb, info )
420  IF( info.NE.0 ) THEN
421  info = n + info
422  RETURN
423  END IF
424 *
425 * Transform problem to standard eigenvalue problem and solve.
426 *
427  CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
428  CALL zheevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol,
429  $ m, w, z, ldz, work, lwork, rwork, iwork, ifail,
430  $ info )
431 *
432  IF( wantz ) THEN
433 *
434 * Backtransform eigenvectors to the original problem.
435 *
436  IF( info.GT.0 )
437  $ m = info - 1
438  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
439 *
440 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
441 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
442 *
443  IF( upper ) THEN
444  trans = 'N'
445  ELSE
446  trans = 'C'
447  END IF
448 *
449  CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, m, cone, b,
450  $ ldb, z, ldz )
451 *
452  ELSE IF( itype.EQ.3 ) THEN
453 *
454 * For B*A*x=(lambda)*x;
455 * backtransform eigenvectors: x = L*y or U**H *y
456 *
457  IF( upper ) THEN
458  trans = 'C'
459  ELSE
460  trans = 'N'
461  END IF
462 *
463  CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, m, cone, b,
464  $ ldb, z, ldz )
465  END IF
466  END IF
467 *
468 * Set WORK(1) to optimal complex workspace size.
469 *
470  work( 1 ) = lwkopt
471 *
472  RETURN
473 *
474 * End of ZHEGVX
475 *
subroutine zpotrf(UPLO, N, A, LDA, INFO)
ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
Definition: zpotrf.f:102
subroutine zhegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
ZHEGST
Definition: zhegst.f:129
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ztrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRMM
Definition: ztrmm.f:179
subroutine zheevx(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
ZHEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: zheevx.f:261
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine ztrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRSM
Definition: ztrsm.f:182
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: