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

◆ chbgvx()

subroutine chbgvx ( character jobz,
character range,
character uplo,
integer n,
integer ka,
integer kb,
complex, dimension( ldab, * ) ab,
integer ldab,
complex, dimension( ldbb, * ) bb,
integer ldbb,
complex, dimension( ldq, * ) q,
integer ldq,
real vl,
real vu,
integer il,
integer iu,
real abstol,
integer m,
real, dimension( * ) w,
complex, dimension( ldz, * ) z,
integer ldz,
complex, dimension( * ) work,
real, dimension( * ) rwork,
integer, dimension( * ) iwork,
integer, dimension( * ) ifail,
integer info )

CHBGVX

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

Purpose:
!>
!> CHBGVX 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 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 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 CPBSTF.
!> 
[in]LDBB
!>          LDBB is INTEGER
!>          The leading dimension of the array BB.  LDBB >= KB+1.
!> 
[out]Q
!>          Q is COMPLEX 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 REAL
!>
!>          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 REAL
!>
!>          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 REAL
!>          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*SLAMCH('S'), not zero.
!>          If this routine returns with INFO>0, indicating that some
!>          eigenvectors did not converge, try setting ABSTOL to
!>          2*SLAMCH('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 REAL array, dimension (N)
!>          If INFO = 0, the eigenvalues in ascending order.
!> 
[out]Z
!>          Z is COMPLEX 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 array, dimension (N)
!> 
[out]RWORK
!>          RWORK is REAL 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 CPBSTF
!>                    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 295 of file chbgvx.f.

298*
299* -- LAPACK driver routine --
300* -- LAPACK is a software package provided by Univ. of Tennessee, --
301* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
302*
303* .. Scalar Arguments ..
304 CHARACTER JOBZ, RANGE, UPLO
305 INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
306 $ N
307 REAL ABSTOL, VL, VU
308* ..
309* .. Array Arguments ..
310 INTEGER IFAIL( * ), IWORK( * )
311 REAL RWORK( * ), W( * )
312 COMPLEX AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
313 $ WORK( * ), Z( LDZ, * )
314* ..
315*
316* =====================================================================
317*
318* .. Parameters ..
319 REAL ZERO
320 parameter( zero = 0.0e+0 )
321 COMPLEX CZERO, CONE
322 parameter( czero = ( 0.0e+0, 0.0e+0 ),
323 $ cone = ( 1.0e+0, 0.0e+0 ) )
324* ..
325* .. Local Scalars ..
326 LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
327 CHARACTER ORDER, VECT
328 INTEGER I, IINFO, INDD, INDE, INDEE, INDISP,
329 $ INDIWK, INDRWK, INDWRK, ITMP1, J, JJ, NSPLIT
330 REAL TMP1
331* ..
332* .. External Functions ..
333 LOGICAL LSAME
334 EXTERNAL lsame
335* ..
336* .. External Subroutines ..
337 EXTERNAL ccopy, cgemv, chbgst, chbtrd, clacpy,
338 $ cpbstf,
340 $ xerbla
341* ..
342* .. Intrinsic Functions ..
343 INTRINSIC min
344* ..
345* .. Executable Statements ..
346*
347* Test the input parameters.
348*
349 wantz = lsame( jobz, 'V' )
350 upper = lsame( uplo, 'U' )
351 alleig = lsame( range, 'A' )
352 valeig = lsame( range, 'V' )
353 indeig = lsame( range, 'I' )
354*
355 info = 0
356 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
357 info = -1
358 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
359 info = -2
360 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
361 info = -3
362 ELSE IF( n.LT.0 ) THEN
363 info = -4
364 ELSE IF( ka.LT.0 ) THEN
365 info = -5
366 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
367 info = -6
368 ELSE IF( ldab.LT.ka+1 ) THEN
369 info = -8
370 ELSE IF( ldbb.LT.kb+1 ) THEN
371 info = -10
372 ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
373 info = -12
374 ELSE
375 IF( valeig ) THEN
376 IF( n.GT.0 .AND. vu.LE.vl )
377 $ info = -14
378 ELSE IF( indeig ) THEN
379 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
380 info = -15
381 ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
382 info = -16
383 END IF
384 END IF
385 END IF
386 IF( info.EQ.0) THEN
387 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
388 info = -21
389 END IF
390 END IF
391*
392 IF( info.NE.0 ) THEN
393 CALL xerbla( 'CHBGVX', -info )
394 RETURN
395 END IF
396*
397* Quick return if possible
398*
399 m = 0
400 IF( n.EQ.0 )
401 $ RETURN
402*
403* Form a split Cholesky factorization of B.
404*
405 CALL cpbstf( uplo, n, kb, bb, ldbb, info )
406 IF( info.NE.0 ) THEN
407 info = n + info
408 RETURN
409 END IF
410*
411* Transform problem to standard eigenvalue problem.
412*
413 CALL chbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
414 $ work, rwork, iinfo )
415*
416* Solve the standard eigenvalue problem.
417* Reduce Hermitian band matrix to tridiagonal form.
418*
419 indd = 1
420 inde = indd + n
421 indrwk = inde + n
422 indwrk = 1
423 IF( wantz ) THEN
424 vect = 'U'
425 ELSE
426 vect = 'N'
427 END IF
428 CALL chbtrd( vect, uplo, n, ka, ab, ldab, rwork( indd ),
429 $ rwork( inde ), q, ldq, work( indwrk ), iinfo )
430*
431* If all eigenvalues are desired and ABSTOL is less than or equal
432* to zero, then call SSTERF or CSTEQR. If this fails for some
433* eigenvalue, then try SSTEBZ.
434*
435 test = .false.
436 IF( indeig ) THEN
437 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
438 test = .true.
439 END IF
440 END IF
441 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
442 CALL scopy( n, rwork( indd ), 1, w, 1 )
443 indee = indrwk + 2*n
444 CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
445 IF( .NOT.wantz ) THEN
446 CALL ssterf( n, w, rwork( indee ), info )
447 ELSE
448 CALL clacpy( 'A', n, n, q, ldq, z, ldz )
449 CALL csteqr( jobz, n, w, rwork( indee ), z, ldz,
450 $ rwork( indrwk ), info )
451 IF( info.EQ.0 ) THEN
452 DO 10 i = 1, n
453 ifail( i ) = 0
454 10 CONTINUE
455 END IF
456 END IF
457 IF( info.EQ.0 ) THEN
458 m = n
459 GO TO 30
460 END IF
461 info = 0
462 END IF
463*
464* Otherwise, call SSTEBZ and, if eigenvectors are desired,
465* call CSTEIN.
466*
467 IF( wantz ) THEN
468 order = 'B'
469 ELSE
470 order = 'E'
471 END IF
472 indisp = 1 + n
473 indiwk = indisp + n
474 CALL sstebz( range, order, n, vl, vu, il, iu, abstol,
475 $ rwork( indd ), rwork( inde ), m, nsplit, w,
476 $ iwork( 1 ), iwork( indisp ), rwork( indrwk ),
477 $ iwork( indiwk ), info )
478*
479 IF( wantz ) THEN
480 CALL cstein( n, rwork( indd ), rwork( inde ), m, w,
481 $ iwork( 1 ), iwork( indisp ), z, ldz,
482 $ rwork( indrwk ), iwork( indiwk ), ifail, info )
483*
484* Apply unitary matrix used in reduction to tridiagonal
485* form to eigenvectors returned by CSTEIN.
486*
487 DO 20 j = 1, m
488 CALL ccopy( n, z( 1, j ), 1, work( 1 ), 1 )
489 CALL cgemv( 'N', n, n, cone, q, ldq, work, 1, czero,
490 $ z( 1, j ), 1 )
491 20 CONTINUE
492 END IF
493*
494 30 CONTINUE
495*
496* If eigenvalues are not in order, then sort them, along with
497* eigenvectors.
498*
499 IF( wantz ) THEN
500 DO 50 j = 1, m - 1
501 i = 0
502 tmp1 = w( j )
503 DO 40 jj = j + 1, m
504 IF( w( jj ).LT.tmp1 ) THEN
505 i = jj
506 tmp1 = w( jj )
507 END IF
508 40 CONTINUE
509*
510 IF( i.NE.0 ) THEN
511 itmp1 = iwork( 1 + i-1 )
512 w( i ) = w( j )
513 iwork( 1 + i-1 ) = iwork( 1 + j-1 )
514 w( j ) = tmp1
515 iwork( 1 + j-1 ) = itmp1
516 CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
517 IF( info.NE.0 ) THEN
518 itmp1 = ifail( i )
519 ifail( i ) = ifail( j )
520 ifail( j ) = itmp1
521 END IF
522 END IF
523 50 CONTINUE
524 END IF
525*
526 RETURN
527*
528* End of CHBGVX
529*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
CGEMV
Definition cgemv.f:160
subroutine chbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
CHBGST
Definition chbgst.f:164
subroutine chbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
CHBTRD
Definition chbtrd.f:161
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine cpbstf(uplo, n, kd, ab, ldab, info)
CPBSTF
Definition cpbstf.f:151
subroutine sstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
SSTEBZ
Definition sstebz.f:272
subroutine cstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
CSTEIN
Definition cstein.f:180
subroutine csteqr(compz, n, d, e, z, ldz, work, info)
CSTEQR
Definition csteqr.f:130
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: