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

◆ chbevx_2stage()

subroutine chbevx_2stage ( character jobz,
character range,
character uplo,
integer n,
integer kd,
complex, dimension( ldab, * ) ab,
integer ldab,
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,
integer lwork,
real, dimension( * ) rwork,
integer, dimension( * ) iwork,
integer, dimension( * ) ifail,
integer info )

CHBEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

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

Purpose:
!>
!> CHBEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
!> of a complex Hermitian band matrix A using the 2stage technique for
!> the reduction to tridiagonal.  Eigenvalues and eigenvectors
!> can be selected by specifying either 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.
!>                  Not available in this release.
!> 
[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 triangle of A is stored;
!>          = 'L':  Lower triangle of A is stored.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in]KD
!>          KD is INTEGER
!>          The number of superdiagonals of the matrix A if UPLO = 'U',
!>          or the number of subdiagonals if UPLO = 'L'.  KD >= 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 KD+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(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
!>          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).
!>
!>          On exit, AB is overwritten by values generated during the
!>          reduction to tridiagonal form.
!> 
[in]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= KD + 1.
!> 
[out]Q
!>          Q is COMPLEX array, dimension (LDQ, N)
!>          If JOBZ = 'V', the N-by-N unitary matrix used in the
!>                          reduction 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 = 'V', then
!>          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 AB 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').
!>
!>          See  by Demmel and
!>          Kahan, LAPACK Working Note #3.
!> 
[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)
!>          The first M elements contain the selected eigenvalues in
!>          ascending order.
!> 
[out]Z
!>          Z is COMPLEX array, dimension (LDZ, max(1,M))
!>          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
!>          contain the orthonormal eigenvectors of the matrix A
!>          corresponding to the selected eigenvalues, with the i-th
!>          column of Z holding the eigenvector associated with W(i).
!>          If an eigenvector fails to converge, then that column of Z
!>          contains the latest approximation to the eigenvector, and the
!>          index of the eigenvector is returned in IFAIL.
!>          If JOBZ = 'N', then Z is not referenced.
!>          Note: the user must ensure that at least max(1,M) columns are
!>          supplied in the array Z; if RANGE = 'V', the exact value of M
!>          is not known in advance and an upper bound must be used.
!> 
[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 COMPLEX array, dimension (LWORK)
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of the array WORK. LWORK >= 1, when N <= 1;
!>          otherwise
!>          If JOBZ = 'N' and N > 1, LWORK must be queried.
!>                                   LWORK = MAX(1, dimension) where
!>                                   dimension = (2KD+1)*N + KD*NTHREADS
!>                                   where KD is the size of the band.
!>                                   NTHREADS is the number of threads used when
!>                                   openMP compilation is enabled, otherwise =1.
!>          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal sizes of the WORK, RWORK and
!>          IWORK arrays, returns these values as the first entries of
!>          the WORK, RWORK and IWORK arrays, and no error message
!>          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
!> 
[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, then i eigenvectors failed to converge.
!>                Their indices are stored in array IFAIL.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  All details about the 2stage techniques are available in:
!>
!>  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
!>  Parallel reduction to condensed forms for symmetric eigenvalue problems
!>  using aggregated fine-grained and memory-aware kernels. In Proceedings
!>  of 2011 International Conference for High Performance Computing,
!>  Networking, Storage and Analysis (SC '11), New York, NY, USA,
!>  Article 8 , 11 pages.
!>  http://doi.acm.org/10.1145/2063384.2063394
!>
!>  A. Haidar, J. Kurzak, P. Luszczek, 2013.
!>  An improved parallel singular value algorithm and its implementation
!>  for multicore hardware, In Proceedings of 2013 International Conference
!>  for High Performance Computing, Networking, Storage and Analysis (SC '13).
!>  Denver, Colorado, USA, 2013.
!>  Article 90, 12 pages.
!>  http://doi.acm.org/10.1145/2503210.2503292
!>
!>  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
!>  A novel hybrid CPU-GPU generalized eigensolver for electronic structure
!>  calculations based on fine-grained memory aware tasks.
!>  International Journal of High Performance Computing Applications.
!>  Volume 28 Issue 2, Pages 196-209, May 2014.
!>  http://hpc.sagepub.com/content/28/2/196
!>
!> 

Definition at line 321 of file chbevx_2stage.f.

325*
326 IMPLICIT NONE
327*
328* -- LAPACK driver routine --
329* -- LAPACK is a software package provided by Univ. of Tennessee, --
330* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
331*
332* .. Scalar Arguments ..
333 CHARACTER JOBZ, RANGE, UPLO
334 INTEGER IL, INFO, IU, KD, LDAB, LDQ, LDZ, M, N, LWORK
335 REAL ABSTOL, VL, VU
336* ..
337* .. Array Arguments ..
338 INTEGER IFAIL( * ), IWORK( * )
339 REAL RWORK( * ), W( * )
340 COMPLEX AB( LDAB, * ), Q( LDQ, * ), WORK( * ),
341 $ Z( LDZ, * )
342* ..
343*
344* =====================================================================
345*
346* .. Parameters ..
347 REAL ZERO, ONE
348 parameter( zero = 0.0e0, one = 1.0e0 )
349 COMPLEX CZERO, CONE
350 parameter( czero = ( 0.0e0, 0.0e0 ),
351 $ cone = ( 1.0e0, 0.0e0 ) )
352* ..
353* .. Local Scalars ..
354 LOGICAL ALLEIG, INDEIG, LOWER, TEST, VALEIG, WANTZ,
355 $ LQUERY
356 CHARACTER ORDER
357 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
358 $ INDISP, INDIWK, INDRWK, INDWRK, ISCALE, ITMP1,
359 $ LLWORK, LWMIN, LHTRD, LWTRD, IB, INDHOUS,
360 $ J, JJ, NSPLIT
361 REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
362 $ SIGMA, SMLNUM, TMP1, VLL, VUU
363 COMPLEX CTMP1
364* ..
365* .. External Functions ..
366 LOGICAL LSAME
367 INTEGER ILAENV2STAGE
368 REAL SLAMCH, CLANHB, SROUNDUP_LWORK
369 EXTERNAL lsame, slamch, clanhb, ilaenv2stage,
371* ..
372* .. External Subroutines ..
373 EXTERNAL scopy, sscal, sstebz, ssterf, xerbla,
374 $ ccopy,
377* ..
378* .. Intrinsic Functions ..
379 INTRINSIC real, max, min, sqrt
380* ..
381* .. Executable Statements ..
382*
383* Test the input parameters.
384*
385 wantz = lsame( jobz, 'V' )
386 alleig = lsame( range, 'A' )
387 valeig = lsame( range, 'V' )
388 indeig = lsame( range, 'I' )
389 lower = lsame( uplo, 'L' )
390 lquery = ( lwork.EQ.-1 )
391*
392 info = 0
393 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
394 info = -1
395 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
396 info = -2
397 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
398 info = -3
399 ELSE IF( n.LT.0 ) THEN
400 info = -4
401 ELSE IF( kd.LT.0 ) THEN
402 info = -5
403 ELSE IF( ldab.LT.kd+1 ) THEN
404 info = -7
405 ELSE IF( wantz .AND. ldq.LT.max( 1, n ) ) THEN
406 info = -9
407 ELSE
408 IF( valeig ) THEN
409 IF( n.GT.0 .AND. vu.LE.vl )
410 $ info = -11
411 ELSE IF( indeig ) THEN
412 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
413 info = -12
414 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
415 info = -13
416 END IF
417 END IF
418 END IF
419 IF( info.EQ.0 ) THEN
420 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) )
421 $ info = -18
422 END IF
423*
424 IF( info.EQ.0 ) THEN
425 IF( n.LE.1 ) THEN
426 lwmin = 1
427 work( 1 ) = sroundup_lwork(lwmin)
428 ELSE
429 ib = ilaenv2stage( 2, 'CHETRD_HB2ST', jobz,
430 $ n, kd, -1, -1 )
431 lhtrd = ilaenv2stage( 3, 'CHETRD_HB2ST', jobz,
432 $ n, kd, ib, -1 )
433 lwtrd = ilaenv2stage( 4, 'CHETRD_HB2ST', jobz,
434 $ n, kd, ib, -1 )
435 lwmin = lhtrd + lwtrd
436 work( 1 ) = sroundup_lwork(lwmin)
437 ENDIF
438*
439 IF( lwork.LT.lwmin .AND. .NOT.lquery )
440 $ info = -20
441 END IF
442*
443 IF( info.NE.0 ) THEN
444 CALL xerbla( 'CHBEVX_2STAGE', -info )
445 RETURN
446 ELSE IF( lquery ) THEN
447 RETURN
448 END IF
449*
450* Quick return if possible
451*
452 m = 0
453 IF( n.EQ.0 )
454 $ RETURN
455*
456 IF( n.EQ.1 ) THEN
457 m = 1
458 IF( lower ) THEN
459 ctmp1 = ab( 1, 1 )
460 ELSE
461 ctmp1 = ab( kd+1, 1 )
462 END IF
463 tmp1 = real( ctmp1 )
464 IF( valeig ) THEN
465 IF( .NOT.( vl.LT.tmp1 .AND. vu.GE.tmp1 ) )
466 $ m = 0
467 END IF
468 IF( m.EQ.1 ) THEN
469 w( 1 ) = real( ctmp1 )
470 IF( wantz )
471 $ z( 1, 1 ) = cone
472 END IF
473 RETURN
474 END IF
475*
476* Get machine constants.
477*
478 safmin = slamch( 'Safe minimum' )
479 eps = slamch( 'Precision' )
480 smlnum = safmin / eps
481 bignum = one / smlnum
482 rmin = sqrt( smlnum )
483 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
484*
485* Scale matrix to allowable range, if necessary.
486*
487 iscale = 0
488 abstll = abstol
489 IF( valeig ) THEN
490 vll = vl
491 vuu = vu
492 ELSE
493 vll = zero
494 vuu = zero
495 END IF
496 anrm = clanhb( 'M', uplo, n, kd, ab, ldab, rwork )
497 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
498 iscale = 1
499 sigma = rmin / anrm
500 ELSE IF( anrm.GT.rmax ) THEN
501 iscale = 1
502 sigma = rmax / anrm
503 END IF
504 IF( iscale.EQ.1 ) THEN
505 IF( lower ) THEN
506 CALL clascl( 'B', kd, kd, one, sigma, n, n, ab, ldab,
507 $ info )
508 ELSE
509 CALL clascl( 'Q', kd, kd, one, sigma, n, n, ab, ldab,
510 $ info )
511 END IF
512 IF( abstol.GT.0 )
513 $ abstll = abstol*sigma
514 IF( valeig ) THEN
515 vll = vl*sigma
516 vuu = vu*sigma
517 END IF
518 END IF
519*
520* Call CHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form.
521*
522 indd = 1
523 inde = indd + n
524 indrwk = inde + n
525*
526 indhous = 1
527 indwrk = indhous + lhtrd
528 llwork = lwork - indwrk + 1
529*
530 CALL chetrd_hb2st( 'N', jobz, uplo, n, kd, ab, ldab,
531 $ rwork( indd ), rwork( inde ), work( indhous ),
532 $ lhtrd, work( indwrk ), llwork, iinfo )
533*
534* If all eigenvalues are desired and ABSTOL is less than or equal
535* to zero, then call SSTERF or CSTEQR. If this fails for some
536* eigenvalue, then try SSTEBZ.
537*
538 test = .false.
539 IF (indeig) THEN
540 IF (il.EQ.1 .AND. iu.EQ.n) THEN
541 test = .true.
542 END IF
543 END IF
544 IF ((alleig .OR. test) .AND. (abstol.LE.zero)) THEN
545 CALL scopy( n, rwork( indd ), 1, w, 1 )
546 indee = indrwk + 2*n
547 IF( .NOT.wantz ) THEN
548 CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
549 CALL ssterf( n, w, rwork( indee ), info )
550 ELSE
551 CALL clacpy( 'A', n, n, q, ldq, z, ldz )
552 CALL scopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
553 CALL csteqr( jobz, n, w, rwork( indee ), z, ldz,
554 $ rwork( indrwk ), info )
555 IF( info.EQ.0 ) THEN
556 DO 10 i = 1, n
557 ifail( i ) = 0
558 10 CONTINUE
559 END IF
560 END IF
561 IF( info.EQ.0 ) THEN
562 m = n
563 GO TO 30
564 END IF
565 info = 0
566 END IF
567*
568* Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
569*
570 IF( wantz ) THEN
571 order = 'B'
572 ELSE
573 order = 'E'
574 END IF
575 indibl = 1
576 indisp = indibl + n
577 indiwk = indisp + n
578 CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
579 $ rwork( indd ), rwork( inde ), m, nsplit, w,
580 $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
581 $ iwork( indiwk ), info )
582*
583 IF( wantz ) THEN
584 CALL cstein( n, rwork( indd ), rwork( inde ), m, w,
585 $ iwork( indibl ), iwork( indisp ), z, ldz,
586 $ rwork( indrwk ), iwork( indiwk ), ifail, info )
587*
588* Apply unitary matrix used in reduction to tridiagonal
589* form to eigenvectors returned by CSTEIN.
590*
591 DO 20 j = 1, m
592 CALL ccopy( n, z( 1, j ), 1, work( 1 ), 1 )
593 CALL cgemv( 'N', n, n, cone, q, ldq, work, 1, czero,
594 $ z( 1, j ), 1 )
595 20 CONTINUE
596 END IF
597*
598* If matrix was scaled, then rescale eigenvalues appropriately.
599*
600 30 CONTINUE
601 IF( iscale.EQ.1 ) THEN
602 IF( info.EQ.0 ) THEN
603 imax = m
604 ELSE
605 imax = info - 1
606 END IF
607 CALL sscal( imax, one / sigma, w, 1 )
608 END IF
609*
610* If eigenvalues are not in order, then sort them, along with
611* eigenvectors.
612*
613 IF( wantz ) THEN
614 DO 50 j = 1, m - 1
615 i = 0
616 tmp1 = w( j )
617 DO 40 jj = j + 1, m
618 IF( w( jj ).LT.tmp1 ) THEN
619 i = jj
620 tmp1 = w( jj )
621 END IF
622 40 CONTINUE
623*
624 IF( i.NE.0 ) THEN
625 itmp1 = iwork( indibl+i-1 )
626 w( i ) = w( j )
627 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
628 w( j ) = tmp1
629 iwork( indibl+j-1 ) = itmp1
630 CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
631 IF( info.NE.0 ) THEN
632 itmp1 = ifail( i )
633 ifail( i ) = ifail( j )
634 ifail( j ) = itmp1
635 END IF
636 END IF
637 50 CONTINUE
638 END IF
639*
640* Set WORK(1) to optimal workspace size.
641*
642 work( 1 ) = sroundup_lwork(lwmin)
643*
644 RETURN
645*
646* End of CHBEVX_2STAGE
647*
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 chetrd_hb2st(stage1, vect, uplo, n, kd, ab, ldab, d, e, hous, lhous, work, lwork, info)
CHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
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
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clanhb(norm, uplo, n, k, ab, ldab, work)
CLANHB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clanhb.f:130
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition clascl.f:142
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
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: