LAPACK  3.6.1
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 ).
          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 334 of file dsyevr.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: