LAPACK  3.10.1
LAPACK: Linear Algebra PACKage

◆ ssyevr_2stage()

subroutine ssyevr_2stage ( character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
real  VL,
real  VU,
integer  IL,
integer  IU,
real  ABSTOL,
integer  M,
real, dimension( * )  W,
real, dimension( ldz, * )  Z,
integer  LDZ,
integer, dimension( * )  ISUPPZ,
real, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

SSYEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices

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

Purpose:
 SSYEVR_2STAGE computes selected eigenvalues and, optionally, eigenvectors
 of a real symmetric 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.

 SSYEVR_2STAGE first reduces the matrix A to tridiagonal form T with a call
 to SSYTRD.  Then, whenever possible, SSYEVR_2STAGE calls SSTEMR to compute
 the eigenspectrum using Relatively Robust Representations.  SSTEMR
 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 SSTEMR'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 : SSYEVR_2STAGE calls SSTEMR when the full spectrum is requested
 on machines which conform to the ieee-754 floating point standard.
 SSYEVR_2STAGE calls SSTEBZ and SSTEIN on non-ieee machines and
 when partial spectrum requests are made.

 Normal execution of SSTEMR 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
          SSTEIN 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 REAL array, dimension (LDA, N)
          On entry, the symmetric 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
          future 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 REAL 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 SSTEMR (tridiagonal
          matrix). The support of the eigenvectors of A is typically 
          1:N because of the orthogonal transformations applied by SORMTR.
          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
[out]WORK
          WORK is REAL 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 + 5*N
                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                                               + max(2*KD*KD, KD*NTHREADS) 
                                               + (KD+1)*N + 5*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 size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LWORK.
[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 size of the IWORK array,
          returns this value as the first entry of the IWORK array, and
          no error message related to 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 \n
Osni Marques, LBNL/NERSC, USA \n
Ken Stanley, Computer Science Division, University of
  California at Berkeley, USA \n
Jason Riedy, Computer Science Division, University of
  California at Berkeley, USA \n
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 378 of file ssyevr_2stage.f.

381 *
382  IMPLICIT NONE
383 *
384 * -- LAPACK driver routine --
385 * -- LAPACK is a software package provided by Univ. of Tennessee, --
386 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
387 *
388 * .. Scalar Arguments ..
389  CHARACTER JOBZ, RANGE, UPLO
390  INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
391  REAL ABSTOL, VL, VU
392 * ..
393 * .. Array Arguments ..
394  INTEGER ISUPPZ( * ), IWORK( * )
395  REAL A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
396 * ..
397 *
398 * =====================================================================
399 *
400 * .. Parameters ..
401  REAL ZERO, ONE, TWO
402  parameter( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0 )
403 * ..
404 * .. Local Scalars ..
405  LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, VALEIG, WANTZ,
406  $ TRYRAC, TEST
407  CHARACTER ORDER
408  INTEGER I, IEEEOK, IINFO, IMAX, INDD, INDDD, INDE,
409  $ INDEE, INDIBL, INDIFL, INDISP, INDIWO, INDTAU,
410  $ INDWK, INDWKN, ISCALE, J, JJ, LIWMIN,
411  $ LLWORK, LLWRKN, LWMIN, NSPLIT,
412  $ LHTRD, LWTRD, KD, IB, INDHOUS
413  REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
414  $ SIGMA, SMLNUM, TMP1, VLL, VUU
415 * ..
416 * .. External Functions ..
417  LOGICAL LSAME
418  INTEGER ILAENV, ILAENV2STAGE
419  REAL SLAMCH, SLANSY
420  EXTERNAL lsame, slamch, slansy, ilaenv, ilaenv2stage
421 * ..
422 * .. External Subroutines ..
423  EXTERNAL scopy, sormtr, sscal, sstebz, sstemr, sstein,
425 * ..
426 * .. Intrinsic Functions ..
427  INTRINSIC max, min, sqrt
428 * ..
429 * .. Executable Statements ..
430 *
431 * Test the input parameters.
432 *
433  ieeeok = ilaenv( 10, 'SSYEVR', 'N', 1, 2, 3, 4 )
434 *
435  lower = lsame( uplo, 'L' )
436  wantz = lsame( jobz, 'V' )
437  alleig = lsame( range, 'A' )
438  valeig = lsame( range, 'V' )
439  indeig = lsame( range, 'I' )
440 *
441  lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
442 *
443  kd = ilaenv2stage( 1, 'SSYTRD_2STAGE', jobz, n, -1, -1, -1 )
444  ib = ilaenv2stage( 2, 'SSYTRD_2STAGE', jobz, n, kd, -1, -1 )
445  lhtrd = ilaenv2stage( 3, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
446  lwtrd = ilaenv2stage( 4, 'SSYTRD_2STAGE', jobz, n, kd, ib, -1 )
447  lwmin = max( 26*n, 5*n + lhtrd + lwtrd )
448  liwmin = max( 1, 10*n )
449 *
450  info = 0
451  IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
452  info = -1
453  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
454  info = -2
455  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
456  info = -3
457  ELSE IF( n.LT.0 ) THEN
458  info = -4
459  ELSE IF( lda.LT.max( 1, n ) ) THEN
460  info = -6
461  ELSE
462  IF( valeig ) THEN
463  IF( n.GT.0 .AND. vu.LE.vl )
464  $ info = -8
465  ELSE IF( indeig ) THEN
466  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
467  info = -9
468  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
469  info = -10
470  END IF
471  END IF
472  END IF
473  IF( info.EQ.0 ) THEN
474  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
475  info = -15
476  ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
477  info = -18
478  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
479  info = -20
480  END IF
481  END IF
482 *
483  IF( info.EQ.0 ) THEN
484 * NB = ILAENV( 1, 'SSYTRD', UPLO, N, -1, -1, -1 )
485 * NB = MAX( NB, ILAENV( 1, 'SORMTR', UPLO, N, -1, -1, -1 ) )
486 * LWKOPT = MAX( ( NB+1 )*N, LWMIN )
487  work( 1 ) = lwmin
488  iwork( 1 ) = liwmin
489  END IF
490 *
491  IF( info.NE.0 ) THEN
492  CALL xerbla( 'SSYEVR_2STAGE', -info )
493  RETURN
494  ELSE IF( lquery ) THEN
495  RETURN
496  END IF
497 *
498 * Quick return if possible
499 *
500  m = 0
501  IF( n.EQ.0 ) THEN
502  work( 1 ) = 1
503  RETURN
504  END IF
505 *
506  IF( n.EQ.1 ) THEN
507  work( 1 ) = 26
508  IF( alleig .OR. indeig ) THEN
509  m = 1
510  w( 1 ) = a( 1, 1 )
511  ELSE
512  IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
513  m = 1
514  w( 1 ) = a( 1, 1 )
515  END IF
516  END IF
517  IF( wantz ) THEN
518  z( 1, 1 ) = one
519  isuppz( 1 ) = 1
520  isuppz( 2 ) = 1
521  END IF
522  RETURN
523  END IF
524 *
525 * Get machine constants.
526 *
527  safmin = slamch( 'Safe minimum' )
528  eps = slamch( 'Precision' )
529  smlnum = safmin / eps
530  bignum = one / smlnum
531  rmin = sqrt( smlnum )
532  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
533 *
534 * Scale matrix to allowable range, if necessary.
535 *
536  iscale = 0
537  abstll = abstol
538  IF (valeig) THEN
539  vll = vl
540  vuu = vu
541  END IF
542  anrm = slansy( 'M', uplo, n, a, lda, work )
543  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
544  iscale = 1
545  sigma = rmin / anrm
546  ELSE IF( anrm.GT.rmax ) THEN
547  iscale = 1
548  sigma = rmax / anrm
549  END IF
550  IF( iscale.EQ.1 ) THEN
551  IF( lower ) THEN
552  DO 10 j = 1, n
553  CALL sscal( n-j+1, sigma, a( j, j ), 1 )
554  10 CONTINUE
555  ELSE
556  DO 20 j = 1, n
557  CALL sscal( j, sigma, a( 1, j ), 1 )
558  20 CONTINUE
559  END IF
560  IF( abstol.GT.0 )
561  $ abstll = abstol*sigma
562  IF( valeig ) THEN
563  vll = vl*sigma
564  vuu = vu*sigma
565  END IF
566  END IF
567 
568 * Initialize indices into workspaces. Note: The IWORK indices are
569 * used only if SSTERF or SSTEMR fail.
570 
571 * WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
572 * elementary reflectors used in SSYTRD.
573  indtau = 1
574 * WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
575  indd = indtau + n
576 * WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
577 * tridiagonal matrix from SSYTRD.
578  inde = indd + n
579 * WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
580 * -written by SSTEMR (the SSTERF path copies the diagonal to W).
581  inddd = inde + n
582 * WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
583 * -written while computing the eigenvalues in SSTERF and SSTEMR.
584  indee = inddd + n
585 * INDHOUS is the starting offset Householder storage of stage 2
586  indhous = indee + n
587 * INDWK is the starting offset of the left-over workspace, and
588 * LLWORK is the remaining workspace size.
589  indwk = indhous + lhtrd
590  llwork = lwork - indwk + 1
591 
592 
593 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in SSTEBZ and
594 * stores the block indices of each of the M<=N eigenvalues.
595  indibl = 1
596 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in SSTEBZ and
597 * stores the starting and finishing indices of each block.
598  indisp = indibl + n
599 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
600 * that corresponding to eigenvectors that fail to converge in
601 * SSTEIN. This information is discarded; if any fail, the driver
602 * returns INFO > 0.
603  indifl = indisp + n
604 * INDIWO is the offset of the remaining integer workspace.
605  indiwo = indifl + n
606 
607 *
608 * Call SSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
609 *
610 *
611  CALL ssytrd_2stage( jobz, uplo, n, a, lda, work( indd ),
612  $ work( inde ), work( indtau ), work( indhous ),
613  $ lhtrd, work( indwk ), llwork, iinfo )
614 *
615 * If all eigenvalues are desired
616 * then call SSTERF or SSTEMR and SORMTR.
617 *
618  test = .false.
619  IF( indeig ) THEN
620  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
621  test = .true.
622  END IF
623  END IF
624  IF( ( alleig.OR.test ) .AND. ( ieeeok.EQ.1 ) ) THEN
625  IF( .NOT.wantz ) THEN
626  CALL scopy( n, work( indd ), 1, w, 1 )
627  CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
628  CALL ssterf( n, w, work( indee ), info )
629  ELSE
630  CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
631  CALL scopy( n, work( indd ), 1, work( inddd ), 1 )
632 *
633  IF (abstol .LE. two*n*eps) THEN
634  tryrac = .true.
635  ELSE
636  tryrac = .false.
637  END IF
638  CALL sstemr( jobz, 'A', n, work( inddd ), work( indee ),
639  $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
640  $ tryrac, work( indwk ), lwork, iwork, liwork,
641  $ info )
642 *
643 *
644 *
645 * Apply orthogonal matrix used in reduction to tridiagonal
646 * form to eigenvectors returned by SSTEMR.
647 *
648  IF( wantz .AND. info.EQ.0 ) THEN
649  indwkn = inde
650  llwrkn = lwork - indwkn + 1
651  CALL sormtr( 'L', uplo, 'N', n, m, a, lda,
652  $ work( indtau ), z, ldz, work( indwkn ),
653  $ llwrkn, iinfo )
654  END IF
655  END IF
656 *
657 *
658  IF( info.EQ.0 ) THEN
659 * Everything worked. Skip SSTEBZ/SSTEIN. IWORK(:) are
660 * undefined.
661  m = n
662  GO TO 30
663  END IF
664  info = 0
665  END IF
666 *
667 * Otherwise, call SSTEBZ and, if eigenvectors are desired, SSTEIN.
668 * Also call SSTEBZ and SSTEIN if SSTEMR fails.
669 *
670  IF( wantz ) THEN
671  order = 'B'
672  ELSE
673  order = 'E'
674  END IF
675 
676  CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
677  $ work( indd ), work( inde ), m, nsplit, w,
678  $ iwork( indibl ), iwork( indisp ), work( indwk ),
679  $ iwork( indiwo ), info )
680 *
681  IF( wantz ) THEN
682  CALL sstein( n, work( indd ), work( inde ), m, w,
683  $ iwork( indibl ), iwork( indisp ), z, ldz,
684  $ work( indwk ), iwork( indiwo ), iwork( indifl ),
685  $ info )
686 *
687 * Apply orthogonal matrix used in reduction to tridiagonal
688 * form to eigenvectors returned by SSTEIN.
689 *
690  indwkn = inde
691  llwrkn = lwork - indwkn + 1
692  CALL sormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
693  $ ldz, work( indwkn ), llwrkn, iinfo )
694  END IF
695 *
696 * If matrix was scaled, then rescale eigenvalues appropriately.
697 *
698 * Jump here if SSTEMR/SSTEIN succeeded.
699  30 CONTINUE
700  IF( iscale.EQ.1 ) THEN
701  IF( info.EQ.0 ) THEN
702  imax = m
703  ELSE
704  imax = info - 1
705  END IF
706  CALL sscal( imax, one / sigma, w, 1 )
707  END IF
708 *
709 * If eigenvalues are not in order, then sort them, along with
710 * eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
711 * It may not be initialized (if SSTEMR/SSTEIN succeeded), and we do
712 * not return this detailed information to the user.
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  w( i ) = w( j )
727  w( j ) = tmp1
728  CALL sswap( n, z( 1, i ), 1, z( 1, j ), 1 )
729  END IF
730  50 CONTINUE
731  END IF
732 *
733 * Set WORK(1) to optimal workspace size.
734 *
735  work( 1 ) = lwmin
736  iwork( 1 ) = liwmin
737 *
738  RETURN
739 *
740 * End of SSYEVR_2STAGE
741 *
integer function ilaenv2stage(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV2STAGE
Definition: ilaenv2stage.f:149
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:86
subroutine sstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
SSTEBZ
Definition: sstebz.f:273
subroutine sormtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
SORMTR
Definition: sormtr.f:172
subroutine sstemr(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)
SSTEMR
Definition: sstemr.f:321
subroutine sstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
SSTEIN
Definition: sstein.f:174
real function slansy(NORM, UPLO, N, A, LDA, WORK)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition: slansy.f:122
subroutine ssytrd_2stage(VECT, UPLO, N, A, LDA, D, E, TAU, HOUS2, LHOUS2, WORK, LWORK, INFO)
SSYTRD_2STAGE
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP
Definition: sswap.f:82
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:82
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:79
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:68
Here is the call graph for this function:
Here is the caller graph for this function: