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

◆ dsyevr_2stage()

subroutine dsyevr_2stage ( character  jobz,
character  range,
character  uplo,
integer  n,
double precision, dimension( lda, * )  a,
integer  lda,
double precision  vl,
double precision  vu,
integer  il,
integer  iu,
double precision  abstol,
integer  m,
double precision, dimension( * )  w,
double precision, dimension( ldz, * )  z,
integer  ldz,
integer, dimension( * )  isuppz,
double precision, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

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

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

Purpose:
 DSYEVR_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.

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

 Normal execution of DSTEMR 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, DSTEBZ and
          DSTEIN 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 DOUBLE PRECISION 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 DOUBLE PRECISION
          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 DOUBLE PRECISION
          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 DOUBLE PRECISION
          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
          DLAMCH( '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 DOUBLE PRECISION array, dimension (N)
          The first M elements contain the selected eigenvalues in
          ascending order.
[out]Z
          Z is DOUBLE PRECISION 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 DSTEMR (tridiagonal
          matrix). The support of the eigenvectors of A is typically
          1:N because of the orthogonal transformations applied by DORMTR.
          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
[out]WORK
          WORK is DOUBLE PRECISION 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 dsyevr_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 DOUBLE PRECISION ABSTOL, VL, VU
392* ..
393* .. Array Arguments ..
394 INTEGER ISUPPZ( * ), IWORK( * )
395 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
396* ..
397*
398* =====================================================================
399*
400* .. Parameters ..
401 DOUBLE PRECISION ZERO, ONE, TWO
402 parameter( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
403* ..
404* .. Local Scalars ..
405 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, VALEIG, WANTZ,
406 $ TRYRAC
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 DOUBLE PRECISION 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 DOUBLE PRECISION DLAMCH, DLANSY
421* ..
422* .. External Subroutines ..
423 EXTERNAL dcopy, dormtr, dscal, dstebz, dstemr, dstein,
425* ..
426* .. Intrinsic Functions ..
427 INTRINSIC max, min, sqrt
428* ..
429* .. Executable Statements ..
430*
431* Test the input parameters.
432*
433 ieeeok = ilaenv( 10, 'DSYEVR', '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, 'DSYTRD_2STAGE', jobz, n, -1, -1, -1 )
444 ib = ilaenv2stage( 2, 'DSYTRD_2STAGE', jobz, n, kd, -1, -1 )
445 lhtrd = ilaenv2stage( 3, 'DSYTRD_2STAGE', jobz, n, kd, ib, -1 )
446 lwtrd = ilaenv2stage( 4, 'DSYTRD_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, 'DSYTRD', UPLO, N, -1, -1, -1 )
485* NB = MAX( NB, ILAENV( 1, 'DORMTR', 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( 'DSYEVR_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 ) = 7
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 = dlamch( 'Safe minimum' )
528 eps = dlamch( '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 = dlansy( '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 dscal( n-j+1, sigma, a( j, j ), 1 )
554 10 CONTINUE
555 ELSE
556 DO 20 j = 1, n
557 CALL dscal( 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 DSTERF or DSTEMR fail.
570
571* WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
572* elementary reflectors used in DSYTRD.
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 DSYTRD.
578 inde = indd + n
579* WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
580* -written by DSTEMR (the DSTERF 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 DSTERF and DSTEMR.
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 DSTEBZ 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 DSTEBZ 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* DSTEIN. 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 DSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
609*
610*
611 CALL dsytrd_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 DSTERF or DSTEMR and DORMTR.
617*
618 IF( ( alleig .OR. ( indeig .AND. il.EQ.1 .AND. iu.EQ.n ) ) .AND.
619 $ ieeeok.EQ.1 ) THEN
620 IF( .NOT.wantz ) THEN
621 CALL dcopy( n, work( indd ), 1, w, 1 )
622 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
623 CALL dsterf( n, w, work( indee ), info )
624 ELSE
625 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
626 CALL dcopy( n, work( indd ), 1, work( inddd ), 1 )
627*
628 IF (abstol .LE. two*n*eps) THEN
629 tryrac = .true.
630 ELSE
631 tryrac = .false.
632 END IF
633 CALL dstemr( jobz, 'A', n, work( inddd ), work( indee ),
634 $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
635 $ tryrac, work( indwk ), lwork, iwork, liwork,
636 $ info )
637*
638*
639*
640* Apply orthogonal matrix used in reduction to tridiagonal
641* form to eigenvectors returned by DSTEMR.
642*
643 IF( wantz .AND. info.EQ.0 ) THEN
644 indwkn = inde
645 llwrkn = lwork - indwkn + 1
646 CALL dormtr( 'L', uplo, 'N', n, m, a, lda,
647 $ work( indtau ), z, ldz, work( indwkn ),
648 $ llwrkn, iinfo )
649 END IF
650 END IF
651*
652*
653 IF( info.EQ.0 ) THEN
654* Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
655* undefined.
656 m = n
657 GO TO 30
658 END IF
659 info = 0
660 END IF
661*
662* Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
663* Also call DSTEBZ and DSTEIN if DSTEMR fails.
664*
665 IF( wantz ) THEN
666 order = 'B'
667 ELSE
668 order = 'E'
669 END IF
670
671 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
672 $ work( indd ), work( inde ), m, nsplit, w,
673 $ iwork( indibl ), iwork( indisp ), work( indwk ),
674 $ iwork( indiwo ), info )
675*
676 IF( wantz ) THEN
677 CALL dstein( n, work( indd ), work( inde ), m, w,
678 $ iwork( indibl ), iwork( indisp ), z, ldz,
679 $ work( indwk ), iwork( indiwo ), iwork( indifl ),
680 $ info )
681*
682* Apply orthogonal matrix used in reduction to tridiagonal
683* form to eigenvectors returned by DSTEIN.
684*
685 indwkn = inde
686 llwrkn = lwork - indwkn + 1
687 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
688 $ ldz, work( indwkn ), llwrkn, iinfo )
689 END IF
690*
691* If matrix was scaled, then rescale eigenvalues appropriately.
692*
693* Jump here if DSTEMR/DSTEIN succeeded.
694 30 CONTINUE
695 IF( iscale.EQ.1 ) THEN
696 IF( info.EQ.0 ) THEN
697 imax = m
698 ELSE
699 imax = info - 1
700 END IF
701 CALL dscal( imax, one / sigma, w, 1 )
702 END IF
703*
704* If eigenvalues are not in order, then sort them, along with
705* eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
706* It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
707* not return this detailed information to the user.
708*
709 IF( wantz ) THEN
710 DO 50 j = 1, m - 1
711 i = 0
712 tmp1 = w( j )
713 DO 40 jj = j + 1, m
714 IF( w( jj ).LT.tmp1 ) THEN
715 i = jj
716 tmp1 = w( jj )
717 END IF
718 40 CONTINUE
719*
720 IF( i.NE.0 ) THEN
721 w( i ) = w( j )
722 w( j ) = tmp1
723 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
724 END IF
725 50 CONTINUE
726 END IF
727*
728* Set WORK(1) to optimal workspace size.
729*
730 work( 1 ) = lwmin
731 iwork( 1 ) = liwmin
732*
733 RETURN
734*
735* End of DSYEVR_2STAGE
736*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dsytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
DSYTRD_2STAGE
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlansy(norm, uplo, n, a, lda, work)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansy.f:122
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
DSTEBZ
Definition dstebz.f:273
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:174
subroutine dstemr(jobz, range, n, d, e, vl, vu, il, iu, m, w, z, ldz, nzc, isuppz, tryrac, work, lwork, iwork, liwork, info)
DSTEMR
Definition dstemr.f:322
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
DORMTR
Definition dormtr.f:171
Here is the call graph for this function:
Here is the caller graph for this function: