LAPACK 3.12.1
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 302 of file zhegvx.f.

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