LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.

 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 "good" 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 DSTEMR's documentation and:
 - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
   to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
 - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
   Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
   2004.  Also LAPACK Working Note 154.
 - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
   tridiagonal eigenvalue/eigenvector problem",
   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.
[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
[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 "Computing Small Singular Values of Bidiagonal Matrices
          with Guaranteed High Relative Accuracy," 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
          furutre releases will. See J. Barlow and J. Demmel,
          "Computing Accurate Eigensystems of Scaled Diagonally
          Dominant Matrices", 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.
[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 ).
          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.  LWORK >= max(1,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.  LRWORK >= max(1,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.  LIWORK >= max(1,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.
Date
June 2016
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 357 of file cheevr.f.

357 *
358 * -- LAPACK driver routine (version 3.6.1) --
359 * -- LAPACK is a software package provided by Univ. of Tennessee, --
360 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
361 * June 2016
362 *
363 * .. Scalar Arguments ..
364  CHARACTER jobz, range, uplo
365  INTEGER il, info, iu, lda, ldz, liwork, lrwork, lwork,
366  $ m, n
367  REAL abstol, vl, vu
368 * ..
369 * .. Array Arguments ..
370  INTEGER isuppz( * ), iwork( * )
371  REAL rwork( * ), w( * )
372  COMPLEX a( lda, * ), work( * ), z( ldz, * )
373 * ..
374 *
375 * =====================================================================
376 *
377 * .. Parameters ..
378  REAL zero, one, two
379  parameter ( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0 )
380 * ..
381 * .. Local Scalars ..
382  LOGICAL alleig, indeig, lower, lquery, test, valeig,
383  $ wantz, tryrac
384  CHARACTER order
385  INTEGER i, ieeeok, iinfo, imax, indibl, indifl, indisp,
386  $ indiwo, indrd, indrdd, indre, indree, indrwk,
387  $ indtau, indwk, indwkn, iscale, itmp1, j, jj,
388  $ liwmin, llwork, llrwork, llwrkn, lrwmin,
389  $ lwkopt, lwmin, nb, nsplit
390  REAL abstll, anrm, bignum, eps, rmax, rmin, safmin,
391  $ sigma, smlnum, tmp1, vll, vuu
392 * ..
393 * .. External Functions ..
394  LOGICAL lsame
395  INTEGER ilaenv
396  REAL clansy, slamch
397  EXTERNAL lsame, ilaenv, clansy, slamch
398 * ..
399 * .. External Subroutines ..
400  EXTERNAL chetrd, csscal, cstemr, cstein, cswap, cunmtr,
402 * ..
403 * .. Intrinsic Functions ..
404  INTRINSIC max, min, REAL, sqrt
405 * ..
406 * .. Executable Statements ..
407 *
408 * Test the input parameters.
409 *
410  ieeeok = ilaenv( 10, 'CHEEVR', 'N', 1, 2, 3, 4 )
411 *
412  lower = lsame( uplo, 'L' )
413  wantz = lsame( jobz, 'V' )
414  alleig = lsame( range, 'A' )
415  valeig = lsame( range, 'V' )
416  indeig = lsame( range, 'I' )
417 *
418  lquery = ( ( lwork.EQ.-1 ) .OR. ( lrwork.EQ.-1 ) .OR.
419  $ ( liwork.EQ.-1 ) )
420 *
421  lrwmin = max( 1, 24*n )
422  liwmin = max( 1, 10*n )
423  lwmin = max( 1, 2*n )
424 *
425  info = 0
426  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
427  info = -1
428  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
429  info = -2
430  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
431  info = -3
432  ELSE IF( n.LT.0 ) THEN
433  info = -4
434  ELSE IF( lda.LT.max( 1, n ) ) THEN
435  info = -6
436  ELSE
437  IF( valeig ) THEN
438  IF( n.GT.0 .AND. vu.LE.vl )
439  $ info = -8
440  ELSE IF( indeig ) THEN
441  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
442  info = -9
443  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
444  info = -10
445  END IF
446  END IF
447  END IF
448  IF( info.EQ.0 ) THEN
449  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
450  info = -15
451  END IF
452  END IF
453 *
454  IF( info.EQ.0 ) THEN
455  nb = ilaenv( 1, 'CHETRD', uplo, n, -1, -1, -1 )
456  nb = max( nb, ilaenv( 1, 'CUNMTR', uplo, n, -1, -1, -1 ) )
457  lwkopt = max( ( nb+1 )*n, lwmin )
458  work( 1 ) = lwkopt
459  rwork( 1 ) = lrwmin
460  iwork( 1 ) = liwmin
461 *
462  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
463  info = -18
464  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
465  info = -20
466  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
467  info = -22
468  END IF
469  END IF
470 *
471  IF( info.NE.0 ) THEN
472  CALL xerbla( 'CHEEVR', -info )
473  RETURN
474  ELSE IF( lquery ) THEN
475  RETURN
476  END IF
477 *
478 * Quick return if possible
479 *
480  m = 0
481  IF( n.EQ.0 ) THEN
482  work( 1 ) = 1
483  RETURN
484  END IF
485 *
486  IF( n.EQ.1 ) THEN
487  work( 1 ) = 2
488  IF( alleig .OR. indeig ) THEN
489  m = 1
490  w( 1 ) = REAL( A( 1, 1 ) )
491  ELSE
492  IF( vl.LT.REAL( A( 1, 1 ) ) .AND. vu.GE.REAL( A( 1, 1 ) ) )
493  $ THEN
494  m = 1
495  w( 1 ) = REAL( A( 1, 1 ) )
496  END IF
497  END IF
498  IF( wantz ) THEN
499  z( 1, 1 ) = one
500  isuppz( 1 ) = 1
501  isuppz( 2 ) = 1
502  END IF
503  RETURN
504  END IF
505 *
506 * Get machine constants.
507 *
508  safmin = slamch( 'Safe minimum' )
509  eps = slamch( 'Precision' )
510  smlnum = safmin / eps
511  bignum = one / smlnum
512  rmin = sqrt( smlnum )
513  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
514 *
515 * Scale matrix to allowable range, if necessary.
516 *
517  iscale = 0
518  abstll = abstol
519  IF (valeig) THEN
520  vll = vl
521  vuu = vu
522  END IF
523  anrm = clansy( 'M', uplo, n, a, lda, rwork )
524  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
525  iscale = 1
526  sigma = rmin / anrm
527  ELSE IF( anrm.GT.rmax ) THEN
528  iscale = 1
529  sigma = rmax / anrm
530  END IF
531  IF( iscale.EQ.1 ) THEN
532  IF( lower ) THEN
533  DO 10 j = 1, n
534  CALL csscal( n-j+1, sigma, a( j, j ), 1 )
535  10 CONTINUE
536  ELSE
537  DO 20 j = 1, n
538  CALL csscal( j, sigma, a( 1, j ), 1 )
539  20 CONTINUE
540  END IF
541  IF( abstol.GT.0 )
542  $ abstll = abstol*sigma
543  IF( valeig ) THEN
544  vll = vl*sigma
545  vuu = vu*sigma
546  END IF
547  END IF
548 
549 * Initialize indices into workspaces. Note: The IWORK indices are
550 * used only if SSTERF or CSTEMR fail.
551 
552 * WORK(INDTAU:INDTAU+N-1) stores the complex scalar factors of the
553 * elementary reflectors used in CHETRD.
554  indtau = 1
555 * INDWK is the starting offset of the remaining complex workspace,
556 * and LLWORK is the remaining complex workspace size.
557  indwk = indtau + n
558  llwork = lwork - indwk + 1
559 
560 * RWORK(INDRD:INDRD+N-1) stores the real tridiagonal's diagonal
561 * entries.
562  indrd = 1
563 * RWORK(INDRE:INDRE+N-1) stores the off-diagonal entries of the
564 * tridiagonal matrix from CHETRD.
565  indre = indrd + n
566 * RWORK(INDRDD:INDRDD+N-1) is a copy of the diagonal entries over
567 * -written by CSTEMR (the SSTERF path copies the diagonal to W).
568  indrdd = indre + n
569 * RWORK(INDREE:INDREE+N-1) is a copy of the off-diagonal entries over
570 * -written while computing the eigenvalues in SSTERF and CSTEMR.
571  indree = indrdd + n
572 * INDRWK is the starting offset of the left-over real workspace, and
573 * LLRWORK is the remaining workspace size.
574  indrwk = indree + n
575  llrwork = lrwork - indrwk + 1
576 
577 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in SSTEBZ and
578 * stores the block indices of each of the M<=N eigenvalues.
579  indibl = 1
580 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in SSTEBZ and
581 * stores the starting and finishing indices of each block.
582  indisp = indibl + n
583 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
584 * that corresponding to eigenvectors that fail to converge in
585 * SSTEIN. This information is discarded; if any fail, the driver
586 * returns INFO > 0.
587  indifl = indisp + n
588 * INDIWO is the offset of the remaining integer workspace.
589  indiwo = indifl + n
590 
591 *
592 * Call CHETRD to reduce Hermitian matrix to tridiagonal form.
593 *
594  CALL chetrd( uplo, n, a, lda, rwork( indrd ), rwork( indre ),
595  $ work( indtau ), work( indwk ), llwork, iinfo )
596 *
597 * If all eigenvalues are desired
598 * then call SSTERF or CSTEMR and CUNMTR.
599 *
600  test = .false.
601  IF( indeig ) THEN
602  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
603  test = .true.
604  END IF
605  END IF
606  IF( ( alleig.OR.test ) .AND. ( ieeeok.EQ.1 ) ) THEN
607  IF( .NOT.wantz ) THEN
608  CALL scopy( n, rwork( indrd ), 1, w, 1 )
609  CALL scopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
610  CALL ssterf( n, w, rwork( indree ), info )
611  ELSE
612  CALL scopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
613  CALL scopy( n, rwork( indrd ), 1, rwork( indrdd ), 1 )
614 *
615  IF (abstol .LE. two*n*eps) THEN
616  tryrac = .true.
617  ELSE
618  tryrac = .false.
619  END IF
620  CALL cstemr( jobz, 'A', n, rwork( indrdd ),
621  $ rwork( indree ), vl, vu, il, iu, m, w,
622  $ z, ldz, n, isuppz, tryrac,
623  $ rwork( indrwk ), llrwork,
624  $ iwork, liwork, info )
625 *
626 * Apply unitary matrix used in reduction to tridiagonal
627 * form to eigenvectors returned by CSTEIN.
628 *
629  IF( wantz .AND. info.EQ.0 ) THEN
630  indwkn = indwk
631  llwrkn = lwork - indwkn + 1
632  CALL cunmtr( 'L', uplo, 'N', n, m, a, lda,
633  $ work( indtau ), z, ldz, work( indwkn ),
634  $ llwrkn, iinfo )
635  END IF
636  END IF
637 *
638 *
639  IF( info.EQ.0 ) THEN
640  m = n
641  GO TO 30
642  END IF
643  info = 0
644  END IF
645 *
646 * Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
647 * Also call SSTEBZ and CSTEIN if CSTEMR fails.
648 *
649  IF( wantz ) THEN
650  order = 'B'
651  ELSE
652  order = 'E'
653  END IF
654 
655  CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
656  $ rwork( indrd ), rwork( indre ), m, nsplit, w,
657  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
658  $ iwork( indiwo ), info )
659 *
660  IF( wantz ) THEN
661  CALL cstein( n, rwork( indrd ), rwork( indre ), m, w,
662  $ iwork( indibl ), iwork( indisp ), z, ldz,
663  $ rwork( indrwk ), iwork( indiwo ), iwork( indifl ),
664  $ info )
665 *
666 * Apply unitary matrix used in reduction to tridiagonal
667 * form to eigenvectors returned by CSTEIN.
668 *
669  indwkn = indwk
670  llwrkn = lwork - indwkn + 1
671  CALL cunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
672  $ ldz, work( indwkn ), llwrkn, iinfo )
673  END IF
674 *
675 * If matrix was scaled, then rescale eigenvalues appropriately.
676 *
677  30 CONTINUE
678  IF( iscale.EQ.1 ) THEN
679  IF( info.EQ.0 ) THEN
680  imax = m
681  ELSE
682  imax = info - 1
683  END IF
684  CALL sscal( imax, one / sigma, w, 1 )
685  END IF
686 *
687 * If eigenvalues are not in order, then sort them, along with
688 * eigenvectors.
689 *
690  IF( wantz ) THEN
691  DO 50 j = 1, m - 1
692  i = 0
693  tmp1 = w( j )
694  DO 40 jj = j + 1, m
695  IF( w( jj ).LT.tmp1 ) THEN
696  i = jj
697  tmp1 = w( jj )
698  END IF
699  40 CONTINUE
700 *
701  IF( i.NE.0 ) THEN
702  itmp1 = iwork( indibl+i-1 )
703  w( i ) = w( j )
704  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
705  w( j ) = tmp1
706  iwork( indibl+j-1 ) = itmp1
707  CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
708  END IF
709  50 CONTINUE
710  END IF
711 *
712 * Set WORK(1) to optimal workspace size.
713 *
714  work( 1 ) = lwkopt
715  rwork( 1 ) = lrwmin
716  iwork( 1 ) = liwmin
717 *
718  RETURN
719 *
720 * End of CHEEVR
721 *
subroutine sstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
SSTEBZ
Definition: sstebz.f:275
subroutine cunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMTR
Definition: cunmtr.f:174
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CSTEIN
Definition: cstein.f:184
subroutine chetrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
CHETRD
Definition: chetrd.f:194
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, or the element of largest absolute value of a complex symmetric matrix.
Definition: clansy.f:125
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
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:340
subroutine cswap(N, CX, INCX, CY, INCY)
CSWAP
Definition: cswap.f:52
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53
subroutine csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function: