LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
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.
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 383 of file ssyevr_2stage.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: