LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine cheevr_2stage ( 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_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices

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

Purpose:
 CHEEVR_2STAGE computes selected eigenvalues and, optionally, eigenvectors
 of a complex Hermitian 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.

 CHEEVR_2STAGE first reduces the matrix A to tridiagonal form T with a call
 to CHETRD.  Then, whenever possible, CHEEVR_2STAGE calls CSTEMR to compute
 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_2STAGE calls CSTEMR when the full spectrum is requested
 on machines which conform to the ieee-754 floating point standard.
 CHEEVR_2STAGE 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.
                  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.
          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 ). 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 dimension of the array WORK.  
          If JOBZ = 'N' and N > 1, LWORK must be queried.
                                   LWORK = MAX(1, 26*N, dimension) where
                                   dimension = max(stage1,stage2) + (KD+1)*N + N
                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                                               + max(2*KD*KD, KD*NTHREADS) 
                                               + (KD+1)*N + N
                                   where KD is the blocking size of the reduction,
                                   FACTOPTNB is the blocking used by the QR or LQ
                                   algorithm, usually FACTOPTNB=128 is a good choice
                                   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 (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
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 408 of file cheevr_2stage.f.

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