File:  [local] / rpl / lapack / lapack / dsyevx.f
Revision 1.13: download - view: text, annotated - select for diffs - revision graph
Sat Aug 27 15:27:11 2016 UTC (7 years, 9 months ago) by bertrand
Branches: MAIN
CVS tags: HEAD

Mise à jour de lapack.

    1: *> \brief <b> DSYEVX 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 DSYEVX + dependencies 
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyevx.f"> 
   11: *> [TGZ]</a> 
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyevx.f"> 
   13: *> [ZIP]</a> 
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyevx.f"> 
   15: *> [TXT]</a>
   16: *> \endhtmlonly 
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DSYEVX( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
   22: *                          ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK,
   23: *                          IFAIL, INFO )
   24:    25: *       .. Scalar Arguments ..
   26: *       CHARACTER          JOBZ, RANGE, UPLO
   27: *       INTEGER            IL, INFO, IU, LDA, LDZ, LWORK, M, N
   28: *       DOUBLE PRECISION   ABSTOL, VL, VU
   29: *       ..
   30: *       .. Array Arguments ..
   31: *       INTEGER            IFAIL( * ), IWORK( * )
   32: *       DOUBLE PRECISION   A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
   33: *       ..
   34: *  
   35: *
   36: *> \par Purpose:
   37: *  =============
   38: *>
   39: *> \verbatim
   40: *>
   41: *> DSYEVX 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 indices
   44: *> for the desired eigenvalues.
   45: *> \endverbatim
   46: *
   47: *  Arguments:
   48: *  ==========
   49: *
   50: *> \param[in] JOBZ
   51: *> \verbatim
   52: *>          JOBZ is CHARACTER*1
   53: *>          = 'N':  Compute eigenvalues only;
   54: *>          = 'V':  Compute eigenvalues and eigenvectors.
   55: *> \endverbatim
   56: *>
   57: *> \param[in] RANGE
   58: *> \verbatim
   59: *>          RANGE is CHARACTER*1
   60: *>          = 'A': all eigenvalues will be found.
   61: *>          = 'V': all eigenvalues in the half-open interval (VL,VU]
   62: *>                 will be found.
   63: *>          = 'I': the IL-th through IU-th eigenvalues will be found.
   64: *> \endverbatim
   65: *>
   66: *> \param[in] UPLO
   67: *> \verbatim
   68: *>          UPLO is CHARACTER*1
   69: *>          = 'U':  Upper triangle of A is stored;
   70: *>          = 'L':  Lower triangle of A is stored.
   71: *> \endverbatim
   72: *>
   73: *> \param[in] N
   74: *> \verbatim
   75: *>          N is INTEGER
   76: *>          The order of the matrix A.  N >= 0.
   77: *> \endverbatim
   78: *>
   79: *> \param[in,out] A
   80: *> \verbatim
   81: *>          A is DOUBLE PRECISION array, dimension (LDA, N)
   82: *>          On entry, the symmetric matrix A.  If UPLO = 'U', the
   83: *>          leading N-by-N upper triangular part of A contains the
   84: *>          upper triangular part of the matrix A.  If UPLO = 'L',
   85: *>          the leading N-by-N lower triangular part of A contains
   86: *>          the lower triangular part of the matrix A.
   87: *>          On exit, the lower triangle (if UPLO='L') or the upper
   88: *>          triangle (if UPLO='U') of A, including the diagonal, is
   89: *>          destroyed.
   90: *> \endverbatim
   91: *>
   92: *> \param[in] LDA
   93: *> \verbatim
   94: *>          LDA is INTEGER
   95: *>          The leading dimension of the array A.  LDA >= max(1,N).
   96: *> \endverbatim
   97: *>
   98: *> \param[in] VL
   99: *> \verbatim
  100: *>          VL is DOUBLE PRECISION
  101: *>          If RANGE='V', the lower bound of the interval to
  102: *>          be searched for eigenvalues. VL < VU.
  103: *>          Not referenced if RANGE = 'A' or 'I'.
  104: *> \endverbatim
  105: *>
  106: *> \param[in] VU
  107: *> \verbatim
  108: *>          VU is DOUBLE PRECISION
  109: *>          If RANGE='V', the upper bound of the interval to
  110: *>          be searched for eigenvalues. VL < VU.
  111: *>          Not referenced if RANGE = 'A' or 'I'.
  112: *> \endverbatim
  113: *>
  114: *> \param[in] IL
  115: *> \verbatim
  116: *>          IL is INTEGER
  117: *>          If RANGE='I', the index of the
  118: *>          smallest eigenvalue to be returned.
  119: *>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
  120: *>          Not referenced if RANGE = 'A' or 'V'.
  121: *> \endverbatim
  122: *>
  123: *> \param[in] IU
  124: *> \verbatim
  125: *>          IU is INTEGER
  126: *>          If RANGE='I', the index of the
  127: *>          largest eigenvalue to be returned.
  128: *>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
  129: *>          Not referenced if RANGE = 'A' or 'V'.
  130: *> \endverbatim
  131: *>
  132: *> \param[in] ABSTOL
  133: *> \verbatim
  134: *>          ABSTOL is DOUBLE PRECISION
  135: *>          The absolute error tolerance for the eigenvalues.
  136: *>          An approximate eigenvalue is accepted as converged
  137: *>          when it is determined to lie in an interval [a,b]
  138: *>          of width less than or equal to
  139: *>
  140: *>                  ABSTOL + EPS *   max( |a|,|b| ) ,
  141: *>
  142: *>          where EPS is the machine precision.  If ABSTOL is less than
  143: *>          or equal to zero, then  EPS*|T|  will be used in its place,
  144: *>          where |T| is the 1-norm of the tridiagonal matrix obtained
  145: *>          by reducing A to tridiagonal form.
  146: *>
  147: *>          Eigenvalues will be computed most accurately when ABSTOL is
  148: *>          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
  149: *>          If this routine returns with INFO>0, indicating that some
  150: *>          eigenvectors did not converge, try setting ABSTOL to
  151: *>          2*DLAMCH('S').
  152: *>
  153: *>          See "Computing Small Singular Values of Bidiagonal Matrices
  154: *>          with Guaranteed High Relative Accuracy," by Demmel and
  155: *>          Kahan, LAPACK Working Note #3.
  156: *> \endverbatim
  157: *>
  158: *> \param[out] M
  159: *> \verbatim
  160: *>          M is INTEGER
  161: *>          The total number of eigenvalues found.  0 <= M <= N.
  162: *>          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
  163: *> \endverbatim
  164: *>
  165: *> \param[out] W
  166: *> \verbatim
  167: *>          W is DOUBLE PRECISION array, dimension (N)
  168: *>          On normal exit, the first M elements contain the selected
  169: *>          eigenvalues in ascending order.
  170: *> \endverbatim
  171: *>
  172: *> \param[out] Z
  173: *> \verbatim
  174: *>          Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
  175: *>          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
  176: *>          contain the orthonormal eigenvectors of the matrix A
  177: *>          corresponding to the selected eigenvalues, with the i-th
  178: *>          column of Z holding the eigenvector associated with W(i).
  179: *>          If an eigenvector fails to converge, then that column of Z
  180: *>          contains the latest approximation to the eigenvector, and the
  181: *>          index of the eigenvector is returned in IFAIL.
  182: *>          If JOBZ = 'N', then Z is not referenced.
  183: *>          Note: the user must ensure that at least max(1,M) columns are
  184: *>          supplied in the array Z; if RANGE = 'V', the exact value of M
  185: *>          is not known in advance and an upper bound must be used.
  186: *> \endverbatim
  187: *>
  188: *> \param[in] LDZ
  189: *> \verbatim
  190: *>          LDZ is INTEGER
  191: *>          The leading dimension of the array Z.  LDZ >= 1, and if
  192: *>          JOBZ = 'V', LDZ >= max(1,N).
  193: *> \endverbatim
  194: *>
  195: *> \param[out] WORK
  196: *> \verbatim
  197: *>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  198: *>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  199: *> \endverbatim
  200: *>
  201: *> \param[in] LWORK
  202: *> \verbatim
  203: *>          LWORK is INTEGER
  204: *>          The length of the array WORK.  LWORK >= 1, when N <= 1;
  205: *>          otherwise 8*N.
  206: *>          For optimal efficiency, LWORK >= (NB+3)*N,
  207: *>          where NB is the max of the blocksize for DSYTRD and DORMTR
  208: *>          returned by ILAENV.
  209: *>
  210: *>          If LWORK = -1, then a workspace query is assumed; the routine
  211: *>          only calculates the optimal size of the WORK array, returns
  212: *>          this value as the first entry of the WORK array, and no error
  213: *>          message related to LWORK is issued by XERBLA.
  214: *> \endverbatim
  215: *>
  216: *> \param[out] IWORK
  217: *> \verbatim
  218: *>          IWORK is INTEGER array, dimension (5*N)
  219: *> \endverbatim
  220: *>
  221: *> \param[out] IFAIL
  222: *> \verbatim
  223: *>          IFAIL is INTEGER array, dimension (N)
  224: *>          If JOBZ = 'V', then if INFO = 0, the first M elements of
  225: *>          IFAIL are zero.  If INFO > 0, then IFAIL contains the
  226: *>          indices of the eigenvectors that failed to converge.
  227: *>          If JOBZ = 'N', then IFAIL is not referenced.
  228: *> \endverbatim
  229: *>
  230: *> \param[out] INFO
  231: *> \verbatim
  232: *>          INFO is INTEGER
  233: *>          = 0:  successful exit
  234: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
  235: *>          > 0:  if INFO = i, then i eigenvectors failed to converge.
  236: *>                Their indices are stored in array IFAIL.
  237: *> \endverbatim
  238: *
  239: *  Authors:
  240: *  ========
  241: *
  242: *> \author Univ. of Tennessee 
  243: *> \author Univ. of California Berkeley 
  244: *> \author Univ. of Colorado Denver 
  245: *> \author NAG Ltd. 
  246: *
  247: *> \date June 2016
  248: *
  249: *> \ingroup doubleSYeigen
  250: *
  251: *  =====================================================================
  252:       SUBROUTINE DSYEVX( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
  253:      $                   ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK,
  254:      $                   IFAIL, INFO )
  255: *
  256: *  -- LAPACK driver routine (version 3.6.1) --
  257: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  258: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  259: *     June 2016
  260: *
  261: *     .. Scalar Arguments ..
  262:       CHARACTER          JOBZ, RANGE, UPLO
  263:       INTEGER            IL, INFO, IU, LDA, LDZ, LWORK, M, N
  264:       DOUBLE PRECISION   ABSTOL, VL, VU
  265: *     ..
  266: *     .. Array Arguments ..
  267:       INTEGER            IFAIL( * ), IWORK( * )
  268:       DOUBLE PRECISION   A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
  269: *     ..
  270: *
  271: * =====================================================================
  272: *
  273: *     .. Parameters ..
  274:       DOUBLE PRECISION   ZERO, ONE
  275:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
  276: *     ..
  277: *     .. Local Scalars ..
  278:       LOGICAL            ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
  279:      $                   WANTZ
  280:       CHARACTER          ORDER
  281:       INTEGER            I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
  282:      $                   INDISP, INDIWO, INDTAU, INDWKN, INDWRK, ISCALE,
  283:      $                   ITMP1, J, JJ, LLWORK, LLWRKN, LWKMIN,
  284:      $                   LWKOPT, NB, NSPLIT
  285:       DOUBLE PRECISION   ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
  286:      $                   SIGMA, SMLNUM, TMP1, VLL, VUU
  287: *     ..
  288: *     .. External Functions ..
  289:       LOGICAL            LSAME
  290:       INTEGER            ILAENV
  291:       DOUBLE PRECISION   DLAMCH, DLANSY
  292:       EXTERNAL           LSAME, ILAENV, DLAMCH, DLANSY
  293: *     ..
  294: *     .. External Subroutines ..
  295:       EXTERNAL           DCOPY, DLACPY, DORGTR, DORMTR, DSCAL, DSTEBZ,
  296:      $                   DSTEIN, DSTEQR, DSTERF, DSWAP, DSYTRD, XERBLA
  297: *     ..
  298: *     .. Intrinsic Functions ..
  299:       INTRINSIC          MAX, MIN, SQRT
  300: *     ..
  301: *     .. Executable Statements ..
  302: *
  303: *     Test the input parameters.
  304: *
  305:       LOWER = LSAME( UPLO, 'L' )
  306:       WANTZ = LSAME( JOBZ, 'V' )
  307:       ALLEIG = LSAME( RANGE, 'A' )
  308:       VALEIG = LSAME( RANGE, 'V' )
  309:       INDEIG = LSAME( RANGE, 'I' )
  310:       LQUERY = ( LWORK.EQ.-1 )
  311: *
  312:       INFO = 0
  313:       IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
  314:          INFO = -1
  315:       ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
  316:          INFO = -2
  317:       ELSE IF( .NOT.( LOWER .OR. LSAME( UPLO, 'U' ) ) ) THEN
  318:          INFO = -3
  319:       ELSE IF( N.LT.0 ) THEN
  320:          INFO = -4
  321:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
  322:          INFO = -6
  323:       ELSE
  324:          IF( VALEIG ) THEN
  325:             IF( N.GT.0 .AND. VU.LE.VL )
  326:      $         INFO = -8
  327:          ELSE IF( INDEIG ) THEN
  328:             IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
  329:                INFO = -9
  330:             ELSE IF( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
  331:                INFO = -10
  332:             END IF
  333:          END IF
  334:       END IF
  335:       IF( INFO.EQ.0 ) THEN
  336:          IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
  337:             INFO = -15
  338:          END IF
  339:       END IF
  340: *
  341:       IF( INFO.EQ.0 ) THEN
  342:          IF( N.LE.1 ) THEN
  343:             LWKMIN = 1
  344:             WORK( 1 ) = LWKMIN
  345:          ELSE
  346:             LWKMIN = 8*N
  347:             NB = ILAENV( 1, 'DSYTRD', UPLO, N, -1, -1, -1 )
  348:             NB = MAX( NB, ILAENV( 1, 'DORMTR', UPLO, N, -1, -1, -1 ) )
  349:             LWKOPT = MAX( LWKMIN, ( NB + 3 )*N )
  350:             WORK( 1 ) = LWKOPT
  351:          END IF
  352: *
  353:          IF( LWORK.LT.LWKMIN .AND. .NOT.LQUERY )
  354:      $      INFO = -17
  355:       END IF
  356: *
  357:       IF( INFO.NE.0 ) THEN
  358:          CALL XERBLA( 'DSYEVX', -INFO )
  359:          RETURN
  360:       ELSE IF( LQUERY ) THEN
  361:          RETURN
  362:       END IF
  363: *
  364: *     Quick return if possible
  365: *
  366:       M = 0
  367:       IF( N.EQ.0 ) THEN
  368:          RETURN
  369:       END IF
  370: *
  371:       IF( N.EQ.1 ) THEN
  372:          IF( ALLEIG .OR. INDEIG ) THEN
  373:             M = 1
  374:             W( 1 ) = A( 1, 1 )
  375:          ELSE
  376:             IF( VL.LT.A( 1, 1 ) .AND. VU.GE.A( 1, 1 ) ) THEN
  377:                M = 1
  378:                W( 1 ) = A( 1, 1 )
  379:             END IF
  380:          END IF
  381:          IF( WANTZ )
  382:      $      Z( 1, 1 ) = ONE
  383:          RETURN
  384:       END IF
  385: *
  386: *     Get machine constants.
  387: *
  388:       SAFMIN = DLAMCH( 'Safe minimum' )
  389:       EPS = DLAMCH( 'Precision' )
  390:       SMLNUM = SAFMIN / EPS
  391:       BIGNUM = ONE / SMLNUM
  392:       RMIN = SQRT( SMLNUM )
  393:       RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) )
  394: *
  395: *     Scale matrix to allowable range, if necessary.
  396: *
  397:       ISCALE = 0
  398:       ABSTLL = ABSTOL
  399:       IF( VALEIG ) THEN
  400:          VLL = VL
  401:          VUU = VU
  402:       END IF
  403:       ANRM = DLANSY( 'M', UPLO, N, A, LDA, WORK )
  404:       IF( ANRM.GT.ZERO .AND. ANRM.LT.RMIN ) THEN
  405:          ISCALE = 1
  406:          SIGMA = RMIN / ANRM
  407:       ELSE IF( ANRM.GT.RMAX ) THEN
  408:          ISCALE = 1
  409:          SIGMA = RMAX / ANRM
  410:       END IF
  411:       IF( ISCALE.EQ.1 ) THEN
  412:          IF( LOWER ) THEN
  413:             DO 10 J = 1, N
  414:                CALL DSCAL( N-J+1, SIGMA, A( J, J ), 1 )
  415:    10       CONTINUE
  416:          ELSE
  417:             DO 20 J = 1, N
  418:                CALL DSCAL( J, SIGMA, A( 1, J ), 1 )
  419:    20       CONTINUE
  420:          END IF
  421:          IF( ABSTOL.GT.0 )
  422:      $      ABSTLL = ABSTOL*SIGMA
  423:          IF( VALEIG ) THEN
  424:             VLL = VL*SIGMA
  425:             VUU = VU*SIGMA
  426:          END IF
  427:       END IF
  428: *
  429: *     Call DSYTRD to reduce symmetric matrix to tridiagonal form.
  430: *
  431:       INDTAU = 1
  432:       INDE = INDTAU + N
  433:       INDD = INDE + N
  434:       INDWRK = INDD + N
  435:       LLWORK = LWORK - INDWRK + 1
  436:       CALL DSYTRD( UPLO, N, A, LDA, WORK( INDD ), WORK( INDE ),
  437:      $             WORK( INDTAU ), WORK( INDWRK ), LLWORK, IINFO )
  438: *
  439: *     If all eigenvalues are desired and ABSTOL is less than or equal to
  440: *     zero, then call DSTERF or DORGTR and SSTEQR.  If this fails for
  441: *     some eigenvalue, then try DSTEBZ.
  442: *
  443:       TEST = .FALSE.
  444:       IF( INDEIG ) THEN
  445:          IF( IL.EQ.1 .AND. IU.EQ.N ) THEN
  446:             TEST = .TRUE.
  447:          END IF
  448:       END IF
  449:       IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN
  450:          CALL DCOPY( N, WORK( INDD ), 1, W, 1 )
  451:          INDEE = INDWRK + 2*N
  452:          IF( .NOT.WANTZ ) THEN
  453:             CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
  454:             CALL DSTERF( N, W, WORK( INDEE ), INFO )
  455:          ELSE
  456:             CALL DLACPY( 'A', N, N, A, LDA, Z, LDZ )
  457:             CALL DORGTR( UPLO, N, Z, LDZ, WORK( INDTAU ),
  458:      $                   WORK( INDWRK ), LLWORK, IINFO )
  459:             CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
  460:             CALL DSTEQR( JOBZ, N, W, WORK( INDEE ), Z, LDZ,
  461:      $                   WORK( INDWRK ), INFO )
  462:             IF( INFO.EQ.0 ) THEN
  463:                DO 30 I = 1, N
  464:                   IFAIL( I ) = 0
  465:    30          CONTINUE
  466:             END IF
  467:          END IF
  468:          IF( INFO.EQ.0 ) THEN
  469:             M = N
  470:             GO TO 40
  471:          END IF
  472:          INFO = 0
  473:       END IF
  474: *
  475: *     Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
  476: *
  477:       IF( WANTZ ) THEN
  478:          ORDER = 'B'
  479:       ELSE
  480:          ORDER = 'E'
  481:       END IF
  482:       INDIBL = 1
  483:       INDISP = INDIBL + N
  484:       INDIWO = INDISP + N
  485:       CALL DSTEBZ( RANGE, ORDER, N, VLL, VUU, IL, IU, ABSTLL,
  486:      $             WORK( INDD ), WORK( INDE ), M, NSPLIT, W,
  487:      $             IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWRK ),
  488:      $             IWORK( INDIWO ), INFO )
  489: *
  490:       IF( WANTZ ) THEN
  491:          CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W,
  492:      $                IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
  493:      $                WORK( INDWRK ), IWORK( INDIWO ), IFAIL, INFO )
  494: *
  495: *        Apply orthogonal matrix used in reduction to tridiagonal
  496: *        form to eigenvectors returned by DSTEIN.
  497: *
  498:          INDWKN = INDE
  499:          LLWRKN = LWORK - INDWKN + 1
  500:          CALL DORMTR( 'L', UPLO, 'N', N, M, A, LDA, WORK( INDTAU ), Z,
  501:      $                LDZ, WORK( INDWKN ), LLWRKN, IINFO )
  502:       END IF
  503: *
  504: *     If matrix was scaled, then rescale eigenvalues appropriately.
  505: *
  506:    40 CONTINUE
  507:       IF( ISCALE.EQ.1 ) THEN
  508:          IF( INFO.EQ.0 ) THEN
  509:             IMAX = M
  510:          ELSE
  511:             IMAX = INFO - 1
  512:          END IF
  513:          CALL DSCAL( IMAX, ONE / SIGMA, W, 1 )
  514:       END IF
  515: *
  516: *     If eigenvalues are not in order, then sort them, along with
  517: *     eigenvectors.
  518: *
  519:       IF( WANTZ ) THEN
  520:          DO 60 J = 1, M - 1
  521:             I = 0
  522:             TMP1 = W( J )
  523:             DO 50 JJ = J + 1, M
  524:                IF( W( JJ ).LT.TMP1 ) THEN
  525:                   I = JJ
  526:                   TMP1 = W( JJ )
  527:                END IF
  528:    50       CONTINUE
  529: *
  530:             IF( I.NE.0 ) THEN
  531:                ITMP1 = IWORK( INDIBL+I-1 )
  532:                W( I ) = W( J )
  533:                IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
  534:                W( J ) = TMP1
  535:                IWORK( INDIBL+J-1 ) = ITMP1
  536:                CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
  537:                IF( INFO.NE.0 ) THEN
  538:                   ITMP1 = IFAIL( I )
  539:                   IFAIL( I ) = IFAIL( J )
  540:                   IFAIL( J ) = ITMP1
  541:                END IF
  542:             END IF
  543:    60    CONTINUE
  544:       END IF
  545: *
  546: *     Set WORK(1) to optimal workspace size.
  547: *
  548:       WORK( 1 ) = LWKOPT
  549: *
  550:       RETURN
  551: *
  552: *     End of DSYEVX
  553: *
  554:       END

CVSweb interface <joel.bertrand@systella.fr>