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

◆ cheevr()

subroutine cheevr ( character jobz,
character range,
character uplo,
integer n,
complex, dimension( lda, * ) a,
integer lda,
real vl,
real vu,
integer il,
integer iu,
real abstol,
integer m,
real, dimension( * ) w,
complex, dimension( ldz, * ) z,
integer ldz,
integer, dimension( * ) isuppz,
complex, dimension( * ) work,
integer lwork,
real, dimension( * ) rwork,
integer lrwork,
integer, dimension( * ) iwork,
integer liwork,
integer info )

CHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices

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

Purpose:
!>
!> CHEEVR computes selected eigenvalues and, optionally, eigenvectors
!> of a complex Hermitian matrix A. Eigenvalues and eigenvectors can be
!> selected by specifying either a range of values or a range of indices
!> for the desired eigenvalues. Invocations with different choices for
!> these parameters may result in the computation of slightly different
!> eigenvalues and/or eigenvectors for the same matrix. The reason for
!> this behavior is that there exists a variety of algorithms (each
!> performing best for a particular set of options) with CHEEVR
!> attempting to select the best based on the various parameters. In all
!> cases, the computed values are accurate within the limits of finite
!> precision arithmetic.
!>
!> CHEEVR first reduces the matrix A to tridiagonal form T with a call
!> to CHETRD.  Then, whenever possible, CHEEVR calls CSTEMR to compute
!> the eigenspectrum using Relatively Robust Representations.  CSTEMR
!> computes eigenvalues by the dqds algorithm, while orthogonal
!> eigenvectors are computed from various  L D L^T representations
!> (also known as Relatively Robust Representations). Gram-Schmidt
!> orthogonalization is avoided as far as possible. More specifically,
!> the various steps of the algorithm are as follows.
!>
!> For each unreduced block (submatrix) of T,
!>    (a) Compute T - sigma I  = L D L^T, so that L and D
!>        define all the wanted eigenvalues to high relative accuracy.
!>        This means that small relative changes in the entries of D and L
!>        cause only small relative changes in the eigenvalues and
!>        eigenvectors. The standard (unfactored) representation of the
!>        tridiagonal matrix T does not have this property in general.
!>    (b) Compute the eigenvalues to suitable accuracy.
!>        If the eigenvectors are desired, the algorithm attains full
!>        accuracy of the computed eigenvalues only right before
!>        the corresponding vectors have to be computed, see steps c) and d).
!>    (c) For each cluster of close eigenvalues, select a new
!>        shift close to the cluster, find a new factorization, and refine
!>        the shifted eigenvalues to suitable accuracy.
!>    (d) For each eigenvalue with a large enough relative separation compute
!>        the corresponding eigenvector by forming a rank revealing twisted
!>        factorization. Go back to (c) for any clusters that remain.
!>
!> The desired accuracy of the output can be specified by the input
!> parameter ABSTOL.
!>
!> For more details, see CSTEMR's documentation and:
!> - Inderjit S. Dhillon and Beresford N. Parlett: 
!>   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
!> - Inderjit Dhillon and Beresford Parlett:  SIAM Journal on Matrix Analysis and Applications, Vol. 25,
!>   2004.  Also LAPACK Working Note 154.
!> - Inderjit Dhillon: ,
!>   Computer Science Division Technical Report No. UCB/CSD-97-971,
!>   UC Berkeley, May 1997.
!>
!>
!> Note 1 : CHEEVR calls CSTEMR when the full spectrum is requested
!> on machines which conform to the ieee-754 floating point standard.
!> CHEEVR calls SSTEBZ and CSTEIN on non-ieee machines and
!> when partial spectrum requests are made.
!>
!> Normal execution of CSTEMR may create NaNs and infinities and
!> hence may abort due to a floating point exception in environments
!> which do not handle NaNs and infinities in the ieee standard default
!> manner.
!> 
Parameters
[in]JOBZ
!>          JOBZ is CHARACTER*1
!>          = 'N':  Compute eigenvalues only;
!>          = 'V':  Compute eigenvalues and eigenvectors.
!>
!>          This parameter influences the choice of the algorithm and
!>          may alter the computed values.
!> 
[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.
!>          For RANGE = 'V' or 'I' and IU - IL < N - 1, SSTEBZ and
!>          CSTEIN are called
!>
!>          This parameter influences the choice of the algorithm and
!>          may alter the computed values.
!> 
[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,out]A
!>          A is COMPLEX array, dimension (LDA, N)
!>          On entry, the Hermitian matrix A.  If UPLO = 'U', the
!>          leading N-by-N upper triangular part of A contains the
!>          upper triangular part of the matrix A.  If UPLO = 'L',
!>          the leading N-by-N lower triangular part of A contains
!>          the lower triangular part of the matrix A.
!>          On exit, the lower triangle (if UPLO='L') or the upper
!>          triangle (if UPLO='U') of A, including the diagonal, is
!>          destroyed.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= 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 A to tridiagonal form.
!>
!>          See  by Demmel and
!>          Kahan, LAPACK Working Note #3.
!>
!>          If high relative accuracy is important, set ABSTOL to
!>          SLAMCH( 'Safe minimum' ).  Doing so will guarantee that
!>          eigenvalues are computed to high relative accuracy when
!>          possible in future releases.  The current code does not
!>          make any guarantees about high relative accuracy, but
!>          future releases will. See J. Barlow and J. Demmel,
!>          , LAPACK Working Note #7, for a discussion
!>          of which matrices define their eigenvalues to high relative
!>          accuracy.
!> 
[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 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.
!>          Supplying N columns is always safe.
!> 
[in]LDZ
!>          LDZ is INTEGER
!>          The leading dimension of the array Z.  LDZ >= 1, and if
!>          JOBZ = 'V', LDZ >= max(1,N).
!> 
[out]ISUPPZ
!>          ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
!>          The support of the eigenvectors in Z, i.e., the indices
!>          indicating the nonzero elements in Z. The i-th eigenvector
!>          is nonzero only in elements ISUPPZ( 2*i-1 ) through
!>          ISUPPZ( 2*i ). This is an output of CSTEMR (tridiagonal
!>          matrix). The support of the eigenvectors of A is typically
!>          1:N because of the unitary transformations applied by CUNMTR.
!>          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
!> 
[out]WORK
!>          WORK is COMPLEX array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of the array WORK.
!>          If N <= 1, LWORK >= 1, else LWORK >= 2*N.
!>          For optimal efficiency, LWORK >= (NB+1)*N,
!>          where NB is the max of the blocksize for CHETRD and for
!>          CUNMTR as returned by ILAENV.
!>
!>          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 (MAX(1,LRWORK))
!>          On exit, if INFO = 0, RWORK(1) returns the optimal
!>          (and minimal) LRWORK.
!> 
[in]LRWORK
!>          LRWORK is INTEGER
!>          The length of the array RWORK.
!>          If N <= 1, LRWORK >= 1, else LRWORK >= 24*N.
!>
!>          If LRWORK = -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]IWORK
!>          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
!>          On exit, if INFO = 0, IWORK(1) returns the optimal
!>          (and minimal) LIWORK.
!> 
[in]LIWORK
!>          LIWORK is INTEGER
!>          The dimension of the array IWORK.
!>          If N <= 1, LIWORK >= 1, else LIWORK >= 10*N.
!>
!>          If LIWORK = -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]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  Internal error
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Inderjit Dhillon, IBM Almaden, USA
Osni Marques, LBNL/NERSC, USA
Ken Stanley, Computer Science Division, University of California at Berkeley, USA
Jason Riedy, Computer Science Division, University of California at Berkeley, USA

Definition at line 369 of file cheevr.f.

373*
374* -- LAPACK driver routine --
375* -- LAPACK is a software package provided by Univ. of Tennessee, --
376* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
377*
378* .. Scalar Arguments ..
379 CHARACTER JOBZ, RANGE, UPLO
380 INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LRWORK, LWORK,
381 $ M, N
382 REAL ABSTOL, VL, VU
383* ..
384* .. Array Arguments ..
385 INTEGER ISUPPZ( * ), IWORK( * )
386 REAL RWORK( * ), W( * )
387 COMPLEX A( LDA, * ), WORK( * ), Z( LDZ, * )
388* ..
389*
390* =====================================================================
391*
392* .. Parameters ..
393 REAL ZERO, ONE, TWO
394 parameter( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0 )
395* ..
396* .. Local Scalars ..
397 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
398 $ WANTZ, TRYRAC
399 CHARACTER ORDER
400 INTEGER I, IEEEOK, IINFO, IMAX, INDIBL, INDIFL, INDISP,
401 $ INDIWO, INDRD, INDRDD, INDRE, INDREE, INDRWK,
402 $ INDTAU, INDWK, INDWKN, ISCALE, ITMP1, J, JJ,
403 $ LIWMIN, LLWORK, LLRWORK, LLWRKN, LRWMIN,
404 $ LWKOPT, LWMIN, NB, NSPLIT
405 REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
406 $ SIGMA, SMLNUM, TMP1, VLL, VUU
407* ..
408* .. External Functions ..
409 LOGICAL LSAME
410 INTEGER ILAENV
411 REAL CLANSY, SLAMCH, SROUNDUP_LWORK
412 EXTERNAL lsame, ilaenv, clansy, slamch,
414* ..
415* .. External Subroutines ..
416 EXTERNAL chetrd, csscal, cstemr, cstein, cswap,
417 $ cunmtr,
419* ..
420* .. Intrinsic Functions ..
421 INTRINSIC max, min, real, sqrt
422* ..
423* .. Executable Statements ..
424*
425* Test the input parameters.
426*
427 ieeeok = ilaenv( 10, 'CHEEVR', 'N', 1, 2, 3, 4 )
428*
429 lower = lsame( uplo, 'L' )
430 wantz = lsame( jobz, 'V' )
431 alleig = lsame( range, 'A' )
432 valeig = lsame( range, 'V' )
433 indeig = lsame( range, 'I' )
434*
435 lquery = ( ( lwork.EQ.-1 ) .OR. ( lrwork.EQ.-1 ) .OR.
436 $ ( liwork.EQ.-1 ) )
437*
438 IF( n.LE.1 ) THEN
439 lwmin = 1
440 lrwmin = 1
441 liwmin = 1
442 ELSE
443 lwmin = 2*n
444 lrwmin = 24*n
445 liwmin = 10*n
446 END IF
447*
448 info = 0
449 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
450 info = -1
451 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
452 info = -2
453 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
454 info = -3
455 ELSE IF( n.LT.0 ) THEN
456 info = -4
457 ELSE IF( lda.LT.max( 1, n ) ) THEN
458 info = -6
459 ELSE
460 IF( valeig ) THEN
461 IF( n.GT.0 .AND. vu.LE.vl )
462 $ info = -8
463 ELSE IF( indeig ) THEN
464 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
465 info = -9
466 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
467 info = -10
468 END IF
469 END IF
470 END IF
471 IF( info.EQ.0 ) THEN
472 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
473 info = -15
474 END IF
475 END IF
476*
477 IF( info.EQ.0 ) THEN
478 nb = ilaenv( 1, 'CHETRD', uplo, n, -1, -1, -1 )
479 nb = max( nb, ilaenv( 1, 'CUNMTR', uplo, n, -1, -1, -1 ) )
480 lwkopt = max( ( nb+1 )*n, lwmin )
481 work( 1 ) = sroundup_lwork( lwkopt )
482 rwork( 1 ) = sroundup_lwork( lrwmin )
483 iwork( 1 ) = liwmin
484*
485 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
486 info = -18
487 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
488 info = -20
489 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
490 info = -22
491 END IF
492 END IF
493*
494 IF( info.NE.0 ) THEN
495 CALL xerbla( 'CHEEVR', -info )
496 RETURN
497 ELSE IF( lquery ) THEN
498 RETURN
499 END IF
500*
501* Quick return if possible
502*
503 m = 0
504 IF( n.EQ.0 ) THEN
505 work( 1 ) = 1
506 RETURN
507 END IF
508*
509 IF( n.EQ.1 ) THEN
510 work( 1 ) = 1
511 IF( alleig .OR. indeig ) THEN
512 m = 1
513 w( 1 ) = real( a( 1, 1 ) )
514 ELSE
515 IF( vl.LT.real( a( 1, 1 ) ) .AND. vu.GE.real( a( 1, 1 ) ) )
516 $ THEN
517 m = 1
518 w( 1 ) = real( a( 1, 1 ) )
519 END IF
520 END IF
521 IF( wantz ) THEN
522 z( 1, 1 ) = one
523 isuppz( 1 ) = 1
524 isuppz( 2 ) = 1
525 END IF
526 RETURN
527 END IF
528*
529* Get machine constants.
530*
531 safmin = slamch( 'Safe minimum' )
532 eps = slamch( 'Precision' )
533 smlnum = safmin / eps
534 bignum = one / smlnum
535 rmin = sqrt( smlnum )
536 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
537*
538* Scale matrix to allowable range, if necessary.
539*
540 iscale = 0
541 abstll = abstol
542 IF (valeig) THEN
543 vll = vl
544 vuu = vu
545 END IF
546 anrm = clansy( 'M', uplo, n, a, lda, rwork )
547 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
548 iscale = 1
549 sigma = rmin / anrm
550 ELSE IF( anrm.GT.rmax ) THEN
551 iscale = 1
552 sigma = rmax / anrm
553 END IF
554 IF( iscale.EQ.1 ) THEN
555 IF( lower ) THEN
556 DO 10 j = 1, n
557 CALL csscal( n-j+1, sigma, a( j, j ), 1 )
558 10 CONTINUE
559 ELSE
560 DO 20 j = 1, n
561 CALL csscal( j, sigma, a( 1, j ), 1 )
562 20 CONTINUE
563 END IF
564 IF( abstol.GT.0 )
565 $ abstll = abstol*sigma
566 IF( valeig ) THEN
567 vll = vl*sigma
568 vuu = vu*sigma
569 END IF
570 END IF
571
572* Initialize indices into workspaces. Note: The IWORK indices are
573* used only if SSTERF or CSTEMR fail.
574
575* WORK(INDTAU:INDTAU+N-1) stores the complex scalar factors of the
576* elementary reflectors used in CHETRD.
577 indtau = 1
578* INDWK is the starting offset of the remaining complex workspace,
579* and LLWORK is the remaining complex workspace size.
580 indwk = indtau + n
581 llwork = lwork - indwk + 1
582
583* RWORK(INDRD:INDRD+N-1) stores the real tridiagonal's diagonal
584* entries.
585 indrd = 1
586* RWORK(INDRE:INDRE+N-1) stores the off-diagonal entries of the
587* tridiagonal matrix from CHETRD.
588 indre = indrd + n
589* RWORK(INDRDD:INDRDD+N-1) is a copy of the diagonal entries over
590* -written by CSTEMR (the SSTERF path copies the diagonal to W).
591 indrdd = indre + n
592* RWORK(INDREE:INDREE+N-1) is a copy of the off-diagonal entries over
593* -written while computing the eigenvalues in SSTERF and CSTEMR.
594 indree = indrdd + n
595* INDRWK is the starting offset of the left-over real workspace, and
596* LLRWORK is the remaining workspace size.
597 indrwk = indree + n
598 llrwork = lrwork - indrwk + 1
599
600* IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in SSTEBZ and
601* stores the block indices of each of the M<=N eigenvalues.
602 indibl = 1
603* IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in SSTEBZ and
604* stores the starting and finishing indices of each block.
605 indisp = indibl + n
606* IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
607* that corresponding to eigenvectors that fail to converge in
608* SSTEIN. This information is discarded; if any fail, the driver
609* returns INFO > 0.
610 indifl = indisp + n
611* INDIWO is the offset of the remaining integer workspace.
612 indiwo = indifl + n
613
614*
615* Call CHETRD to reduce Hermitian matrix to tridiagonal form.
616*
617 CALL chetrd( uplo, n, a, lda, rwork( indrd ), rwork( indre ),
618 $ work( indtau ), work( indwk ), llwork, iinfo )
619*
620* If all eigenvalues are desired
621* then call SSTERF or CSTEMR and CUNMTR.
622*
623 test = .false.
624 IF( indeig ) THEN
625 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
626 test = .true.
627 END IF
628 END IF
629 IF( ( alleig.OR.test ) .AND. ( ieeeok.EQ.1 ) ) THEN
630 IF( .NOT.wantz ) THEN
631 CALL scopy( n, rwork( indrd ), 1, w, 1 )
632 CALL scopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
633 CALL ssterf( n, w, rwork( indree ), info )
634 ELSE
635 CALL scopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
636 CALL scopy( n, rwork( indrd ), 1, rwork( indrdd ), 1 )
637*
638 IF (abstol .LE. two*real( n )*eps) THEN
639 tryrac = .true.
640 ELSE
641 tryrac = .false.
642 END IF
643 CALL cstemr( jobz, 'A', n, rwork( indrdd ),
644 $ rwork( indree ), vl, vu, il, iu, m, w,
645 $ z, ldz, n, isuppz, tryrac,
646 $ rwork( indrwk ), llrwork,
647 $ iwork, liwork, info )
648*
649* Apply unitary matrix used in reduction to tridiagonal
650* form to eigenvectors returned by CSTEMR.
651*
652 IF( wantz .AND. info.EQ.0 ) THEN
653 indwkn = indwk
654 llwrkn = lwork - indwkn + 1
655 CALL cunmtr( 'L', uplo, 'N', n, m, a, lda,
656 $ work( indtau ), z, ldz, work( indwkn ),
657 $ llwrkn, iinfo )
658 END IF
659 END IF
660*
661*
662 IF( info.EQ.0 ) THEN
663 m = n
664 GO TO 30
665 END IF
666 info = 0
667 END IF
668*
669* Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
670* Also call SSTEBZ and CSTEIN if CSTEMR fails.
671*
672 IF( wantz ) THEN
673 order = 'B'
674 ELSE
675 order = 'E'
676 END IF
677
678 CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
679 $ rwork( indrd ), rwork( indre ), m, nsplit, w,
680 $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
681 $ iwork( indiwo ), info )
682*
683 IF( wantz ) THEN
684 CALL cstein( n, rwork( indrd ), rwork( indre ), m, w,
685 $ iwork( indibl ), iwork( indisp ), z, ldz,
686 $ rwork( indrwk ), iwork( indiwo ), iwork( indifl ),
687 $ info )
688*
689* Apply unitary matrix used in reduction to tridiagonal
690* form to eigenvectors returned by CSTEIN.
691*
692 indwkn = indwk
693 llwrkn = lwork - indwkn + 1
694 CALL cunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
695 $ z,
696 $ ldz, work( indwkn ), llwrkn, iinfo )
697 END IF
698*
699* If matrix was scaled, then rescale eigenvalues appropriately.
700*
701 30 CONTINUE
702 IF( iscale.EQ.1 ) THEN
703 IF( info.EQ.0 ) THEN
704 imax = m
705 ELSE
706 imax = info - 1
707 END IF
708 CALL sscal( imax, one / sigma, w, 1 )
709 END IF
710*
711* If eigenvalues are not in order, then sort them, along with
712* eigenvectors.
713*
714 IF( wantz ) THEN
715 DO 50 j = 1, m - 1
716 i = 0
717 tmp1 = w( j )
718 DO 40 jj = j + 1, m
719 IF( w( jj ).LT.tmp1 ) THEN
720 i = jj
721 tmp1 = w( jj )
722 END IF
723 40 CONTINUE
724*
725 IF( i.NE.0 ) THEN
726 itmp1 = iwork( indibl+i-1 )
727 w( i ) = w( j )
728 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
729 w( j ) = tmp1
730 iwork( indibl+j-1 ) = itmp1
731 CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
732 END IF
733 50 CONTINUE
734 END IF
735*
736* Set WORK(1) to optimal workspace size.
737*
738 work( 1 ) = sroundup_lwork( lwkopt )
739 rwork( 1 ) = sroundup_lwork( lrwmin )
740 iwork( 1 ) = liwmin
741*
742 RETURN
743*
744* End of CHEEVR
745*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine chetrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
CHETRD
Definition chetrd.f:191
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clansy(norm, uplo, n, a, lda, work)
CLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clansy.f:121
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine csscal(n, sa, cx, incx)
CSSCAL
Definition csscal.f:78
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 cstemr(jobz, range, n, d, e, vl, vu, il, iu, m, w, z, ldz, nzc, isuppz, tryrac, work, lwork, iwork, liwork, info)
CSTEMR
Definition cstemr.f:337
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
subroutine cunmtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
CUNMTR
Definition cunmtr.f:171
Here is the call graph for this function:
Here is the caller graph for this function: