File:  [local] / rpl / lapack / lapack / dsbgvx.f
Revision 1.20: download - view: text, annotated - select for diffs - revision graph
Mon Aug 7 08:39:06 2023 UTC (9 months ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_35, rpl-4_1_34, HEAD
Première mise à jour de lapack et blas.

    1: *> \brief \b DSBGVX
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DSBGVX + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsbgvx.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsbgvx.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsbgvx.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
   22: *                          LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
   23: *                          LDZ, WORK, IWORK, IFAIL, INFO )
   24: *
   25: *       .. Scalar Arguments ..
   26: *       CHARACTER          JOBZ, RANGE, UPLO
   27: *       INTEGER            IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
   28: *      $                   N
   29: *       DOUBLE PRECISION   ABSTOL, VL, VU
   30: *       ..
   31: *       .. Array Arguments ..
   32: *       INTEGER            IFAIL( * ), IWORK( * )
   33: *       DOUBLE PRECISION   AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
   34: *      $                   W( * ), WORK( * ), Z( LDZ, * )
   35: *       ..
   36: *
   37: *
   38: *> \par Purpose:
   39: *  =============
   40: *>
   41: *> \verbatim
   42: *>
   43: *> DSBGVX computes selected eigenvalues, and optionally, eigenvectors
   44: *> of a real generalized symmetric-definite banded eigenproblem, of
   45: *> the form A*x=(lambda)*B*x.  Here A and B are assumed to be symmetric
   46: *> and banded, and B is also positive definite.  Eigenvalues and
   47: *> eigenvectors can be selected by specifying either all eigenvalues,
   48: *> a range of values or a range of indices for the desired eigenvalues.
   49: *> \endverbatim
   50: *
   51: *  Arguments:
   52: *  ==========
   53: *
   54: *> \param[in] JOBZ
   55: *> \verbatim
   56: *>          JOBZ is CHARACTER*1
   57: *>          = 'N':  Compute eigenvalues only;
   58: *>          = 'V':  Compute eigenvalues and eigenvectors.
   59: *> \endverbatim
   60: *>
   61: *> \param[in] RANGE
   62: *> \verbatim
   63: *>          RANGE is CHARACTER*1
   64: *>          = 'A': all eigenvalues will be found.
   65: *>          = 'V': all eigenvalues in the half-open interval (VL,VU]
   66: *>                 will be found.
   67: *>          = 'I': the IL-th through IU-th eigenvalues will be found.
   68: *> \endverbatim
   69: *>
   70: *> \param[in] UPLO
   71: *> \verbatim
   72: *>          UPLO is CHARACTER*1
   73: *>          = 'U':  Upper triangles of A and B are stored;
   74: *>          = 'L':  Lower triangles of A and B are stored.
   75: *> \endverbatim
   76: *>
   77: *> \param[in] N
   78: *> \verbatim
   79: *>          N is INTEGER
   80: *>          The order of the matrices A and B.  N >= 0.
   81: *> \endverbatim
   82: *>
   83: *> \param[in] KA
   84: *> \verbatim
   85: *>          KA is INTEGER
   86: *>          The number of superdiagonals of the matrix A if UPLO = 'U',
   87: *>          or the number of subdiagonals if UPLO = 'L'.  KA >= 0.
   88: *> \endverbatim
   89: *>
   90: *> \param[in] KB
   91: *> \verbatim
   92: *>          KB is INTEGER
   93: *>          The number of superdiagonals of the matrix B if UPLO = 'U',
   94: *>          or the number of subdiagonals if UPLO = 'L'.  KB >= 0.
   95: *> \endverbatim
   96: *>
   97: *> \param[in,out] AB
   98: *> \verbatim
   99: *>          AB is DOUBLE PRECISION array, dimension (LDAB, N)
  100: *>          On entry, the upper or lower triangle of the symmetric band
  101: *>          matrix A, stored in the first ka+1 rows of the array.  The
  102: *>          j-th column of A is stored in the j-th column of the array AB
  103: *>          as follows:
  104: *>          if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
  105: *>          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+ka).
  106: *>
  107: *>          On exit, the contents of AB are destroyed.
  108: *> \endverbatim
  109: *>
  110: *> \param[in] LDAB
  111: *> \verbatim
  112: *>          LDAB is INTEGER
  113: *>          The leading dimension of the array AB.  LDAB >= KA+1.
  114: *> \endverbatim
  115: *>
  116: *> \param[in,out] BB
  117: *> \verbatim
  118: *>          BB is DOUBLE PRECISION array, dimension (LDBB, N)
  119: *>          On entry, the upper or lower triangle of the symmetric band
  120: *>          matrix B, stored in the first kb+1 rows of the array.  The
  121: *>          j-th column of B is stored in the j-th column of the array BB
  122: *>          as follows:
  123: *>          if UPLO = 'U', BB(ka+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
  124: *>          if UPLO = 'L', BB(1+i-j,j)    = B(i,j) for j<=i<=min(n,j+kb).
  125: *>
  126: *>          On exit, the factor S from the split Cholesky factorization
  127: *>          B = S**T*S, as returned by DPBSTF.
  128: *> \endverbatim
  129: *>
  130: *> \param[in] LDBB
  131: *> \verbatim
  132: *>          LDBB is INTEGER
  133: *>          The leading dimension of the array BB.  LDBB >= KB+1.
  134: *> \endverbatim
  135: *>
  136: *> \param[out] Q
  137: *> \verbatim
  138: *>          Q is DOUBLE PRECISION array, dimension (LDQ, N)
  139: *>          If JOBZ = 'V', the n-by-n matrix used in the reduction of
  140: *>          A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
  141: *>          and consequently C to tridiagonal form.
  142: *>          If JOBZ = 'N', the array Q is not referenced.
  143: *> \endverbatim
  144: *>
  145: *> \param[in] LDQ
  146: *> \verbatim
  147: *>          LDQ is INTEGER
  148: *>          The leading dimension of the array Q.  If JOBZ = 'N',
  149: *>          LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).
  150: *> \endverbatim
  151: *>
  152: *> \param[in] VL
  153: *> \verbatim
  154: *>          VL is DOUBLE PRECISION
  155: *>
  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: *>
  165: *>          If RANGE='V', the upper bound of the interval to
  166: *>          be searched for eigenvalues. VL < VU.
  167: *>          Not referenced if RANGE = 'A' or 'I'.
  168: *> \endverbatim
  169: *>
  170: *> \param[in] IL
  171: *> \verbatim
  172: *>          IL is INTEGER
  173: *>
  174: *>          If RANGE='I', the index of the
  175: *>          smallest eigenvalue to be returned.
  176: *>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
  177: *>          Not referenced if RANGE = 'A' or 'V'.
  178: *> \endverbatim
  179: *>
  180: *> \param[in] IU
  181: *> \verbatim
  182: *>          IU is INTEGER
  183: *>
  184: *>          If RANGE='I', the index of the
  185: *>          largest eigenvalue to be returned.
  186: *>          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
  187: *>          Not referenced if RANGE = 'A' or 'V'.
  188: *> \endverbatim
  189: *>
  190: *> \param[in] ABSTOL
  191: *> \verbatim
  192: *>          ABSTOL is DOUBLE PRECISION
  193: *>          The absolute error tolerance for the eigenvalues.
  194: *>          An approximate eigenvalue is accepted as converged
  195: *>          when it is determined to lie in an interval [a,b]
  196: *>          of width less than or equal to
  197: *>
  198: *>                  ABSTOL + EPS *   max( |a|,|b| ) ,
  199: *>
  200: *>          where EPS is the machine precision.  If ABSTOL is less than
  201: *>          or equal to zero, then  EPS*|T|  will be used in its place,
  202: *>          where |T| is the 1-norm of the tridiagonal matrix obtained
  203: *>          by reducing A to tridiagonal form.
  204: *>
  205: *>          Eigenvalues will be computed most accurately when ABSTOL is
  206: *>          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
  207: *>          If this routine returns with INFO>0, indicating that some
  208: *>          eigenvectors did not converge, try setting ABSTOL to
  209: *>          2*DLAMCH('S').
  210: *> \endverbatim
  211: *>
  212: *> \param[out] M
  213: *> \verbatim
  214: *>          M is INTEGER
  215: *>          The total number of eigenvalues found.  0 <= M <= N.
  216: *>          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
  217: *> \endverbatim
  218: *>
  219: *> \param[out] W
  220: *> \verbatim
  221: *>          W is DOUBLE PRECISION array, dimension (N)
  222: *>          If INFO = 0, the eigenvalues in ascending order.
  223: *> \endverbatim
  224: *>
  225: *> \param[out] Z
  226: *> \verbatim
  227: *>          Z is DOUBLE PRECISION array, dimension (LDZ, N)
  228: *>          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
  229: *>          eigenvectors, with the i-th column of Z holding the
  230: *>          eigenvector associated with W(i).  The eigenvectors are
  231: *>          normalized so Z**T*B*Z = I.
  232: *>          If JOBZ = 'N', then Z is not referenced.
  233: *> \endverbatim
  234: *>
  235: *> \param[in] LDZ
  236: *> \verbatim
  237: *>          LDZ is INTEGER
  238: *>          The leading dimension of the array Z.  LDZ >= 1, and if
  239: *>          JOBZ = 'V', LDZ >= max(1,N).
  240: *> \endverbatim
  241: *>
  242: *> \param[out] WORK
  243: *> \verbatim
  244: *>          WORK is DOUBLE PRECISION array, dimension (7*N)
  245: *> \endverbatim
  246: *>
  247: *> \param[out] IWORK
  248: *> \verbatim
  249: *>          IWORK is INTEGER array, dimension (5*N)
  250: *> \endverbatim
  251: *>
  252: *> \param[out] IFAIL
  253: *> \verbatim
  254: *>          IFAIL is INTEGER array, dimension (M)
  255: *>          If JOBZ = 'V', then if INFO = 0, the first M elements of
  256: *>          IFAIL are zero.  If INFO > 0, then IFAIL contains the
  257: *>          indices of the eigenvalues that failed to converge.
  258: *>          If JOBZ = 'N', then IFAIL is not referenced.
  259: *> \endverbatim
  260: *>
  261: *> \param[out] INFO
  262: *> \verbatim
  263: *>          INFO is INTEGER
  264: *>          = 0:  successful exit
  265: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
  266: *>          <= N: if INFO = i, then i eigenvectors failed to converge.
  267: *>                  Their indices are stored in IFAIL.
  268: *>          > N:  DPBSTF returned an error code; i.e.,
  269: *>                if INFO = N + i, for 1 <= i <= N, then the leading
  270: *>                minor of order i of B is not positive definite.
  271: *>                The factorization of B could not be completed and
  272: *>                no eigenvalues or eigenvectors were computed.
  273: *> \endverbatim
  274: *
  275: *  Authors:
  276: *  ========
  277: *
  278: *> \author Univ. of Tennessee
  279: *> \author Univ. of California Berkeley
  280: *> \author Univ. of Colorado Denver
  281: *> \author NAG Ltd.
  282: *
  283: *> \ingroup doubleOTHEReigen
  284: *
  285: *> \par Contributors:
  286: *  ==================
  287: *>
  288: *>     Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
  289: *
  290: *  =====================================================================
  291:       SUBROUTINE DSBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
  292:      $                   LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
  293:      $                   LDZ, WORK, IWORK, IFAIL, INFO )
  294: *
  295: *  -- LAPACK driver routine --
  296: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  297: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  298: *
  299: *     .. Scalar Arguments ..
  300:       CHARACTER          JOBZ, RANGE, UPLO
  301:       INTEGER            IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
  302:      $                   N
  303:       DOUBLE PRECISION   ABSTOL, VL, VU
  304: *     ..
  305: *     .. Array Arguments ..
  306:       INTEGER            IFAIL( * ), IWORK( * )
  307:       DOUBLE PRECISION   AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
  308:      $                   W( * ), WORK( * ), Z( LDZ, * )
  309: *     ..
  310: *
  311: *  =====================================================================
  312: *
  313: *     .. Parameters ..
  314:       DOUBLE PRECISION   ZERO, ONE
  315:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
  316: *     ..
  317: *     .. Local Scalars ..
  318:       LOGICAL            ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
  319:       CHARACTER          ORDER, VECT
  320:       INTEGER            I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP,
  321:      $                   INDIWO, INDWRK, ITMP1, J, JJ, NSPLIT
  322:       DOUBLE PRECISION   TMP1
  323: *     ..
  324: *     .. External Functions ..
  325:       LOGICAL            LSAME
  326:       EXTERNAL           LSAME
  327: *     ..
  328: *     .. External Subroutines ..
  329:       EXTERNAL           DCOPY, DGEMV, DLACPY, DPBSTF, DSBGST, DSBTRD,
  330:      $                   DSTEBZ, DSTEIN, DSTEQR, DSTERF, DSWAP, XERBLA
  331: *     ..
  332: *     .. Intrinsic Functions ..
  333:       INTRINSIC          MIN
  334: *     ..
  335: *     .. Executable Statements ..
  336: *
  337: *     Test the input parameters.
  338: *
  339:       WANTZ = LSAME( JOBZ, 'V' )
  340:       UPPER = LSAME( UPLO, 'U' )
  341:       ALLEIG = LSAME( RANGE, 'A' )
  342:       VALEIG = LSAME( RANGE, 'V' )
  343:       INDEIG = LSAME( RANGE, 'I' )
  344: *
  345:       INFO = 0
  346:       IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
  347:          INFO = -1
  348:       ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN
  349:          INFO = -2
  350:       ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
  351:          INFO = -3
  352:       ELSE IF( N.LT.0 ) THEN
  353:          INFO = -4
  354:       ELSE IF( KA.LT.0 ) THEN
  355:          INFO = -5
  356:       ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN
  357:          INFO = -6
  358:       ELSE IF( LDAB.LT.KA+1 ) THEN
  359:          INFO = -8
  360:       ELSE IF( LDBB.LT.KB+1 ) THEN
  361:          INFO = -10
  362:       ELSE IF( LDQ.LT.1 .OR. ( WANTZ .AND. LDQ.LT.N ) ) THEN
  363:          INFO = -12
  364:       ELSE
  365:          IF( VALEIG ) THEN
  366:             IF( N.GT.0 .AND. VU.LE.VL )
  367:      $         INFO = -14
  368:          ELSE IF( INDEIG ) THEN
  369:             IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN
  370:                INFO = -15
  371:             ELSE IF ( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN
  372:                INFO = -16
  373:             END IF
  374:          END IF
  375:       END IF
  376:       IF( INFO.EQ.0) THEN
  377:          IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
  378:             INFO = -21
  379:          END IF
  380:       END IF
  381: *
  382:       IF( INFO.NE.0 ) THEN
  383:          CALL XERBLA( 'DSBGVX', -INFO )
  384:          RETURN
  385:       END IF
  386: *
  387: *     Quick return if possible
  388: *
  389:       M = 0
  390:       IF( N.EQ.0 )
  391:      $   RETURN
  392: *
  393: *     Form a split Cholesky factorization of B.
  394: *
  395:       CALL DPBSTF( UPLO, N, KB, BB, LDBB, INFO )
  396:       IF( INFO.NE.0 ) THEN
  397:          INFO = N + INFO
  398:          RETURN
  399:       END IF
  400: *
  401: *     Transform problem to standard eigenvalue problem.
  402: *
  403:       CALL DSBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ,
  404:      $             WORK, IINFO )
  405: *
  406: *     Reduce symmetric band matrix to tridiagonal form.
  407: *
  408:       INDD = 1
  409:       INDE = INDD + N
  410:       INDWRK = INDE + N
  411:       IF( WANTZ ) THEN
  412:          VECT = 'U'
  413:       ELSE
  414:          VECT = 'N'
  415:       END IF
  416:       CALL DSBTRD( VECT, UPLO, N, KA, AB, LDAB, WORK( INDD ),
  417:      $             WORK( INDE ), Q, LDQ, WORK( INDWRK ), IINFO )
  418: *
  419: *     If all eigenvalues are desired and ABSTOL is less than or equal
  420: *     to zero, then call DSTERF or SSTEQR.  If this fails for some
  421: *     eigenvalue, then try DSTEBZ.
  422: *
  423:       TEST = .FALSE.
  424:       IF( INDEIG ) THEN
  425:          IF( IL.EQ.1 .AND. IU.EQ.N ) THEN
  426:             TEST = .TRUE.
  427:          END IF
  428:       END IF
  429:       IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN
  430:          CALL DCOPY( N, WORK( INDD ), 1, W, 1 )
  431:          INDEE = INDWRK + 2*N
  432:          CALL DCOPY( N-1, WORK( INDE ), 1, WORK( INDEE ), 1 )
  433:          IF( .NOT.WANTZ ) THEN
  434:             CALL DSTERF( N, W, WORK( INDEE ), INFO )
  435:          ELSE
  436:             CALL DLACPY( 'A', N, N, Q, LDQ, Z, LDZ )
  437:             CALL DSTEQR( JOBZ, N, W, WORK( INDEE ), Z, LDZ,
  438:      $                   WORK( INDWRK ), INFO )
  439:             IF( INFO.EQ.0 ) THEN
  440:                DO 10 I = 1, N
  441:                   IFAIL( I ) = 0
  442:    10          CONTINUE
  443:             END IF
  444:          END IF
  445:          IF( INFO.EQ.0 ) THEN
  446:             M = N
  447:             GO TO 30
  448:          END IF
  449:          INFO = 0
  450:       END IF
  451: *
  452: *     Otherwise, call DSTEBZ and, if eigenvectors are desired,
  453: *     call DSTEIN.
  454: *
  455:       IF( WANTZ ) THEN
  456:          ORDER = 'B'
  457:       ELSE
  458:          ORDER = 'E'
  459:       END IF
  460:       INDIBL = 1
  461:       INDISP = INDIBL + N
  462:       INDIWO = INDISP + N
  463:       CALL DSTEBZ( RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL,
  464:      $             WORK( INDD ), WORK( INDE ), M, NSPLIT, W,
  465:      $             IWORK( INDIBL ), IWORK( INDISP ), WORK( INDWRK ),
  466:      $             IWORK( INDIWO ), INFO )
  467: *
  468:       IF( WANTZ ) THEN
  469:          CALL DSTEIN( N, WORK( INDD ), WORK( INDE ), M, W,
  470:      $                IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ,
  471:      $                WORK( INDWRK ), IWORK( INDIWO ), IFAIL, INFO )
  472: *
  473: *        Apply transformation matrix used in reduction to tridiagonal
  474: *        form to eigenvectors returned by DSTEIN.
  475: *
  476:          DO 20 J = 1, M
  477:             CALL DCOPY( N, Z( 1, J ), 1, WORK( 1 ), 1 )
  478:             CALL DGEMV( 'N', N, N, ONE, Q, LDQ, WORK, 1, ZERO,
  479:      $                  Z( 1, J ), 1 )
  480:    20    CONTINUE
  481:       END IF
  482: *
  483:    30 CONTINUE
  484: *
  485: *     If eigenvalues are not in order, then sort them, along with
  486: *     eigenvectors.
  487: *
  488:       IF( WANTZ ) THEN
  489:          DO 50 J = 1, M - 1
  490:             I = 0
  491:             TMP1 = W( J )
  492:             DO 40 JJ = J + 1, M
  493:                IF( W( JJ ).LT.TMP1 ) THEN
  494:                   I = JJ
  495:                   TMP1 = W( JJ )
  496:                END IF
  497:    40       CONTINUE
  498: *
  499:             IF( I.NE.0 ) THEN
  500:                ITMP1 = IWORK( INDIBL+I-1 )
  501:                W( I ) = W( J )
  502:                IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 )
  503:                W( J ) = TMP1
  504:                IWORK( INDIBL+J-1 ) = ITMP1
  505:                CALL DSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 )
  506:                IF( INFO.NE.0 ) THEN
  507:                   ITMP1 = IFAIL( I )
  508:                   IFAIL( I ) = IFAIL( J )
  509:                   IFAIL( J ) = ITMP1
  510:                END IF
  511:             END IF
  512:    50    CONTINUE
  513:       END IF
  514: *
  515:       RETURN
  516: *
  517: *     End of DSBGVX
  518: *
  519:       END

CVSweb interface <joel.bertrand@systella.fr>