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

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

Purpose:
 DSYEVR computes selected eigenvalues and, optionally, eigenvectors
 of a real symmetric matrix A.  Eigenvalues and eigenvectors can be
 selected by specifying either a range of values or a range of
 indices for the desired eigenvalues.

 DSYEVR first reduces the matrix A to tridiagonal form T with a call
 to DSYTRD.  Then, whenever possible, DSYEVR 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 calls DSTEMR when the full spectrum is requested
 on machines which conform to the ieee-754 floating point standard.
 DSYEVR 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.
[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.  LWORK >= max(1,26*N).
          For optimal efficiency, LWORK >= (NB+6)*N,
          where NB is the max of the blocksize for DSYTRD and DORMTR
          returned by ILAENV.

          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

Definition at line 336 of file dsyevr.f.

336 *
337 * -- LAPACK driver routine (version 3.7.0) --
338 * -- LAPACK is a software package provided by Univ. of Tennessee, --
339 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
340 * June 2016
341 *
342 * .. Scalar Arguments ..
343  CHARACTER jobz, range, uplo
344  INTEGER il, info, iu, lda, ldz, liwork, lwork, m, n
345  DOUBLE PRECISION abstol, vl, vu
346 * ..
347 * .. Array Arguments ..
348  INTEGER isuppz( * ), iwork( * )
349  DOUBLE PRECISION a( lda, * ), w( * ), work( * ), z( ldz, * )
350 * ..
351 *
352 * =====================================================================
353 *
354 * .. Parameters ..
355  DOUBLE PRECISION zero, one, two
356  parameter ( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
357 * ..
358 * .. Local Scalars ..
359  LOGICAL alleig, indeig, lower, lquery, valeig, wantz,
360  $ tryrac
361  CHARACTER order
362  INTEGER i, ieeeok, iinfo, imax, indd, inddd, inde,
363  $ indee, indibl, indifl, indisp, indiwo, indtau,
364  $ indwk, indwkn, iscale, j, jj, liwmin,
365  $ llwork, llwrkn, lwkopt, lwmin, nb, nsplit
366  DOUBLE PRECISION abstll, anrm, bignum, eps, rmax, rmin, safmin,
367  $ sigma, smlnum, tmp1, vll, vuu
368 * ..
369 * .. External Functions ..
370  LOGICAL lsame
371  INTEGER ilaenv
372  DOUBLE PRECISION dlamch, dlansy
373  EXTERNAL lsame, ilaenv, dlamch, dlansy
374 * ..
375 * .. External Subroutines ..
376  EXTERNAL dcopy, dormtr, dscal, dstebz, dstemr, dstein,
378 * ..
379 * .. Intrinsic Functions ..
380  INTRINSIC max, min, sqrt
381 * ..
382 * .. Executable Statements ..
383 *
384 * Test the input parameters.
385 *
386  ieeeok = ilaenv( 10, 'DSYEVR', 'N', 1, 2, 3, 4 )
387 *
388  lower = lsame( uplo, 'L' )
389  wantz = lsame( jobz, 'V' )
390  alleig = lsame( range, 'A' )
391  valeig = lsame( range, 'V' )
392  indeig = lsame( range, 'I' )
393 *
394  lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
395 *
396  lwmin = max( 1, 26*n )
397  liwmin = max( 1, 10*n )
398 *
399  info = 0
400  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
401  info = -1
402  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
403  info = -2
404  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
405  info = -3
406  ELSE IF( n.LT.0 ) THEN
407  info = -4
408  ELSE IF( lda.LT.max( 1, n ) ) THEN
409  info = -6
410  ELSE
411  IF( valeig ) THEN
412  IF( n.GT.0 .AND. vu.LE.vl )
413  $ info = -8
414  ELSE IF( indeig ) THEN
415  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
416  info = -9
417  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
418  info = -10
419  END IF
420  END IF
421  END IF
422  IF( info.EQ.0 ) THEN
423  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
424  info = -15
425  ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
426  info = -18
427  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
428  info = -20
429  END IF
430  END IF
431 *
432  IF( info.EQ.0 ) THEN
433  nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
434  nb = max( nb, ilaenv( 1, 'DORMTR', uplo, n, -1, -1, -1 ) )
435  lwkopt = max( ( nb+1 )*n, lwmin )
436  work( 1 ) = lwkopt
437  iwork( 1 ) = liwmin
438  END IF
439 *
440  IF( info.NE.0 ) THEN
441  CALL xerbla( 'DSYEVR', -info )
442  RETURN
443  ELSE IF( lquery ) THEN
444  RETURN
445  END IF
446 *
447 * Quick return if possible
448 *
449  m = 0
450  IF( n.EQ.0 ) THEN
451  work( 1 ) = 1
452  RETURN
453  END IF
454 *
455  IF( n.EQ.1 ) THEN
456  work( 1 ) = 7
457  IF( alleig .OR. indeig ) THEN
458  m = 1
459  w( 1 ) = a( 1, 1 )
460  ELSE
461  IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
462  m = 1
463  w( 1 ) = a( 1, 1 )
464  END IF
465  END IF
466  IF( wantz ) THEN
467  z( 1, 1 ) = one
468  isuppz( 1 ) = 1
469  isuppz( 2 ) = 1
470  END IF
471  RETURN
472  END IF
473 *
474 * Get machine constants.
475 *
476  safmin = dlamch( 'Safe minimum' )
477  eps = dlamch( 'Precision' )
478  smlnum = safmin / eps
479  bignum = one / smlnum
480  rmin = sqrt( smlnum )
481  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
482 *
483 * Scale matrix to allowable range, if necessary.
484 *
485  iscale = 0
486  abstll = abstol
487  IF (valeig) THEN
488  vll = vl
489  vuu = vu
490  END IF
491  anrm = dlansy( 'M', uplo, n, a, lda, work )
492  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
493  iscale = 1
494  sigma = rmin / anrm
495  ELSE IF( anrm.GT.rmax ) THEN
496  iscale = 1
497  sigma = rmax / anrm
498  END IF
499  IF( iscale.EQ.1 ) THEN
500  IF( lower ) THEN
501  DO 10 j = 1, n
502  CALL dscal( n-j+1, sigma, a( j, j ), 1 )
503  10 CONTINUE
504  ELSE
505  DO 20 j = 1, n
506  CALL dscal( j, sigma, a( 1, j ), 1 )
507  20 CONTINUE
508  END IF
509  IF( abstol.GT.0 )
510  $ abstll = abstol*sigma
511  IF( valeig ) THEN
512  vll = vl*sigma
513  vuu = vu*sigma
514  END IF
515  END IF
516 
517 * Initialize indices into workspaces. Note: The IWORK indices are
518 * used only if DSTERF or DSTEMR fail.
519 
520 * WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
521 * elementary reflectors used in DSYTRD.
522  indtau = 1
523 * WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
524  indd = indtau + n
525 * WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
526 * tridiagonal matrix from DSYTRD.
527  inde = indd + n
528 * WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
529 * -written by DSTEMR (the DSTERF path copies the diagonal to W).
530  inddd = inde + n
531 * WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
532 * -written while computing the eigenvalues in DSTERF and DSTEMR.
533  indee = inddd + n
534 * INDWK is the starting offset of the left-over workspace, and
535 * LLWORK is the remaining workspace size.
536  indwk = indee + n
537  llwork = lwork - indwk + 1
538 
539 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
540 * stores the block indices of each of the M<=N eigenvalues.
541  indibl = 1
542 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
543 * stores the starting and finishing indices of each block.
544  indisp = indibl + n
545 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
546 * that corresponding to eigenvectors that fail to converge in
547 * DSTEIN. This information is discarded; if any fail, the driver
548 * returns INFO > 0.
549  indifl = indisp + n
550 * INDIWO is the offset of the remaining integer workspace.
551  indiwo = indifl + n
552 
553 *
554 * Call DSYTRD to reduce symmetric matrix to tridiagonal form.
555 *
556  CALL dsytrd( uplo, n, a, lda, work( indd ), work( inde ),
557  $ work( indtau ), work( indwk ), llwork, iinfo )
558 *
559 * If all eigenvalues are desired
560 * then call DSTERF or DSTEMR and DORMTR.
561 *
562  IF( ( alleig .OR. ( indeig .AND. il.EQ.1 .AND. iu.EQ.n ) ) .AND.
563  $ ieeeok.EQ.1 ) THEN
564  IF( .NOT.wantz ) THEN
565  CALL dcopy( n, work( indd ), 1, w, 1 )
566  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
567  CALL dsterf( n, w, work( indee ), info )
568  ELSE
569  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
570  CALL dcopy( n, work( indd ), 1, work( inddd ), 1 )
571 *
572  IF (abstol .LE. two*n*eps) THEN
573  tryrac = .true.
574  ELSE
575  tryrac = .false.
576  END IF
577  CALL dstemr( jobz, 'A', n, work( inddd ), work( indee ),
578  $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
579  $ tryrac, work( indwk ), lwork, iwork, liwork,
580  $ info )
581 *
582 *
583 *
584 * Apply orthogonal matrix used in reduction to tridiagonal
585 * form to eigenvectors returned by DSTEMR.
586 *
587  IF( wantz .AND. info.EQ.0 ) THEN
588  indwkn = inde
589  llwrkn = lwork - indwkn + 1
590  CALL dormtr( 'L', uplo, 'N', n, m, a, lda,
591  $ work( indtau ), z, ldz, work( indwkn ),
592  $ llwrkn, iinfo )
593  END IF
594  END IF
595 *
596 *
597  IF( info.EQ.0 ) THEN
598 * Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
599 * undefined.
600  m = n
601  GO TO 30
602  END IF
603  info = 0
604  END IF
605 *
606 * Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
607 * Also call DSTEBZ and DSTEIN if DSTEMR fails.
608 *
609  IF( wantz ) THEN
610  order = 'B'
611  ELSE
612  order = 'E'
613  END IF
614 
615  CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
616  $ work( indd ), work( inde ), m, nsplit, w,
617  $ iwork( indibl ), iwork( indisp ), work( indwk ),
618  $ iwork( indiwo ), info )
619 *
620  IF( wantz ) THEN
621  CALL dstein( n, work( indd ), work( inde ), m, w,
622  $ iwork( indibl ), iwork( indisp ), z, ldz,
623  $ work( indwk ), iwork( indiwo ), iwork( indifl ),
624  $ info )
625 *
626 * Apply orthogonal matrix used in reduction to tridiagonal
627 * form to eigenvectors returned by DSTEIN.
628 *
629  indwkn = inde
630  llwrkn = lwork - indwkn + 1
631  CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
632  $ ldz, work( indwkn ), llwrkn, iinfo )
633  END IF
634 *
635 * If matrix was scaled, then rescale eigenvalues appropriately.
636 *
637 * Jump here if DSTEMR/DSTEIN succeeded.
638  30 CONTINUE
639  IF( iscale.EQ.1 ) THEN
640  IF( info.EQ.0 ) THEN
641  imax = m
642  ELSE
643  imax = info - 1
644  END IF
645  CALL dscal( imax, one / sigma, w, 1 )
646  END IF
647 *
648 * If eigenvalues are not in order, then sort them, along with
649 * eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
650 * It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
651 * not return this detailed information to the user.
652 *
653  IF( wantz ) THEN
654  DO 50 j = 1, m - 1
655  i = 0
656  tmp1 = w( j )
657  DO 40 jj = j + 1, m
658  IF( w( jj ).LT.tmp1 ) THEN
659  i = jj
660  tmp1 = w( jj )
661  END IF
662  40 CONTINUE
663 *
664  IF( i.NE.0 ) THEN
665  w( i ) = w( j )
666  w( j ) = tmp1
667  CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
668  END IF
669  50 CONTINUE
670  END IF
671 *
672 * Set WORK(1) to optimal workspace size.
673 *
674  work( 1 ) = lwkopt
675  iwork( 1 ) = liwmin
676 *
677  RETURN
678 *
679 * End of DSYEVR
680 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
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, or the element of largest absolute value of a real symmetric matrix.
Definition: dlansy.f:124
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dsytrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
DSYTRD
Definition: dsytrd.f:194
subroutine dstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
DSTEIN
Definition: dstein.f:176
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
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:323
subroutine dormtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMTR
Definition: dormtr.f:173

Here is the call graph for this function:

Here is the caller graph for this function: