LAPACK  3.5.0
LAPACK: Linear Algebra PACKage
zhbgvx.f File Reference

Go to the source code of this file.

Functions/Subroutines

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)
 ZHBGST More...
 

Function/Subroutine Documentation

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 
)

ZHBGST

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
[in]VU
          VU is DOUBLE PRECISION

          If RANGE='V', the lower and upper bounds of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]IL
          IL is INTEGER
[in]IU
          IU is INTEGER

          If RANGE='I', the indices (in ascending order) of the
          smallest and largest eigenvalues 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
November 2011
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 293 of file zhbgvx.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: