LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ dsyevr()

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 DSTEIN 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.
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 331 of file dsyevr.f.

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