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

◆ dsbgvx()

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
                principal minor of order i of B is not positive.
                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 291 of file dsbgvx.f.

294*
295* -- LAPACK driver routine --
296* -- LAPACK is a software package provided by Univ. of Tennessee, --
297* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
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 AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
308 $ W( * ), WORK( * ), Z( LDZ, * )
309* ..
310*
311* =====================================================================
312*
313* .. Parameters ..
314 DOUBLE PRECISION ZERO, ONE
315 parameter( zero = 0.0d+0, one = 1.0d+0 )
316* ..
317* .. Local Scalars ..
318 LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
319 CHARACTER ORDER, VECT
320 INTEGER I, IINFO, INDD, INDE, INDEE, INDISP,
321 $ INDIWO, INDWRK, ITMP1, J, JJ, NSPLIT
322 DOUBLE PRECISION TMP1
323* ..
324* .. External Functions ..
325 LOGICAL LSAME
326 EXTERNAL lsame
327* ..
328* .. External Subroutines ..
329 EXTERNAL dcopy, dgemv, dlacpy, dpbstf, dsbgst, dsbtrd,
331* ..
332* .. Intrinsic Functions ..
333 INTRINSIC min
334* ..
335* .. Executable Statements ..
336*
337* Test the input parameters.
338*
339 wantz = lsame( jobz, 'V' )
340 upper = lsame( uplo, 'U' )
341 alleig = lsame( range, 'A' )
342 valeig = lsame( range, 'V' )
343 indeig = lsame( range, 'I' )
344*
345 info = 0
346 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
347 info = -1
348 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
349 info = -2
350 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
351 info = -3
352 ELSE IF( n.LT.0 ) THEN
353 info = -4
354 ELSE IF( ka.LT.0 ) THEN
355 info = -5
356 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
357 info = -6
358 ELSE IF( ldab.LT.ka+1 ) THEN
359 info = -8
360 ELSE IF( ldbb.LT.kb+1 ) THEN
361 info = -10
362 ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
363 info = -12
364 ELSE
365 IF( valeig ) THEN
366 IF( n.GT.0 .AND. vu.LE.vl )
367 $ info = -14
368 ELSE IF( indeig ) THEN
369 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
370 info = -15
371 ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
372 info = -16
373 END IF
374 END IF
375 END IF
376 IF( info.EQ.0) THEN
377 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
378 info = -21
379 END IF
380 END IF
381*
382 IF( info.NE.0 ) THEN
383 CALL xerbla( 'DSBGVX', -info )
384 RETURN
385 END IF
386*
387* Quick return if possible
388*
389 m = 0
390 IF( n.EQ.0 )
391 $ RETURN
392*
393* Form a split Cholesky factorization of B.
394*
395 CALL dpbstf( uplo, n, kb, bb, ldbb, info )
396 IF( info.NE.0 ) THEN
397 info = n + info
398 RETURN
399 END IF
400*
401* Transform problem to standard eigenvalue problem.
402*
403 CALL dsbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
404 $ work, iinfo )
405*
406* Reduce symmetric band matrix to tridiagonal form.
407*
408 indd = 1
409 inde = indd + n
410 indwrk = inde + n
411 IF( wantz ) THEN
412 vect = 'U'
413 ELSE
414 vect = 'N'
415 END IF
416 CALL dsbtrd( vect, uplo, n, ka, ab, ldab, work( indd ),
417 $ work( inde ), q, ldq, work( indwrk ), iinfo )
418*
419* If all eigenvalues are desired and ABSTOL is less than or equal
420* to zero, then call DSTERF or SSTEQR. If this fails for some
421* eigenvalue, then try DSTEBZ.
422*
423 test = .false.
424 IF( indeig ) THEN
425 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
426 test = .true.
427 END IF
428 END IF
429 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
430 CALL dcopy( n, work( indd ), 1, w, 1 )
431 indee = indwrk + 2*n
432 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
433 IF( .NOT.wantz ) THEN
434 CALL dsterf( n, w, work( indee ), info )
435 ELSE
436 CALL dlacpy( 'A', n, n, q, ldq, z, ldz )
437 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
438 $ work( indwrk ), info )
439 IF( info.EQ.0 ) THEN
440 DO 10 i = 1, n
441 ifail( i ) = 0
442 10 CONTINUE
443 END IF
444 END IF
445 IF( info.EQ.0 ) THEN
446 m = n
447 GO TO 30
448 END IF
449 info = 0
450 END IF
451*
452* Otherwise, call DSTEBZ and, if eigenvectors are desired,
453* call DSTEIN.
454*
455 IF( wantz ) THEN
456 order = 'B'
457 ELSE
458 order = 'E'
459 END IF
460 indisp = 1 + n
461 indiwo = indisp + n
462 CALL dstebz( range, order, n, vl, vu, il, iu, abstol,
463 $ work( indd ), work( inde ), m, nsplit, w,
464 $ iwork( 1 ), iwork( indisp ), work( indwrk ),
465 $ iwork( indiwo ), info )
466*
467 IF( wantz ) THEN
468 CALL dstein( n, work( indd ), work( inde ), m, w,
469 $ iwork( 1 ), iwork( indisp ), z, ldz,
470 $ work( indwrk ), iwork( indiwo ), ifail, info )
471*
472* Apply transformation matrix used in reduction to tridiagonal
473* form to eigenvectors returned by DSTEIN.
474*
475 DO 20 j = 1, m
476 CALL dcopy( n, z( 1, j ), 1, work( 1 ), 1 )
477 CALL dgemv( 'N', n, n, one, q, ldq, work, 1, zero,
478 $ z( 1, j ), 1 )
479 20 CONTINUE
480 END IF
481*
482 30 CONTINUE
483*
484* If eigenvalues are not in order, then sort them, along with
485* eigenvectors.
486*
487 IF( wantz ) THEN
488 DO 50 j = 1, m - 1
489 i = 0
490 tmp1 = w( j )
491 DO 40 jj = j + 1, m
492 IF( w( jj ).LT.tmp1 ) THEN
493 i = jj
494 tmp1 = w( jj )
495 END IF
496 40 CONTINUE
497*
498 IF( i.NE.0 ) THEN
499 itmp1 = iwork( 1 + i-1 )
500 w( i ) = w( j )
501 iwork( 1 + i-1 ) = iwork( 1 + j-1 )
502 w( j ) = tmp1
503 iwork( 1 + j-1 ) = itmp1
504 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
505 IF( info.NE.0 ) THEN
506 itmp1 = ifail( i )
507 ifail( i ) = ifail( j )
508 ifail( j ) = itmp1
509 END IF
510 END IF
511 50 CONTINUE
512 END IF
513*
514 RETURN
515*
516* End of DSBGVX
517*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
subroutine dsbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, info)
DSBGST
Definition dsbgst.f:159
subroutine dsbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
DSBTRD
Definition dsbtrd.f:163
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpbstf(uplo, n, kd, ab, ldab, info)
DPBSTF
Definition dpbstf.f:152
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 dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:174
subroutine dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:131
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
Here is the call graph for this function:
Here is the caller graph for this function: