297 SUBROUTINE zhbgvx( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
298 $ LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
299 $ LDZ, WORK, RWORK, IWORK, IFAIL, INFO )
306 CHARACTER JOBZ, RANGE, UPLO
307 INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
309 DOUBLE PRECISION ABSTOL, VL, VU
312 INTEGER IFAIL( * ), IWORK( * )
313 DOUBLE PRECISION RWORK( * ), W( * )
314 COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
315 $ work( * ), z( ldz, * )
321 DOUBLE PRECISION ZERO
322 PARAMETER ( ZERO = 0.0d+0 )
323 COMPLEX*16 CZERO, CONE
324 parameter( czero = ( 0.0d+0, 0.0d+0 ),
325 $ cone = ( 1.0d+0, 0.0d+0 ) )
328 LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
329 CHARACTER ORDER, VECT
330 INTEGER I, IINFO, INDD, INDE, INDEE, INDISP,
331 $ indiwk, indrwk, indwrk, itmp1, j, jj, nsplit
332 DOUBLE PRECISION TMP1
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' )
357 IF( .NOT.( wantz .OR. lsame( jobz,
'N' ) ) )
THEN
359 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) )
THEN
361 ELSE IF( .NOT.( upper .OR. lsame( uplo,
'L' ) ) )
THEN
363 ELSE IF( n.LT.0 )
THEN
365 ELSE IF( ka.LT.0 )
THEN
367 ELSE IF( kb.LT.0 .OR. kb.GT.ka )
THEN
369 ELSE IF( ldab.LT.ka+1 )
THEN
371 ELSE IF( ldbb.LT.kb+1 )
THEN
373 ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) )
THEN
377 IF( n.GT.0 .AND. vu.LE.vl )
379 ELSE IF( indeig )
THEN
380 IF( il.LT.1 .OR. il.GT.max( 1, n ) )
THEN
382 ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n )
THEN
388 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
THEN
394 CALL xerbla(
'ZHBGVX', -info )
406 CALL zpbstf( uplo, n, kb, bb, ldbb, info )
414 CALL zhbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
415 $ work, rwork, iinfo )
429 CALL zhbtrd( vect, uplo, n, ka, ab, ldab, rwork( indd ),
430 $ rwork( inde ), q, ldq, work( indwrk ), iinfo )
438 IF( il.EQ.1 .AND. iu.EQ.n )
THEN
442 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) )
THEN
443 CALL dcopy( n, rwork( indd ), 1, w, 1 )
445 CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
446 IF( .NOT.wantz )
THEN
447 CALL dsterf( n, w, rwork( indee ), info )
449 CALL zlacpy(
'A', n, n, q, ldq, z, ldz )
450 CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
451 $ rwork( indrwk ), info )
475 CALL dstebz( range, order, n, vl, vu, il, iu, abstol,
476 $ rwork( indd ), rwork( inde ), m, nsplit, w,
477 $ iwork( 1 ), iwork( indisp ), rwork( indrwk ),
478 $ iwork( indiwk ), info )
481 CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
482 $ iwork( 1 ), iwork( indisp ), z, ldz,
483 $ rwork( indrwk ), iwork( indiwk ), ifail, info )
489 CALL zcopy( n, z( 1, j ), 1, work( 1 ), 1 )
490 CALL zgemv(
'N', n, n, cone, q, ldq, work, 1, czero,
505 IF( w( jj ).LT.tmp1 )
THEN
512 itmp1 = iwork( 1 + i-1 )
514 iwork( 1 + i-1 ) = iwork( 1 + j-1 )
516 iwork( 1 + j-1 ) = itmp1
517 CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
520 ifail( i ) = ifail( j )
subroutine xerbla(srname, info)
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
subroutine zhbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
ZHBGST
subroutine zhbgvx(jobz, range, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
ZHBGVX
subroutine zhbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
ZHBTRD
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
subroutine zpbstf(uplo, n, kd, ab, ldab, info)
ZPBSTF
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 zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
subroutine dsterf(n, d, e, info)
DSTERF
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP