LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zhegvx()

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
                    principal minor of order i of B is not positive.
                    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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 304 of file zhegvx.f.

307*
308* -- LAPACK driver routine --
309* -- LAPACK is a software package provided by Univ. of Tennessee, --
310* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
311*
312* .. Scalar Arguments ..
313 CHARACTER JOBZ, RANGE, UPLO
314 INTEGER IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N
315 DOUBLE PRECISION ABSTOL, VL, VU
316* ..
317* .. Array Arguments ..
318 INTEGER IFAIL( * ), IWORK( * )
319 DOUBLE PRECISION RWORK( * ), W( * )
320 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * ),
321 $ Z( LDZ, * )
322* ..
323*
324* =====================================================================
325*
326* .. Parameters ..
327 COMPLEX*16 CONE
328 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
329* ..
330* .. Local Scalars ..
331 LOGICAL ALLEIG, INDEIG, LQUERY, UPPER, VALEIG, WANTZ
332 CHARACTER TRANS
333 INTEGER LWKOPT, NB
334* ..
335* .. External Functions ..
336 LOGICAL LSAME
337 INTEGER ILAENV
338 EXTERNAL lsame, ilaenv
339* ..
340* .. External Subroutines ..
341 EXTERNAL xerbla, zheevx, zhegst, zpotrf, ztrmm, ztrsm
342* ..
343* .. Intrinsic Functions ..
344 INTRINSIC max, min
345* ..
346* .. Executable Statements ..
347*
348* Test the input parameters.
349*
350 wantz = lsame( jobz, 'V' )
351 upper = lsame( uplo, 'U' )
352 alleig = lsame( range, 'A' )
353 valeig = lsame( range, 'V' )
354 indeig = lsame( range, 'I' )
355 lquery = ( lwork.EQ.-1 )
356*
357 info = 0
358 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
359 info = -1
360 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
361 info = -2
362 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
363 info = -3
364 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
365 info = -4
366 ELSE IF( n.LT.0 ) THEN
367 info = -5
368 ELSE IF( lda.LT.max( 1, n ) ) THEN
369 info = -7
370 ELSE IF( ldb.LT.max( 1, n ) ) THEN
371 info = -9
372 ELSE
373 IF( valeig ) THEN
374 IF( n.GT.0 .AND. vu.LE.vl )
375 $ info = -11
376 ELSE IF( indeig ) THEN
377 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
378 info = -12
379 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
380 info = -13
381 END IF
382 END IF
383 END IF
384 IF (info.EQ.0) THEN
385 IF (ldz.LT.1 .OR. (wantz .AND. ldz.LT.n)) THEN
386 info = -18
387 END IF
388 END IF
389*
390 IF( info.EQ.0 ) THEN
391 nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
392 lwkopt = max( 1, ( nb + 1 )*n )
393 work( 1 ) = lwkopt
394*
395 IF( lwork.LT.max( 1, 2*n ) .AND. .NOT.lquery ) THEN
396 info = -20
397 END IF
398 END IF
399*
400 IF( info.NE.0 ) THEN
401 CALL xerbla( 'ZHEGVX', -info )
402 RETURN
403 ELSE IF( lquery ) THEN
404 RETURN
405 END IF
406*
407* Quick return if possible
408*
409 m = 0
410 IF( n.EQ.0 ) THEN
411 RETURN
412 END IF
413*
414* Form a Cholesky factorization of B.
415*
416 CALL zpotrf( uplo, n, b, ldb, info )
417 IF( info.NE.0 ) THEN
418 info = n + info
419 RETURN
420 END IF
421*
422* Transform problem to standard eigenvalue problem and solve.
423*
424 CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
425 CALL zheevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol,
426 $ m, w, z, ldz, work, lwork, rwork, iwork, ifail,
427 $ info )
428*
429 IF( wantz ) THEN
430*
431* Backtransform eigenvectors to the original problem.
432*
433 IF( info.GT.0 )
434 $ m = info - 1
435 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
436*
437* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
438* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
439*
440 IF( upper ) THEN
441 trans = 'N'
442 ELSE
443 trans = 'C'
444 END IF
445*
446 CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, m, cone, b,
447 $ ldb, z, ldz )
448*
449 ELSE IF( itype.EQ.3 ) THEN
450*
451* For B*A*x=(lambda)*x;
452* backtransform eigenvectors: x = L*y or U**H *y
453*
454 IF( upper ) THEN
455 trans = 'C'
456 ELSE
457 trans = 'N'
458 END IF
459*
460 CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, m, cone, b,
461 $ ldb, z, ldz )
462 END IF
463 END IF
464*
465* Set WORK(1) to optimal complex workspace size.
466*
467 work( 1 ) = lwkopt
468*
469 RETURN
470*
471* End of ZHEGVX
472*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:259
subroutine zhegst(itype, uplo, n, a, lda, b, ldb, info)
ZHEGST
Definition zhegst.f:128
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:107
subroutine ztrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRMM
Definition ztrmm.f:177
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180
Here is the call graph for this function:
Here is the caller graph for this function: