LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
June 2016
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 302 of file zhbgvx.f.

302 *
303 * -- LAPACK driver routine (version 3.6.1) --
304 * -- LAPACK is a software package provided by Univ. of Tennessee, --
305 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
306 * June 2016
307 *
308 * .. Scalar Arguments ..
309  CHARACTER jobz, range, uplo
310  INTEGER il, info, iu, ka, kb, ldab, ldbb, ldq, ldz, m,
311  $ n
312  DOUBLE PRECISION abstol, vl, vu
313 * ..
314 * .. Array Arguments ..
315  INTEGER ifail( * ), iwork( * )
316  DOUBLE PRECISION rwork( * ), w( * )
317  COMPLEX*16 ab( ldab, * ), bb( ldbb, * ), q( ldq, * ),
318  $ work( * ), z( ldz, * )
319 * ..
320 *
321 * =====================================================================
322 *
323 * .. Parameters ..
324  DOUBLE PRECISION zero
325  parameter ( zero = 0.0d+0 )
326  COMPLEX*16 czero, cone
327  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
328  $ cone = ( 1.0d+0, 0.0d+0 ) )
329 * ..
330 * .. Local Scalars ..
331  LOGICAL alleig, indeig, test, upper, valeig, wantz
332  CHARACTER order, vect
333  INTEGER i, iinfo, indd, inde, indee, indibl, indisp,
334  $ indiwk, indrwk, indwrk, itmp1, j, jj, nsplit
335  DOUBLE PRECISION tmp1
336 * ..
337 * .. External Functions ..
338  LOGICAL lsame
339  EXTERNAL lsame
340 * ..
341 * .. External Subroutines ..
342  EXTERNAL dcopy, dstebz, dsterf, xerbla, zcopy, zgemv,
344  $ zswap
345 * ..
346 * .. Intrinsic Functions ..
347  INTRINSIC min
348 * ..
349 * .. Executable Statements ..
350 *
351 * Test the input parameters.
352 *
353  wantz = lsame( jobz, 'V' )
354  upper = lsame( uplo, 'U' )
355  alleig = lsame( range, 'A' )
356  valeig = lsame( range, 'V' )
357  indeig = lsame( range, 'I' )
358 *
359  info = 0
360  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
361  info = -1
362  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
363  info = -2
364  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
365  info = -3
366  ELSE IF( n.LT.0 ) THEN
367  info = -4
368  ELSE IF( ka.LT.0 ) THEN
369  info = -5
370  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
371  info = -6
372  ELSE IF( ldab.LT.ka+1 ) THEN
373  info = -8
374  ELSE IF( ldbb.LT.kb+1 ) THEN
375  info = -10
376  ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
377  info = -12
378  ELSE
379  IF( valeig ) THEN
380  IF( n.GT.0 .AND. vu.LE.vl )
381  $ info = -14
382  ELSE IF( indeig ) THEN
383  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
384  info = -15
385  ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
386  info = -16
387  END IF
388  END IF
389  END IF
390  IF( info.EQ.0) THEN
391  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
392  info = -21
393  END IF
394  END IF
395 *
396  IF( info.NE.0 ) THEN
397  CALL xerbla( 'ZHBGVX', -info )
398  RETURN
399  END IF
400 *
401 * Quick return if possible
402 *
403  m = 0
404  IF( n.EQ.0 )
405  $ RETURN
406 *
407 * Form a split Cholesky factorization of B.
408 *
409  CALL zpbstf( uplo, n, kb, bb, ldbb, info )
410  IF( info.NE.0 ) THEN
411  info = n + info
412  RETURN
413  END IF
414 *
415 * Transform problem to standard eigenvalue problem.
416 *
417  CALL zhbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
418  $ work, rwork, iinfo )
419 *
420 * Solve the standard eigenvalue problem.
421 * Reduce Hermitian band matrix to tridiagonal form.
422 *
423  indd = 1
424  inde = indd + n
425  indrwk = inde + n
426  indwrk = 1
427  IF( wantz ) THEN
428  vect = 'U'
429  ELSE
430  vect = 'N'
431  END IF
432  CALL zhbtrd( vect, uplo, n, ka, ab, ldab, rwork( indd ),
433  $ rwork( inde ), q, ldq, work( indwrk ), iinfo )
434 *
435 * If all eigenvalues are desired and ABSTOL is less than or equal
436 * to zero, then call DSTERF or ZSTEQR. If this fails for some
437 * eigenvalue, then try DSTEBZ.
438 *
439  test = .false.
440  IF( indeig ) THEN
441  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
442  test = .true.
443  END IF
444  END IF
445  IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
446  CALL dcopy( n, rwork( indd ), 1, w, 1 )
447  indee = indrwk + 2*n
448  CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
449  IF( .NOT.wantz ) THEN
450  CALL dsterf( n, w, rwork( indee ), info )
451  ELSE
452  CALL zlacpy( 'A', n, n, q, ldq, z, ldz )
453  CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
454  $ rwork( indrwk ), info )
455  IF( info.EQ.0 ) THEN
456  DO 10 i = 1, n
457  ifail( i ) = 0
458  10 CONTINUE
459  END IF
460  END IF
461  IF( info.EQ.0 ) THEN
462  m = n
463  GO TO 30
464  END IF
465  info = 0
466  END IF
467 *
468 * Otherwise, call DSTEBZ and, if eigenvectors are desired,
469 * call ZSTEIN.
470 *
471  IF( wantz ) THEN
472  order = 'B'
473  ELSE
474  order = 'E'
475  END IF
476  indibl = 1
477  indisp = indibl + n
478  indiwk = indisp + n
479  CALL dstebz( range, order, n, vl, vu, il, iu, abstol,
480  $ rwork( indd ), rwork( inde ), m, nsplit, w,
481  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
482  $ iwork( indiwk ), info )
483 *
484  IF( wantz ) THEN
485  CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
486  $ iwork( indibl ), iwork( indisp ), z, ldz,
487  $ rwork( indrwk ), iwork( indiwk ), ifail, info )
488 *
489 * Apply unitary matrix used in reduction to tridiagonal
490 * form to eigenvectors returned by ZSTEIN.
491 *
492  DO 20 j = 1, m
493  CALL zcopy( n, z( 1, j ), 1, work( 1 ), 1 )
494  CALL zgemv( 'N', n, n, cone, q, ldq, work, 1, czero,
495  $ z( 1, j ), 1 )
496  20 CONTINUE
497  END IF
498 *
499  30 CONTINUE
500 *
501 * If eigenvalues are not in order, then sort them, along with
502 * eigenvectors.
503 *
504  IF( wantz ) THEN
505  DO 50 j = 1, m - 1
506  i = 0
507  tmp1 = w( j )
508  DO 40 jj = j + 1, m
509  IF( w( jj ).LT.tmp1 ) THEN
510  i = jj
511  tmp1 = w( jj )
512  END IF
513  40 CONTINUE
514 *
515  IF( i.NE.0 ) THEN
516  itmp1 = iwork( indibl+i-1 )
517  w( i ) = w( j )
518  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
519  w( j ) = tmp1
520  iwork( indibl+j-1 ) = itmp1
521  CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
522  IF( info.NE.0 ) THEN
523  itmp1 = ifail( i )
524  ifail( i ) = ifail( j )
525  ifail( j ) = itmp1
526  END IF
527  END IF
528  50 CONTINUE
529  END IF
530 *
531  RETURN
532 *
533 * End of ZHBGVX
534 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine zhbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, RWORK, INFO)
ZHBGST
Definition: zhbgst.f:167
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZSTEQR
Definition: zsteqr.f:134
subroutine zstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
ZSTEIN
Definition: zstein.f:184
subroutine zpbstf(UPLO, N, KD, AB, LDAB, INFO)
ZPBSTF
Definition: zpbstf.f:155
subroutine zhbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
ZHBTRD
Definition: zhbtrd.f:165
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: