LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dsbgvx ( character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
integer  KA,
integer  KB,
double precision, dimension( ldab, * )  AB,
integer  LDAB,
double precision, dimension( ldbb, * )  BB,
integer  LDBB,
double precision, dimension( ldq, * )  Q,
integer  LDQ,
double precision  VL,
double precision  VU,
integer  IL,
integer  IU,
double precision  ABSTOL,
integer  M,
double precision, dimension( * )  W,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
double precision, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer, dimension( * )  IFAIL,
integer  INFO 
)

DSBGVX

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

Purpose:
 DSBGVX computes selected eigenvalues, and optionally, eigenvectors
 of a real generalized symmetric-definite banded eigenproblem, of
 the form A*x=(lambda)*B*x.  Here A and B are assumed to be symmetric
 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 DOUBLE PRECISION array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the symmetric 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 DOUBLE PRECISION array, dimension (LDBB, N)
          On entry, the upper or lower triangle of the symmetric 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(ka+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**T*S, as returned by DPBSTF.
[in]LDBB
          LDBB is INTEGER
          The leading dimension of the array BB.  LDBB >= KB+1.
[out]Q
          Q is DOUBLE PRECISION 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 A 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 DOUBLE PRECISION 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 Z**T*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 >= max(1,N).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (7*N)
[out]IWORK
          IWORK is INTEGER array, dimension (5*N)
[out]IFAIL
          IFAIL is INTEGER array, dimension (M)
          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 eigenvalues 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
          <= N: if INFO = i, then i eigenvectors failed to converge.
                  Their indices are stored in IFAIL.
          > N : DPBSTF returned an error code; i.e.,
                if INFO = N + i, for 1 <= i <= N, then the leading
                minor of order i of 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 296 of file dsbgvx.f.

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