LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dstevr ( character  JOBZ,
character  RANGE,
integer  N,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
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 
)

DSTEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

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

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

 Whenever possible, DSTEVR 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 the i-th
 unreduced block of T,
    (a) Compute T - sigma_i = L_i D_i L_i^T, such that L_i D_i L_i^T
         is a relatively robust representation,
    (b) Compute the eigenvalues, lambda_j, of L_i D_i L_i^T to high
        relative accuracy by the dqds algorithm,
    (c) If there is a cluster of close eigenvalues, "choose" sigma_i
        close to the cluster, and go to step (a),
    (d) Given the approximate eigenvalue lambda_j of L_i D_i L_i^T,
        compute the corresponding eigenvector by forming a
        rank-revealing twisted factorization.
 The desired accuracy of the output can be specified by the input
 parameter ABSTOL.

 For more details, see "A new O(n^2) algorithm for the symmetric
 tridiagonal eigenvalue/eigenvector problem", by Inderjit Dhillon,
 Computer Science Division Technical Report No. UCB//CSD-97-971,
 UC Berkeley, May 1997.


 Note 1 : DSTEVR calls DSTEMR when the full spectrum is requested
 on machines which conform to the ieee-754 floating point standard.
 DSTEVR 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]N
          N is INTEGER
          The order of the matrix.  N >= 0.
[in,out]D
          D is DOUBLE PRECISION array, dimension (N)
          On entry, the n diagonal elements of the tridiagonal matrix
          A.
          On exit, D may be multiplied by a constant factor chosen
          to avoid over/underflow in computing the eigenvalues.
[in,out]E
          E is DOUBLE PRECISION array, dimension (max(1,N-1))
          On entry, the (n-1) subdiagonal elements of the tridiagonal
          matrix A in elements 1 to N-1 of E.
          On exit, E may be multiplied by a constant factor chosen
          to avoid over/underflow in computing the eigenvalues.
[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).
          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.
[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 (and
          minimal) LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.  LWORK >= max(1,20*N).

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK or
          LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal (and
          minimal) LIWORK.
[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 sizes of the WORK and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK or 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

Definition at line 306 of file dstevr.f.

306 *
307 * -- LAPACK driver routine (version 3.6.1) --
308 * -- LAPACK is a software package provided by Univ. of Tennessee, --
309 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
310 * June 2016
311 *
312 * .. Scalar Arguments ..
313  CHARACTER jobz, range
314  INTEGER il, info, iu, ldz, liwork, lwork, m, n
315  DOUBLE PRECISION abstol, vl, vu
316 * ..
317 * .. Array Arguments ..
318  INTEGER isuppz( * ), iwork( * )
319  DOUBLE PRECISION d( * ), e( * ), w( * ), work( * ), z( ldz, * )
320 * ..
321 *
322 * =====================================================================
323 *
324 * .. Parameters ..
325  DOUBLE PRECISION zero, one, two
326  parameter ( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
327 * ..
328 * .. Local Scalars ..
329  LOGICAL alleig, indeig, test, lquery, valeig, wantz,
330  $ tryrac
331  CHARACTER order
332  INTEGER i, ieeeok, imax, indibl, indifl, indisp,
333  $ indiwo, iscale, itmp1, j, jj, liwmin, lwmin,
334  $ nsplit
335  DOUBLE PRECISION bignum, eps, rmax, rmin, safmin, sigma, smlnum,
336  $ tmp1, tnrm, vll, vuu
337 * ..
338 * .. External Functions ..
339  LOGICAL lsame
340  INTEGER ilaenv
341  DOUBLE PRECISION dlamch, dlanst
342  EXTERNAL lsame, ilaenv, dlamch, dlanst
343 * ..
344 * .. External Subroutines ..
345  EXTERNAL dcopy, dscal, dstebz, dstemr, dstein, dsterf,
346  $ dswap, xerbla
347 * ..
348 * .. Intrinsic Functions ..
349  INTRINSIC max, min, sqrt
350 * ..
351 * .. Executable Statements ..
352 *
353 *
354 * Test the input parameters.
355 *
356  ieeeok = ilaenv( 10, 'DSTEVR', 'N', 1, 2, 3, 4 )
357 *
358  wantz = lsame( jobz, 'V' )
359  alleig = lsame( range, 'A' )
360  valeig = lsame( range, 'V' )
361  indeig = lsame( range, 'I' )
362 *
363  lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
364  lwmin = max( 1, 20*n )
365  liwmin = max( 1, 10*n )
366 *
367 *
368  info = 0
369  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
370  info = -1
371  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
372  info = -2
373  ELSE IF( n.LT.0 ) THEN
374  info = -3
375  ELSE
376  IF( valeig ) THEN
377  IF( n.GT.0 .AND. vu.LE.vl )
378  $ info = -7
379  ELSE IF( indeig ) THEN
380  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
381  info = -8
382  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
383  info = -9
384  END IF
385  END IF
386  END IF
387  IF( info.EQ.0 ) THEN
388  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
389  info = -14
390  END IF
391  END IF
392 *
393  IF( info.EQ.0 ) THEN
394  work( 1 ) = lwmin
395  iwork( 1 ) = liwmin
396 *
397  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
398  info = -17
399  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
400  info = -19
401  END IF
402  END IF
403 *
404  IF( info.NE.0 ) THEN
405  CALL xerbla( 'DSTEVR', -info )
406  RETURN
407  ELSE IF( lquery ) THEN
408  RETURN
409  END IF
410 *
411 * Quick return if possible
412 *
413  m = 0
414  IF( n.EQ.0 )
415  $ RETURN
416 *
417  IF( n.EQ.1 ) THEN
418  IF( alleig .OR. indeig ) THEN
419  m = 1
420  w( 1 ) = d( 1 )
421  ELSE
422  IF( vl.LT.d( 1 ) .AND. vu.GE.d( 1 ) ) THEN
423  m = 1
424  w( 1 ) = d( 1 )
425  END IF
426  END IF
427  IF( wantz )
428  $ z( 1, 1 ) = one
429  RETURN
430  END IF
431 *
432 * Get machine constants.
433 *
434  safmin = dlamch( 'Safe minimum' )
435  eps = dlamch( 'Precision' )
436  smlnum = safmin / eps
437  bignum = one / smlnum
438  rmin = sqrt( smlnum )
439  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
440 *
441 *
442 * Scale matrix to allowable range, if necessary.
443 *
444  iscale = 0
445  IF( valeig ) THEN
446  vll = vl
447  vuu = vu
448  END IF
449 *
450  tnrm = dlanst( 'M', n, d, e )
451  IF( tnrm.GT.zero .AND. tnrm.LT.rmin ) THEN
452  iscale = 1
453  sigma = rmin / tnrm
454  ELSE IF( tnrm.GT.rmax ) THEN
455  iscale = 1
456  sigma = rmax / tnrm
457  END IF
458  IF( iscale.EQ.1 ) THEN
459  CALL dscal( n, sigma, d, 1 )
460  CALL dscal( n-1, sigma, e( 1 ), 1 )
461  IF( valeig ) THEN
462  vll = vl*sigma
463  vuu = vu*sigma
464  END IF
465  END IF
466 
467 * Initialize indices into workspaces. Note: These indices are used only
468 * if DSTERF or DSTEMR fail.
469 
470 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
471 * stores the block indices of each of the M<=N eigenvalues.
472  indibl = 1
473 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
474 * stores the starting and finishing indices of each block.
475  indisp = indibl + n
476 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
477 * that corresponding to eigenvectors that fail to converge in
478 * DSTEIN. This information is discarded; if any fail, the driver
479 * returns INFO > 0.
480  indifl = indisp + n
481 * INDIWO is the offset of the remaining integer workspace.
482  indiwo = indisp + n
483 *
484 * If all eigenvalues are desired, then
485 * call DSTERF or DSTEMR. If this fails for some eigenvalue, then
486 * try DSTEBZ.
487 *
488 *
489  test = .false.
490  IF( indeig ) THEN
491  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
492  test = .true.
493  END IF
494  END IF
495  IF( ( alleig .OR. test ) .AND. ieeeok.EQ.1 ) THEN
496  CALL dcopy( n-1, e( 1 ), 1, work( 1 ), 1 )
497  IF( .NOT.wantz ) THEN
498  CALL dcopy( n, d, 1, w, 1 )
499  CALL dsterf( n, w, work, info )
500  ELSE
501  CALL dcopy( n, d, 1, work( n+1 ), 1 )
502  IF (abstol .LE. two*n*eps) THEN
503  tryrac = .true.
504  ELSE
505  tryrac = .false.
506  END IF
507  CALL dstemr( jobz, 'A', n, work( n+1 ), work, vl, vu, il,
508  $ iu, m, w, z, ldz, n, isuppz, tryrac,
509  $ work( 2*n+1 ), lwork-2*n, iwork, liwork, info )
510 *
511  END IF
512  IF( info.EQ.0 ) THEN
513  m = n
514  GO TO 10
515  END IF
516  info = 0
517  END IF
518 *
519 * Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
520 *
521  IF( wantz ) THEN
522  order = 'B'
523  ELSE
524  order = 'E'
525  END IF
526 
527  CALL dstebz( range, order, n, vll, vuu, il, iu, abstol, d, e, m,
528  $ nsplit, w, iwork( indibl ), iwork( indisp ), work,
529  $ iwork( indiwo ), info )
530 *
531  IF( wantz ) THEN
532  CALL dstein( n, d, e, m, w, iwork( indibl ), iwork( indisp ),
533  $ z, ldz, work, iwork( indiwo ), iwork( indifl ),
534  $ info )
535  END IF
536 *
537 * If matrix was scaled, then rescale eigenvalues appropriately.
538 *
539  10 CONTINUE
540  IF( iscale.EQ.1 ) THEN
541  IF( info.EQ.0 ) THEN
542  imax = m
543  ELSE
544  imax = info - 1
545  END IF
546  CALL dscal( imax, one / sigma, w, 1 )
547  END IF
548 *
549 * If eigenvalues are not in order, then sort them, along with
550 * eigenvectors.
551 *
552  IF( wantz ) THEN
553  DO 30 j = 1, m - 1
554  i = 0
555  tmp1 = w( j )
556  DO 20 jj = j + 1, m
557  IF( w( jj ).LT.tmp1 ) THEN
558  i = jj
559  tmp1 = w( jj )
560  END IF
561  20 CONTINUE
562 *
563  IF( i.NE.0 ) THEN
564  itmp1 = iwork( i )
565  w( i ) = w( j )
566  iwork( i ) = iwork( j )
567  w( j ) = tmp1
568  iwork( j ) = itmp1
569  CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
570  END IF
571  30 CONTINUE
572  END IF
573 *
574 * Causes problems with tests 19 & 20:
575 * IF (wantz .and. INDEIG ) Z( 1,1) = Z(1,1) / 1.002 + .002
576 *
577 *
578  work( 1 ) = lwmin
579  iwork( 1 ) = liwmin
580  RETURN
581 *
582 * End of DSTEVR
583 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
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 xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
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
double precision function dlanst(NORM, N, D, E)
DLANST 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 tridiagonal matrix.
Definition: dlanst.f:102
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: