File:  [local] / rpl / lapack / lapack / dsyevr.f
Revision 1.17: download - view: text, annotated - select for diffs - revision graph
Sat Jun 17 10:54:04 2017 UTC (6 years, 11 months ago) by bertrand
Branches: MAIN
CVS tags: HEAD
Mise à jour de lapack.

    1: *> \brief <b> DSYEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices</b>
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DSYEVR + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyevr.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyevr.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyevr.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DSYEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
   22: *                          ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK,
   23: *                          IWORK, LIWORK, INFO )
   24: *
   25: *       .. Scalar Arguments ..
   26: *       CHARACTER          JOBZ, RANGE, UPLO
   27: *       INTEGER            IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
   28: *       DOUBLE PRECISION   ABSTOL, VL, VU
   29: *       ..
   30: *       .. Array Arguments ..
   31: *       INTEGER            ISUPPZ( * ), IWORK( * )
   32: *       DOUBLE PRECISION   A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
   33: *       ..
   34: *
   35: *
   36: *> \par Purpose:
   37: *  =============
   38: *>
   39: *> \verbatim
   40: *>
   41: *> DSYEVR computes selected eigenvalues and, optionally, eigenvectors
   42: *> of a real symmetric matrix A.  Eigenvalues and eigenvectors can be
   43: *> selected by specifying either a range of values or a range of
   44: *> indices for the desired eigenvalues.
   45: *>
   46: *> DSYEVR first reduces the matrix A to tridiagonal form T with a call
   47: *> to DSYTRD.  Then, whenever possible, DSYEVR calls DSTEMR to compute
   48: *> the eigenspectrum using Relatively Robust Representations.  DSTEMR
   49: *> computes eigenvalues by the dqds algorithm, while orthogonal
   50: *> eigenvectors are computed from various "good" L D L^T representations
   51: *> (also known as Relatively Robust Representations). Gram-Schmidt
   52: *> orthogonalization is avoided as far as possible. More specifically,
   53: *> the various steps of the algorithm are as follows.
   54: *>
   55: *> For each unreduced block (submatrix) of T,
   56: *>    (a) Compute T - sigma I  = L D L^T, so that L and D
   57: *>        define all the wanted eigenvalues to high relative accuracy.
   58: *>        This means that small relative changes in the entries of D and L
   59: *>        cause only small relative changes in the eigenvalues and
   60: *>        eigenvectors. The standard (unfactored) representation of the
   61: *>        tridiagonal matrix T does not have this property in general.
   62: *>    (b) Compute the eigenvalues to suitable accuracy.
   63: *>        If the eigenvectors are desired, the algorithm attains full
   64: *>        accuracy of the computed eigenvalues only right before
   65: *>        the corresponding vectors have to be computed, see steps c) and d).
   66: *>    (c) For each cluster of close eigenvalues, select a new
   67: *>        shift close to the cluster, find a new factorization, and refine
   68: *>        the shifted eigenvalues to suitable accuracy.
   69: *>    (d) For each eigenvalue with a large enough relative separation compute
   70: *>        the corresponding eigenvector by forming a rank revealing twisted
   71: *>        factorization. Go back to (c) for any clusters that remain.
   72: *>
   73: *> The desired accuracy of the output can be specified by the input
   74: *> parameter ABSTOL.
   75: *>
   76: *> For more details, see DSTEMR's documentation and:
   77: *> - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
   78: *>   to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
   79: *>   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
   80: *> - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
   81: *>   Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
   82: *>   2004.  Also LAPACK Working Note 154.
   83: *> - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
   84: *>   tridiagonal eigenvalue/eigenvector problem",
   85: *>   Computer Science Division Technical Report No. UCB/CSD-97-971,
   86: *>   UC Berkeley, May 1997.
   87: *>
   88: *>
   89: *> Note 1 : DSYEVR calls DSTEMR when the full spectrum is requested
   90: *> on machines which conform to the ieee-754 floating point standard.
   91: *> DSYEVR calls DSTEBZ and SSTEIN on non-ieee machines and
   92: *> when partial spectrum requests are made.
   93: *>
   94: *> Normal execution of DSTEMR may create NaNs and infinities and
   95: *> hence may abort due to a floating point exception in environments
   96: *> which do not handle NaNs and infinities in the ieee standard default
   97: *> manner.
   98: *> \endverbatim
   99: *
  100: *  Arguments:
  101: *  ==========
  102: *
  103: *> \param[in] JOBZ
  104: *> \verbatim
  105: *>          JOBZ is CHARACTER*1
  106: *>          = 'N':  Compute eigenvalues only;
  107: *>          = 'V':  Compute eigenvalues and eigenvectors.
  108: *> \endverbatim
  109: *>
  110: *> \param[in] RANGE
  111: *> \verbatim
  112: *>          RANGE is CHARACTER*1
  113: *>          = 'A': all eigenvalues will be found.
  114: *>          = 'V': all eigenvalues in the half-open interval (VL,VU]
  115: *>                 will be found.
  116: *>          = 'I': the IL-th through IU-th eigenvalues will be found.
  117: *>          For RANGE = 'V' or 'I' and IU - IL < N - 1, DSTEBZ and
  118: *>          DSTEIN are called
  119: *> \endverbatim
  120: *>
  121: *> \param[in] UPLO
  122: *> \verbatim
  123: *>          UPLO is CHARACTER*1
  124: *>          = 'U':  Upper triangle of A is stored;
  125: *>          = 'L':  Lower triangle of A is stored.
  126: *> \endverbatim
  127: *>
  128: *> \param[in] N
  129: *> \verbatim
  130: *>          N is INTEGER
  131: *>          The order of the matrix A.  N >= 0.
  132: *> \endverbatim
  133: *>
  134: *> \param[in,out] A
  135: *> \verbatim
  136: *>          A is DOUBLE PRECISION array, dimension (LDA, N)
  137: *>          On entry, the symmetric matrix A.  If UPLO = 'U', the
  138: *>          leading N-by-N upper triangular part of A contains the
  139: *>          upper triangular part of the matrix A.  If UPLO = 'L',
  140: *>          the leading N-by-N lower triangular part of A contains
  141: *>          the lower triangular part of the matrix A.
  142: *>          On exit, the lower triangle (if UPLO='L') or the upper
  143: *>          triangle (if UPLO='U') of A, including the diagonal, is
  144: *>          destroyed.
  145: *> \endverbatim
  146: *>
  147: *> \param[in] LDA
  148: *> \verbatim
  149: *>          LDA is INTEGER
  150: *>          The leading dimension of the array A.  LDA >= max(1,N).
  151: *> \endverbatim
  152: *>
  153: *> \param[in] VL
  154: *> \verbatim
  155: *>          VL is DOUBLE PRECISION
  156: *>          If RANGE='V', the lower bound of the interval to
  157: *>          be searched for eigenvalues. VL < VU.
  158: *>          Not referenced if RANGE = 'A' or 'I'.
  159: *> \endverbatim
  160: *>
  161: *> \param[in] VU
  162: *> \verbatim
  163: *>          VU is DOUBLE PRECISION
  164: *>          If RANGE='V', the upper bound of the interval to
  165: *>          be searched for eigenvalues. VL < VU.
  166: *>          Not referenced if RANGE = 'A' or 'I'.
  167: *> \endverbatim
  168: *>
  169: *> \param[in] IL
  170: *> \verbatim
  171: *>          IL is INTEGER
  172: *>          If RANGE='I', the index of the
  173: *>          smallest eigenvalue to be returned.
  174: *>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
  175: *>          Not referenced if RANGE = 'A' or 'V'.
  176: *> \endverbatim
  177: *>
  178: *> \param[in] IU
  179: *> \verbatim
  180: *>          IU is INTEGER
  181: *>          If RANGE='I', the index of the
  182: *>          largest eigenvalue to be returned.
  183: *>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
  184: *>          Not referenced if RANGE = 'A' or 'V'.
  185: *> \endverbatim
  186: *>
  187: *> \param[in] ABSTOL
  188: *> \verbatim
  189: *>          ABSTOL is DOUBLE PRECISION
  190: *>          The absolute error tolerance for the eigenvalues.
  191: *>          An approximate eigenvalue is accepted as converged
  192: *>          when it is determined to lie in an interval [a,b]
  193: *>          of width less than or equal to
  194: *>
  195: *>                  ABSTOL + EPS *   max( |a|,|b| ) ,
  196: *>
  197: *>          where EPS is the machine precision.  If ABSTOL is less than
  198: *>          or equal to zero, then  EPS*|T|  will be used in its place,
  199: *>          where |T| is the 1-norm of the tridiagonal matrix obtained
  200: *>          by reducing A to tridiagonal form.
  201: *>
  202: *>          See "Computing Small Singular Values of Bidiagonal Matrices
  203: *>          with Guaranteed High Relative Accuracy," by Demmel and
  204: *>          Kahan, LAPACK Working Note #3.
  205: *>
  206: *>          If high relative accuracy is important, set ABSTOL to
  207: *>          DLAMCH( 'Safe minimum' ).  Doing so will guarantee that
  208: *>          eigenvalues are computed to high relative accuracy when
  209: *>          possible in future releases.  The current code does not
  210: *>          make any guarantees about high relative accuracy, but
  211: *>          future releases will. See J. Barlow and J. Demmel,
  212: *>          "Computing Accurate Eigensystems of Scaled Diagonally
  213: *>          Dominant Matrices", LAPACK Working Note #7, for a discussion
  214: *>          of which matrices define their eigenvalues to high relative
  215: *>          accuracy.
  216: *> \endverbatim
  217: *>
  218: *> \param[out] M
  219: *> \verbatim
  220: *>          M is INTEGER
  221: *>          The total number of eigenvalues found.  0 <= M <= N.
  222: *>          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
  223: *> \endverbatim
  224: *>
  225: *> \param[out] W
  226: *> \verbatim
  227: *>          W is DOUBLE PRECISION array, dimension (N)
  228: *>          The first M elements contain the selected eigenvalues in
  229: *>          ascending order.
  230: *> \endverbatim
  231: *>
  232: *> \param[out] Z
  233: *> \verbatim
  234: *>          Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
  235: *>          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
  236: *>          contain the orthonormal eigenvectors of the matrix A
  237: *>          corresponding to the selected eigenvalues, with the i-th
  238: *>          column of Z holding the eigenvector associated with W(i).
  239: *>          If JOBZ = 'N', then Z is not referenced.
  240: *>          Note: the user must ensure that at least max(1,M) columns are
  241: *>          supplied in the array Z; if RANGE = 'V', the exact value of M
  242: *>          is not known in advance and an upper bound must be used.
  243: *>          Supplying N columns is always safe.
  244: *> \endverbatim
  245: *>
  246: *> \param[in] LDZ
  247: *> \verbatim
  248: *>          LDZ is INTEGER
  249: *>          The leading dimension of the array Z.  LDZ >= 1, and if
  250: *>          JOBZ = 'V', LDZ >= max(1,N).
  251: *> \endverbatim
  252: *>
  253: *> \param[out] ISUPPZ
  254: *> \verbatim
  255: *>          ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
  256: *>          The support of the eigenvectors in Z, i.e., the indices
  257: *>          indicating the nonzero elements in Z. The i-th eigenvector
  258: *>          is nonzero only in elements ISUPPZ( 2*i-1 ) through
  259: *>          ISUPPZ( 2*i ). This is an output of DSTEMR (tridiagonal
  260: *>          matrix). The support of the eigenvectors of A is typically
  261: *>          1:N because of the orthogonal transformations applied by DORMTR.
  262: *>          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
  263: *> \endverbatim
  264: *>
  265: *> \param[out] WORK
  266: *> \verbatim
  267: *>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  268: *>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  269: *> \endverbatim
  270: *>
  271: *> \param[in] LWORK
  272: *> \verbatim
  273: *>          LWORK is INTEGER
  274: *>          The dimension of the array WORK.  LWORK >= max(1,26*N).
  275: *>          For optimal efficiency, LWORK >= (NB+6)*N,
  276: *>          where NB is the max of the blocksize for DSYTRD and DORMTR
  277: *>          returned by ILAENV.
  278: *>
  279: *>          If LWORK = -1, then a workspace query is assumed; the routine
  280: *>          only calculates the optimal size of the WORK array, returns
  281: *>          this value as the first entry of the WORK array, and no error
  282: *>          message related to LWORK is issued by XERBLA.
  283: *> \endverbatim
  284: *>
  285: *> \param[out] IWORK
  286: *> \verbatim
  287: *>          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
  288: *>          On exit, if INFO = 0, IWORK(1) returns the optimal LWORK.
  289: *> \endverbatim
  290: *>
  291: *> \param[in] LIWORK
  292: *> \verbatim
  293: *>          LIWORK is INTEGER
  294: *>          The dimension of the array IWORK.  LIWORK >= max(1,10*N).
  295: *>
  296: *>          If LIWORK = -1, then a workspace query is assumed; the
  297: *>          routine only calculates the optimal size of the IWORK array,
  298: *>          returns this value as the first entry of the IWORK array, and
  299: *>          no error message related to LIWORK is issued by XERBLA.
  300: *> \endverbatim
  301: *>
  302: *> \param[out] INFO
  303: *> \verbatim
  304: *>          INFO is INTEGER
  305: *>          = 0:  successful exit
  306: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
  307: *>          > 0:  Internal error
  308: *> \endverbatim
  309: *
  310: *  Authors:
  311: *  ========
  312: *
  313: *> \author Univ. of Tennessee
  314: *> \author Univ. of California Berkeley
  315: *> \author Univ. of Colorado Denver
  316: *> \author NAG Ltd.
  317: *
  318: *> \date June 2016
  319: *
  320: *> \ingroup doubleSYeigen
  321: *
  322: *> \par Contributors:
  323: *  ==================
  324: *>
  325: *>     Inderjit Dhillon, IBM Almaden, USA \n
  326: *>     Osni Marques, LBNL/NERSC, USA \n
  327: *>     Ken Stanley, Computer Science Division, University of
  328: *>       California at Berkeley, USA \n
  329: *>     Jason Riedy, Computer Science Division, University of
  330: *>       California at Berkeley, USA \n
  331: *>
  332: *  =====================================================================
  333:       SUBROUTINE DSYEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
  334:      $                   ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK,
  335:      $                   IWORK, LIWORK, INFO )
  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,
  377:      $                   DSTERF, DSWAP, DSYTRD, XERBLA
  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: *
  681:       END

CVSweb interface <joel.bertrand@systella.fr>