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

◆ zhbgvx()

subroutine zhbgvx ( character  jobz,
character  range,
character  uplo,
integer  n,
integer  ka,
integer  kb,
complex*16, dimension( ldab, * )  ab,
integer  ldab,
complex*16, dimension( ldbb, * )  bb,
integer  ldbb,
complex*16, dimension( ldq, * )  q,
integer  ldq,
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,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer, dimension( * )  ifail,
integer  info 
)

ZHBGVX

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

Purpose:
 ZHBGVX computes all the eigenvalues, and optionally, the eigenvectors
 of a complex generalized Hermitian-definite banded eigenproblem, of
 the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
 and banded, and B is also positive definite.  Eigenvalues and
 eigenvectors can be selected by specifying either all eigenvalues,
 a range of values or a range of indices for the desired eigenvalues.
Parameters
[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]KA
          KA is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'. KA >= 0.
[in]KB
          KB is INTEGER
          The number of superdiagonals of the matrix B if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'. KB >= 0.
[in,out]AB
          AB is COMPLEX*16 array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix A, stored in the first ka+1 rows of the array.  The
          j-th column of A is stored in the j-th column of the array AB
          as follows:
          if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+ka).

          On exit, the contents of AB are destroyed.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KA+1.
[in,out]BB
          BB is COMPLEX*16 array, dimension (LDBB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix B, stored in the first kb+1 rows of the array.  The
          j-th column of B is stored in the j-th column of the array BB
          as follows:
          if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
          if UPLO = 'L', BB(1+i-j,j)    = B(i,j) for j<=i<=min(n,j+kb).

          On exit, the factor S from the split Cholesky factorization
          B = S**H*S, as returned by ZPBSTF.
[in]LDBB
          LDBB is INTEGER
          The leading dimension of the array BB.  LDBB >= KB+1.
[out]Q
          Q is COMPLEX*16 array, dimension (LDQ, N)
          If JOBZ = 'V', the n-by-n matrix used in the reduction of
          A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
          and consequently C to tridiagonal form.
          If JOBZ = 'N', the array Q is not referenced.
[in]LDQ
          LDQ is INTEGER
          The leading dimension of the array Q.  If JOBZ = 'N',
          LDQ >= 1. If JOBZ = 'V', LDQ >= 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 AP to tridiagonal form.

          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)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
          eigenvectors, with the i-th column of Z holding the
          eigenvector associated with W(i). The eigenvectors are
          normalized so that Z**H*B*Z = I.
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= N.
[out]WORK
          WORK is COMPLEX*16 array, dimension (N)
[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:  if INFO = i, and i is:
             <= N:  then i eigenvectors failed to converge.  Their
                    indices are stored in array IFAIL.
             > N:   if INFO = N + i, for 1 <= i <= N, then ZPBSTF
                    returned INFO = i: 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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 297 of file zhbgvx.f.

300*
301* -- LAPACK driver routine --
302* -- LAPACK is a software package provided by Univ. of Tennessee, --
303* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
304*
305* .. Scalar Arguments ..
306 CHARACTER JOBZ, RANGE, UPLO
307 INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
308 $ N
309 DOUBLE PRECISION ABSTOL, VL, VU
310* ..
311* .. Array Arguments ..
312 INTEGER IFAIL( * ), IWORK( * )
313 DOUBLE PRECISION RWORK( * ), W( * )
314 COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
315 $ WORK( * ), Z( LDZ, * )
316* ..
317*
318* =====================================================================
319*
320* .. Parameters ..
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 ) )
326* ..
327* .. Local Scalars ..
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
333* ..
334* .. External Functions ..
335 LOGICAL LSAME
336 EXTERNAL lsame
337* ..
338* .. External Subroutines ..
339 EXTERNAL dcopy, dstebz, dsterf, xerbla, zcopy, zgemv,
341 $ zswap
342* ..
343* .. Intrinsic Functions ..
344 INTRINSIC 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*
356 info = 0
357 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
358 info = -1
359 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
360 info = -2
361 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
362 info = -3
363 ELSE IF( n.LT.0 ) THEN
364 info = -4
365 ELSE IF( ka.LT.0 ) THEN
366 info = -5
367 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
368 info = -6
369 ELSE IF( ldab.LT.ka+1 ) THEN
370 info = -8
371 ELSE IF( ldbb.LT.kb+1 ) THEN
372 info = -10
373 ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
374 info = -12
375 ELSE
376 IF( valeig ) THEN
377 IF( n.GT.0 .AND. vu.LE.vl )
378 $ info = -14
379 ELSE IF( indeig ) THEN
380 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
381 info = -15
382 ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
383 info = -16
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 = -21
390 END IF
391 END IF
392*
393 IF( info.NE.0 ) THEN
394 CALL xerbla( 'ZHBGVX', -info )
395 RETURN
396 END IF
397*
398* Quick return if possible
399*
400 m = 0
401 IF( n.EQ.0 )
402 $ RETURN
403*
404* Form a split Cholesky factorization of B.
405*
406 CALL zpbstf( uplo, n, kb, bb, ldbb, info )
407 IF( info.NE.0 ) THEN
408 info = n + info
409 RETURN
410 END IF
411*
412* Transform problem to standard eigenvalue problem.
413*
414 CALL zhbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
415 $ work, rwork, iinfo )
416*
417* Solve the standard eigenvalue problem.
418* Reduce Hermitian band matrix to tridiagonal form.
419*
420 indd = 1
421 inde = indd + n
422 indrwk = inde + n
423 indwrk = 1
424 IF( wantz ) THEN
425 vect = 'U'
426 ELSE
427 vect = 'N'
428 END IF
429 CALL zhbtrd( vect, uplo, n, ka, ab, ldab, rwork( indd ),
430 $ rwork( inde ), q, ldq, work( indwrk ), iinfo )
431*
432* If all eigenvalues are desired and ABSTOL is less than or equal
433* to zero, then call DSTERF or ZSTEQR. If this fails for some
434* eigenvalue, then try DSTEBZ.
435*
436 test = .false.
437 IF( indeig ) THEN
438 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
439 test = .true.
440 END IF
441 END IF
442 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
443 CALL dcopy( n, rwork( indd ), 1, w, 1 )
444 indee = indrwk + 2*n
445 CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
446 IF( .NOT.wantz ) THEN
447 CALL dsterf( n, w, rwork( indee ), info )
448 ELSE
449 CALL zlacpy( 'A', n, n, q, ldq, z, ldz )
450 CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
451 $ rwork( indrwk ), info )
452 IF( info.EQ.0 ) THEN
453 DO 10 i = 1, n
454 ifail( i ) = 0
455 10 CONTINUE
456 END IF
457 END IF
458 IF( info.EQ.0 ) THEN
459 m = n
460 GO TO 30
461 END IF
462 info = 0
463 END IF
464*
465* Otherwise, call DSTEBZ and, if eigenvectors are desired,
466* call ZSTEIN.
467*
468 IF( wantz ) THEN
469 order = 'B'
470 ELSE
471 order = 'E'
472 END IF
473 indisp = 1 + n
474 indiwk = indisp + n
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 )
479*
480 IF( wantz ) THEN
481 CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
482 $ iwork( 1 ), iwork( indisp ), z, ldz,
483 $ rwork( indrwk ), iwork( indiwk ), ifail, info )
484*
485* Apply unitary matrix used in reduction to tridiagonal
486* form to eigenvectors returned by ZSTEIN.
487*
488 DO 20 j = 1, m
489 CALL zcopy( n, z( 1, j ), 1, work( 1 ), 1 )
490 CALL zgemv( 'N', n, n, cone, q, ldq, work, 1, czero,
491 $ z( 1, j ), 1 )
492 20 CONTINUE
493 END IF
494*
495 30 CONTINUE
496*
497* If eigenvalues are not in order, then sort them, along with
498* eigenvectors.
499*
500 IF( wantz ) THEN
501 DO 50 j = 1, m - 1
502 i = 0
503 tmp1 = w( j )
504 DO 40 jj = j + 1, m
505 IF( w( jj ).LT.tmp1 ) THEN
506 i = jj
507 tmp1 = w( jj )
508 END IF
509 40 CONTINUE
510*
511 IF( i.NE.0 ) THEN
512 itmp1 = iwork( 1 + i-1 )
513 w( i ) = w( j )
514 iwork( 1 + i-1 ) = iwork( 1 + j-1 )
515 w( j ) = tmp1
516 iwork( 1 + j-1 ) = itmp1
517 CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
518 IF( info.NE.0 ) THEN
519 itmp1 = ifail( i )
520 ifail( i ) = ifail( j )
521 ifail( j ) = itmp1
522 END IF
523 END IF
524 50 CONTINUE
525 END IF
526*
527 RETURN
528*
529* End of ZHBGVX
530*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:160
subroutine zhbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
ZHBGST
Definition zhbgst.f:165
subroutine zhbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
ZHBTRD
Definition zhbtrd.f:163
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpbstf(uplo, n, kd, ab, ldab, info)
ZPBSTF
Definition zpbstf.f:153
subroutine dstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
DSTEBZ
Definition dstebz.f:273
subroutine zstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
ZSTEIN
Definition zstein.f:182
subroutine zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
Definition zsteqr.f:132
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: