305 SUBROUTINE zheevx_2stage( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
306 $ IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
307 $ LWORK, RWORK, IWORK, IFAIL, INFO )
317 CHARACTER JOBZ, RANGE, UPLO
318 INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
319 DOUBLE PRECISION ABSTOL, VL, VU
322 INTEGER IFAIL( * ), IWORK( * )
323 DOUBLE PRECISION RWORK( * ), W( * )
324 COMPLEX*16 A( lda, * ), WORK( * ), Z( ldz, * )
330 DOUBLE PRECISION ZERO, ONE
331 parameter( zero = 0.0d+0, one = 1.0d+0 )
333 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
336 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
339 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
340 $ indisp, indiwk, indrwk, indtau, indwrk, iscale,
341 $ itmp1, j, jj, llwork,
342 $ nsplit, lwmin, lhtrd, lwtrd, kd, ib, indhous
343 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
344 $ sigma, smlnum, tmp1, vll, vuu
349 DOUBLE PRECISION DLAMCH, ZLANHE
350 EXTERNAL lsame, ilaenv, dlamch, zlanhe
358 INTRINSIC dble, max, min, sqrt
364 lower = lsame( uplo,
'L' )
365 wantz = lsame( jobz,
'V' )
366 alleig = lsame( range,
'A' )
367 valeig = lsame( range,
'V' )
368 indeig = lsame( range,
'I' )
369 lquery = ( lwork.EQ.-1 )
372 IF( .NOT.( lsame( jobz,
'N' ) ) )
THEN 374 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) )
THEN 376 ELSE IF( .NOT.( lower .OR. lsame( uplo,
'U' ) ) )
THEN 378 ELSE IF( n.LT.0 )
THEN 380 ELSE IF( lda.LT.max( 1, n ) )
THEN 384 IF( n.GT.0 .AND. vu.LE.vl )
386 ELSE IF( indeig )
THEN 387 IF( il.LT.1 .OR. il.GT.max( 1, n ) )
THEN 389 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n )
THEN 395 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
THEN 405 kd = ilaenv( 17,
'ZHETRD_2STAGE', jobz, n, -1, -1, -1 )
406 ib = ilaenv( 18,
'ZHETRD_2STAGE', jobz, n, kd, -1, -1 )
407 lhtrd = ilaenv( 19,
'ZHETRD_2STAGE', jobz, n, kd, ib, -1 )
408 lwtrd = ilaenv( 20,
'ZHETRD_2STAGE', jobz, n, kd, ib, -1 )
409 lwmin = n + lhtrd + lwtrd
413 IF( lwork.LT.lwmin .AND. .NOT.lquery )
418 CALL xerbla(
'ZHEEVX_2STAGE', -info )
420 ELSE IF( lquery )
THEN 432 IF( alleig .OR. indeig )
THEN 434 w( 1 ) = dble( a( 1, 1 ) )
435 ELSE IF( valeig )
THEN 436 IF( vl.LT.dble( a( 1, 1 ) ) .AND. vu.GE.dble( a( 1, 1 ) ) )
439 w( 1 ) = dble( a( 1, 1 ) )
449 safmin = dlamch(
'Safe minimum' )
450 eps = dlamch(
'Precision' )
451 smlnum = safmin / eps
452 bignum = one / smlnum
453 rmin = sqrt( smlnum )
454 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
464 anrm = zlanhe(
'M', uplo, n, a, lda, rwork )
465 IF( anrm.GT.zero .AND. anrm.LT.rmin )
THEN 468 ELSE IF( anrm.GT.rmax )
THEN 472 IF( iscale.EQ.1 )
THEN 475 CALL zdscal( n-j+1, sigma, a( j, j ), 1 )
479 CALL zdscal( j, sigma, a( 1, j ), 1 )
483 $ abstll = abstol*sigma
497 indwrk = indhous + lhtrd
498 llwork = lwork - indwrk + 1
501 $ rwork( inde ), work( indtau ),
502 $ work( indhous ), lhtrd, work( indwrk ),
511 IF( il.EQ.1 .AND. iu.EQ.n )
THEN 515 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) )
THEN 516 CALL dcopy( n, rwork( indd ), 1, w, 1 )
518 IF( .NOT.wantz )
THEN 519 CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
520 CALL dsterf( n, w, rwork( indee ), info )
522 CALL zlacpy(
'A', n, n, a, lda, z, ldz )
523 CALL zungtr( uplo, n, z, ldz, work( indtau ),
524 $ work( indwrk ), llwork, iinfo )
525 CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
526 CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
527 $ rwork( indrwk ), info )
551 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
552 $ rwork( indd ), rwork( inde ), m, nsplit, w,
553 $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
554 $ iwork( indiwk ), info )
557 CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
558 $ iwork( indibl ), iwork( indisp ), z, ldz,
559 $ rwork( indrwk ), iwork( indiwk ), ifail, info )
564 CALL zunmtr(
'L', uplo,
'N', n, m, a, lda, work( indtau ), z,
565 $ ldz, work( indwrk ), llwork, iinfo )
571 IF( iscale.EQ.1 )
THEN 577 CALL dscal( imax, one / sigma, w, 1 )
588 IF( w( jj ).LT.tmp1 )
THEN 595 itmp1 = iwork( indibl+i-1 )
597 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
599 iwork( indibl+j-1 ) = itmp1
600 CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
603 ifail( i ) = ifail( j )
subroutine zheevx_2stage(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, IFAIL, INFO)
ZHEEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE ma...
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
subroutine zungtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGTR
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
subroutine dsterf(N, D, E, INFO)
DSTERF
subroutine zsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZSTEQR
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine dscal(N, DA, DX, INCX)
DSCAL
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
subroutine zstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
ZSTEIN
subroutine zunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMTR
subroutine zhetrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
ZHETRD_2STAGE