File:  [local] / rpl / lapack / lapack / dlaqr0.f
Revision 1.1: download - view: text, annotated - select for diffs - revision graph
Tue Jan 26 15:22:45 2010 UTC (14 years, 4 months ago) by bertrand
Branches: MAIN
CVS tags: HEAD
Initial revision

    1:       SUBROUTINE DLAQR0( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI,
    2:      $                   ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO )
    3: *
    4: *  -- LAPACK auxiliary routine (version 3.2) --
    5: *     Univ. of Tennessee, Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..
    6: *     November 2006
    7: *
    8: *     .. Scalar Arguments ..
    9:       INTEGER            IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N
   10:       LOGICAL            WANTT, WANTZ
   11: *     ..
   12: *     .. Array Arguments ..
   13:       DOUBLE PRECISION   H( LDH, * ), WI( * ), WORK( * ), WR( * ),
   14:      $                   Z( LDZ, * )
   15: *     ..
   16: *
   17: *     Purpose
   18: *     =======
   19: *
   20: *     DLAQR0 computes the eigenvalues of a Hessenberg matrix H
   21: *     and, optionally, the matrices T and Z from the Schur decomposition
   22: *     H = Z T Z**T, where T is an upper quasi-triangular matrix (the
   23: *     Schur form), and Z is the orthogonal matrix of Schur vectors.
   24: *
   25: *     Optionally Z may be postmultiplied into an input orthogonal
   26: *     matrix Q so that this routine can give the Schur factorization
   27: *     of a matrix A which has been reduced to the Hessenberg form H
   28: *     by the orthogonal matrix Q:  A = Q*H*Q**T = (QZ)*T*(QZ)**T.
   29: *
   30: *     Arguments
   31: *     =========
   32: *
   33: *     WANTT   (input) LOGICAL
   34: *          = .TRUE. : the full Schur form T is required;
   35: *          = .FALSE.: only eigenvalues are required.
   36: *
   37: *     WANTZ   (input) LOGICAL
   38: *          = .TRUE. : the matrix of Schur vectors Z is required;
   39: *          = .FALSE.: Schur vectors are not required.
   40: *
   41: *     N     (input) INTEGER
   42: *           The order of the matrix H.  N .GE. 0.
   43: *
   44: *     ILO   (input) INTEGER
   45: *     IHI   (input) INTEGER
   46: *           It is assumed that H is already upper triangular in rows
   47: *           and columns 1:ILO-1 and IHI+1:N and, if ILO.GT.1,
   48: *           H(ILO,ILO-1) is zero. ILO and IHI are normally set by a
   49: *           previous call to DGEBAL, and then passed to DGEHRD when the
   50: *           matrix output by DGEBAL is reduced to Hessenberg form.
   51: *           Otherwise, ILO and IHI should be set to 1 and N,
   52: *           respectively.  If N.GT.0, then 1.LE.ILO.LE.IHI.LE.N.
   53: *           If N = 0, then ILO = 1 and IHI = 0.
   54: *
   55: *     H     (input/output) DOUBLE PRECISION array, dimension (LDH,N)
   56: *           On entry, the upper Hessenberg matrix H.
   57: *           On exit, if INFO = 0 and WANTT is .TRUE., then H contains
   58: *           the upper quasi-triangular matrix T from the Schur
   59: *           decomposition (the Schur form); 2-by-2 diagonal blocks
   60: *           (corresponding to complex conjugate pairs of eigenvalues)
   61: *           are returned in standard form, with H(i,i) = H(i+1,i+1)
   62: *           and H(i+1,i)*H(i,i+1).LT.0. If INFO = 0 and WANTT is
   63: *           .FALSE., then the contents of H are unspecified on exit.
   64: *           (The output value of H when INFO.GT.0 is given under the
   65: *           description of INFO below.)
   66: *
   67: *           This subroutine may explicitly set H(i,j) = 0 for i.GT.j and
   68: *           j = 1, 2, ... ILO-1 or j = IHI+1, IHI+2, ... N.
   69: *
   70: *     LDH   (input) INTEGER
   71: *           The leading dimension of the array H. LDH .GE. max(1,N).
   72: *
   73: *     WR    (output) DOUBLE PRECISION array, dimension (IHI)
   74: *     WI    (output) DOUBLE PRECISION array, dimension (IHI)
   75: *           The real and imaginary parts, respectively, of the computed
   76: *           eigenvalues of H(ILO:IHI,ILO:IHI) are stored in WR(ILO:IHI)
   77: *           and WI(ILO:IHI). If two eigenvalues are computed as a
   78: *           complex conjugate pair, they are stored in consecutive
   79: *           elements of WR and WI, say the i-th and (i+1)th, with
   80: *           WI(i) .GT. 0 and WI(i+1) .LT. 0. If WANTT is .TRUE., then
   81: *           the eigenvalues are stored in the same order as on the
   82: *           diagonal of the Schur form returned in H, with
   83: *           WR(i) = H(i,i) and, if H(i:i+1,i:i+1) is a 2-by-2 diagonal
   84: *           block, WI(i) = sqrt(-H(i+1,i)*H(i,i+1)) and
   85: *           WI(i+1) = -WI(i).
   86: *
   87: *     ILOZ     (input) INTEGER
   88: *     IHIZ     (input) INTEGER
   89: *           Specify the rows of Z to which transformations must be
   90: *           applied if WANTZ is .TRUE..
   91: *           1 .LE. ILOZ .LE. ILO; IHI .LE. IHIZ .LE. N.
   92: *
   93: *     Z     (input/output) DOUBLE PRECISION array, dimension (LDZ,IHI)
   94: *           If WANTZ is .FALSE., then Z is not referenced.
   95: *           If WANTZ is .TRUE., then Z(ILO:IHI,ILOZ:IHIZ) is
   96: *           replaced by Z(ILO:IHI,ILOZ:IHIZ)*U where U is the
   97: *           orthogonal Schur factor of H(ILO:IHI,ILO:IHI).
   98: *           (The output value of Z when INFO.GT.0 is given under
   99: *           the description of INFO below.)
  100: *
  101: *     LDZ   (input) INTEGER
  102: *           The leading dimension of the array Z.  if WANTZ is .TRUE.
  103: *           then LDZ.GE.MAX(1,IHIZ).  Otherwize, LDZ.GE.1.
  104: *
  105: *     WORK  (workspace/output) DOUBLE PRECISION array, dimension LWORK
  106: *           On exit, if LWORK = -1, WORK(1) returns an estimate of
  107: *           the optimal value for LWORK.
  108: *
  109: *     LWORK (input) INTEGER
  110: *           The dimension of the array WORK.  LWORK .GE. max(1,N)
  111: *           is sufficient, but LWORK typically as large as 6*N may
  112: *           be required for optimal performance.  A workspace query
  113: *           to determine the optimal workspace size is recommended.
  114: *
  115: *           If LWORK = -1, then DLAQR0 does a workspace query.
  116: *           In this case, DLAQR0 checks the input parameters and
  117: *           estimates the optimal workspace size for the given
  118: *           values of N, ILO and IHI.  The estimate is returned
  119: *           in WORK(1).  No error message related to LWORK is
  120: *           issued by XERBLA.  Neither H nor Z are accessed.
  121: *
  122: *
  123: *     INFO  (output) INTEGER
  124: *             =  0:  successful exit
  125: *           .GT. 0:  if INFO = i, DLAQR0 failed to compute all of
  126: *                the eigenvalues.  Elements 1:ilo-1 and i+1:n of WR
  127: *                and WI contain those eigenvalues which have been
  128: *                successfully computed.  (Failures are rare.)
  129: *
  130: *                If INFO .GT. 0 and WANT is .FALSE., then on exit,
  131: *                the remaining unconverged eigenvalues are the eigen-
  132: *                values of the upper Hessenberg matrix rows and
  133: *                columns ILO through INFO of the final, output
  134: *                value of H.
  135: *
  136: *                If INFO .GT. 0 and WANTT is .TRUE., then on exit
  137: *
  138: *           (*)  (initial value of H)*U  = U*(final value of H)
  139: *
  140: *                where U is an orthogonal matrix.  The final
  141: *                value of H is upper Hessenberg and quasi-triangular
  142: *                in rows and columns INFO+1 through IHI.
  143: *
  144: *                If INFO .GT. 0 and WANTZ is .TRUE., then on exit
  145: *
  146: *                  (final value of Z(ILO:IHI,ILOZ:IHIZ)
  147: *                   =  (initial value of Z(ILO:IHI,ILOZ:IHIZ)*U
  148: *
  149: *                where U is the orthogonal matrix in (*) (regard-
  150: *                less of the value of WANTT.)
  151: *
  152: *                If INFO .GT. 0 and WANTZ is .FALSE., then Z is not
  153: *                accessed.
  154: *
  155: *     ================================================================
  156: *     Based on contributions by
  157: *        Karen Braman and Ralph Byers, Department of Mathematics,
  158: *        University of Kansas, USA
  159: *
  160: *     ================================================================
  161: *     References:
  162: *       K. Braman, R. Byers and R. Mathias, The Multi-Shift QR
  163: *       Algorithm Part I: Maintaining Well Focused Shifts, and Level 3
  164: *       Performance, SIAM Journal of Matrix Analysis, volume 23, pages
  165: *       929--947, 2002.
  166: *
  167: *       K. Braman, R. Byers and R. Mathias, The Multi-Shift QR
  168: *       Algorithm Part II: Aggressive Early Deflation, SIAM Journal
  169: *       of Matrix Analysis, volume 23, pages 948--973, 2002.
  170: *
  171: *     ================================================================
  172: *     .. Parameters ..
  173: *
  174: *     ==== Matrices of order NTINY or smaller must be processed by
  175: *     .    DLAHQR because of insufficient subdiagonal scratch space.
  176: *     .    (This is a hard limit.) ====
  177:       INTEGER            NTINY
  178:       PARAMETER          ( NTINY = 11 )
  179: *
  180: *     ==== Exceptional deflation windows:  try to cure rare
  181: *     .    slow convergence by varying the size of the
  182: *     .    deflation window after KEXNW iterations. ====
  183:       INTEGER            KEXNW
  184:       PARAMETER          ( KEXNW = 5 )
  185: *
  186: *     ==== Exceptional shifts: try to cure rare slow convergence
  187: *     .    with ad-hoc exceptional shifts every KEXSH iterations.
  188: *     .    ====
  189:       INTEGER            KEXSH
  190:       PARAMETER          ( KEXSH = 6 )
  191: *
  192: *     ==== The constants WILK1 and WILK2 are used to form the
  193: *     .    exceptional shifts. ====
  194:       DOUBLE PRECISION   WILK1, WILK2
  195:       PARAMETER          ( WILK1 = 0.75d0, WILK2 = -0.4375d0 )
  196:       DOUBLE PRECISION   ZERO, ONE
  197:       PARAMETER          ( ZERO = 0.0d0, ONE = 1.0d0 )
  198: *     ..
  199: *     .. Local Scalars ..
  200:       DOUBLE PRECISION   AA, BB, CC, CS, DD, SN, SS, SWAP
  201:       INTEGER            I, INF, IT, ITMAX, K, KACC22, KBOT, KDU, KS,
  202:      $                   KT, KTOP, KU, KV, KWH, KWTOP, KWV, LD, LS,
  203:      $                   LWKOPT, NDEC, NDFL, NH, NHO, NIBBLE, NMIN, NS,
  204:      $                   NSMAX, NSR, NVE, NW, NWMAX, NWR, NWUPBD
  205:       LOGICAL            SORTED
  206:       CHARACTER          JBCMPZ*2
  207: *     ..
  208: *     .. External Functions ..
  209:       INTEGER            ILAENV
  210:       EXTERNAL           ILAENV
  211: *     ..
  212: *     .. Local Arrays ..
  213:       DOUBLE PRECISION   ZDUM( 1, 1 )
  214: *     ..
  215: *     .. External Subroutines ..
  216:       EXTERNAL           DLACPY, DLAHQR, DLANV2, DLAQR3, DLAQR4, DLAQR5
  217: *     ..
  218: *     .. Intrinsic Functions ..
  219:       INTRINSIC          ABS, DBLE, INT, MAX, MIN, MOD
  220: *     ..
  221: *     .. Executable Statements ..
  222:       INFO = 0
  223: *
  224: *     ==== Quick return for N = 0: nothing to do. ====
  225: *
  226:       IF( N.EQ.0 ) THEN
  227:          WORK( 1 ) = ONE
  228:          RETURN
  229:       END IF
  230: *
  231:       IF( N.LE.NTINY ) THEN
  232: *
  233: *        ==== Tiny matrices must use DLAHQR. ====
  234: *
  235:          LWKOPT = 1
  236:          IF( LWORK.NE.-1 )
  237:      $      CALL DLAHQR( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI,
  238:      $                   ILOZ, IHIZ, Z, LDZ, INFO )
  239:       ELSE
  240: *
  241: *        ==== Use small bulge multi-shift QR with aggressive early
  242: *        .    deflation on larger-than-tiny matrices. ====
  243: *
  244: *        ==== Hope for the best. ====
  245: *
  246:          INFO = 0
  247: *
  248: *        ==== Set up job flags for ILAENV. ====
  249: *
  250:          IF( WANTT ) THEN
  251:             JBCMPZ( 1: 1 ) = 'S'
  252:          ELSE
  253:             JBCMPZ( 1: 1 ) = 'E'
  254:          END IF
  255:          IF( WANTZ ) THEN
  256:             JBCMPZ( 2: 2 ) = 'V'
  257:          ELSE
  258:             JBCMPZ( 2: 2 ) = 'N'
  259:          END IF
  260: *
  261: *        ==== NWR = recommended deflation window size.  At this
  262: *        .    point,  N .GT. NTINY = 11, so there is enough
  263: *        .    subdiagonal workspace for NWR.GE.2 as required.
  264: *        .    (In fact, there is enough subdiagonal space for
  265: *        .    NWR.GE.3.) ====
  266: *
  267:          NWR = ILAENV( 13, 'DLAQR0', JBCMPZ, N, ILO, IHI, LWORK )
  268:          NWR = MAX( 2, NWR )
  269:          NWR = MIN( IHI-ILO+1, ( N-1 ) / 3, NWR )
  270: *
  271: *        ==== NSR = recommended number of simultaneous shifts.
  272: *        .    At this point N .GT. NTINY = 11, so there is at
  273: *        .    enough subdiagonal workspace for NSR to be even
  274: *        .    and greater than or equal to two as required. ====
  275: *
  276:          NSR = ILAENV( 15, 'DLAQR0', JBCMPZ, N, ILO, IHI, LWORK )
  277:          NSR = MIN( NSR, ( N+6 ) / 9, IHI-ILO )
  278:          NSR = MAX( 2, NSR-MOD( NSR, 2 ) )
  279: *
  280: *        ==== Estimate optimal workspace ====
  281: *
  282: *        ==== Workspace query call to DLAQR3 ====
  283: *
  284:          CALL DLAQR3( WANTT, WANTZ, N, ILO, IHI, NWR+1, H, LDH, ILOZ,
  285:      $                IHIZ, Z, LDZ, LS, LD, WR, WI, H, LDH, N, H, LDH,
  286:      $                N, H, LDH, WORK, -1 )
  287: *
  288: *        ==== Optimal workspace = MAX(DLAQR5, DLAQR3) ====
  289: *
  290:          LWKOPT = MAX( 3*NSR / 2, INT( WORK( 1 ) ) )
  291: *
  292: *        ==== Quick return in case of workspace query. ====
  293: *
  294:          IF( LWORK.EQ.-1 ) THEN
  295:             WORK( 1 ) = DBLE( LWKOPT )
  296:             RETURN
  297:          END IF
  298: *
  299: *        ==== DLAHQR/DLAQR0 crossover point ====
  300: *
  301:          NMIN = ILAENV( 12, 'DLAQR0', JBCMPZ, N, ILO, IHI, LWORK )
  302:          NMIN = MAX( NTINY, NMIN )
  303: *
  304: *        ==== Nibble crossover point ====
  305: *
  306:          NIBBLE = ILAENV( 14, 'DLAQR0', JBCMPZ, N, ILO, IHI, LWORK )
  307:          NIBBLE = MAX( 0, NIBBLE )
  308: *
  309: *        ==== Accumulate reflections during ttswp?  Use block
  310: *        .    2-by-2 structure during matrix-matrix multiply? ====
  311: *
  312:          KACC22 = ILAENV( 16, 'DLAQR0', JBCMPZ, N, ILO, IHI, LWORK )
  313:          KACC22 = MAX( 0, KACC22 )
  314:          KACC22 = MIN( 2, KACC22 )
  315: *
  316: *        ==== NWMAX = the largest possible deflation window for
  317: *        .    which there is sufficient workspace. ====
  318: *
  319:          NWMAX = MIN( ( N-1 ) / 3, LWORK / 2 )
  320:          NW = NWMAX
  321: *
  322: *        ==== NSMAX = the Largest number of simultaneous shifts
  323: *        .    for which there is sufficient workspace. ====
  324: *
  325:          NSMAX = MIN( ( N+6 ) / 9, 2*LWORK / 3 )
  326:          NSMAX = NSMAX - MOD( NSMAX, 2 )
  327: *
  328: *        ==== NDFL: an iteration count restarted at deflation. ====
  329: *
  330:          NDFL = 1
  331: *
  332: *        ==== ITMAX = iteration limit ====
  333: *
  334:          ITMAX = MAX( 30, 2*KEXSH )*MAX( 10, ( IHI-ILO+1 ) )
  335: *
  336: *        ==== Last row and column in the active block ====
  337: *
  338:          KBOT = IHI
  339: *
  340: *        ==== Main Loop ====
  341: *
  342:          DO 80 IT = 1, ITMAX
  343: *
  344: *           ==== Done when KBOT falls below ILO ====
  345: *
  346:             IF( KBOT.LT.ILO )
  347:      $         GO TO 90
  348: *
  349: *           ==== Locate active block ====
  350: *
  351:             DO 10 K = KBOT, ILO + 1, -1
  352:                IF( H( K, K-1 ).EQ.ZERO )
  353:      $            GO TO 20
  354:    10       CONTINUE
  355:             K = ILO
  356:    20       CONTINUE
  357:             KTOP = K
  358: *
  359: *           ==== Select deflation window size:
  360: *           .    Typical Case:
  361: *           .      If possible and advisable, nibble the entire
  362: *           .      active block.  If not, use size MIN(NWR,NWMAX)
  363: *           .      or MIN(NWR+1,NWMAX) depending upon which has
  364: *           .      the smaller corresponding subdiagonal entry
  365: *           .      (a heuristic).
  366: *           .
  367: *           .    Exceptional Case:
  368: *           .      If there have been no deflations in KEXNW or
  369: *           .      more iterations, then vary the deflation window
  370: *           .      size.   At first, because, larger windows are,
  371: *           .      in general, more powerful than smaller ones,
  372: *           .      rapidly increase the window to the maximum possible.
  373: *           .      Then, gradually reduce the window size. ====
  374: *
  375:             NH = KBOT - KTOP + 1
  376:             NWUPBD = MIN( NH, NWMAX )
  377:             IF( NDFL.LT.KEXNW ) THEN
  378:                NW = MIN( NWUPBD, NWR )
  379:             ELSE
  380:                NW = MIN( NWUPBD, 2*NW )
  381:             END IF
  382:             IF( NW.LT.NWMAX ) THEN
  383:                IF( NW.GE.NH-1 ) THEN
  384:                   NW = NH
  385:                ELSE
  386:                   KWTOP = KBOT - NW + 1
  387:                   IF( ABS( H( KWTOP, KWTOP-1 ) ).GT.
  388:      $                ABS( H( KWTOP-1, KWTOP-2 ) ) )NW = NW + 1
  389:                END IF
  390:             END IF
  391:             IF( NDFL.LT.KEXNW ) THEN
  392:                NDEC = -1
  393:             ELSE IF( NDEC.GE.0 .OR. NW.GE.NWUPBD ) THEN
  394:                NDEC = NDEC + 1
  395:                IF( NW-NDEC.LT.2 )
  396:      $            NDEC = 0
  397:                NW = NW - NDEC
  398:             END IF
  399: *
  400: *           ==== Aggressive early deflation:
  401: *           .    split workspace under the subdiagonal into
  402: *           .      - an nw-by-nw work array V in the lower
  403: *           .        left-hand-corner,
  404: *           .      - an NW-by-at-least-NW-but-more-is-better
  405: *           .        (NW-by-NHO) horizontal work array along
  406: *           .        the bottom edge,
  407: *           .      - an at-least-NW-but-more-is-better (NHV-by-NW)
  408: *           .        vertical work array along the left-hand-edge.
  409: *           .        ====
  410: *
  411:             KV = N - NW + 1
  412:             KT = NW + 1
  413:             NHO = ( N-NW-1 ) - KT + 1
  414:             KWV = NW + 2
  415:             NVE = ( N-NW ) - KWV + 1
  416: *
  417: *           ==== Aggressive early deflation ====
  418: *
  419:             CALL DLAQR3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
  420:      $                   IHIZ, Z, LDZ, LS, LD, WR, WI, H( KV, 1 ), LDH,
  421:      $                   NHO, H( KV, KT ), LDH, NVE, H( KWV, 1 ), LDH,
  422:      $                   WORK, LWORK )
  423: *
  424: *           ==== Adjust KBOT accounting for new deflations. ====
  425: *
  426:             KBOT = KBOT - LD
  427: *
  428: *           ==== KS points to the shifts. ====
  429: *
  430:             KS = KBOT - LS + 1
  431: *
  432: *           ==== Skip an expensive QR sweep if there is a (partly
  433: *           .    heuristic) reason to expect that many eigenvalues
  434: *           .    will deflate without it.  Here, the QR sweep is
  435: *           .    skipped if many eigenvalues have just been deflated
  436: *           .    or if the remaining active block is small.
  437: *
  438:             IF( ( LD.EQ.0 ) .OR. ( ( 100*LD.LE.NW*NIBBLE ) .AND. ( KBOT-
  439:      $          KTOP+1.GT.MIN( NMIN, NWMAX ) ) ) ) THEN
  440: *
  441: *              ==== NS = nominal number of simultaneous shifts.
  442: *              .    This may be lowered (slightly) if DLAQR3
  443: *              .    did not provide that many shifts. ====
  444: *
  445:                NS = MIN( NSMAX, NSR, MAX( 2, KBOT-KTOP ) )
  446:                NS = NS - MOD( NS, 2 )
  447: *
  448: *              ==== If there have been no deflations
  449: *              .    in a multiple of KEXSH iterations,
  450: *              .    then try exceptional shifts.
  451: *              .    Otherwise use shifts provided by
  452: *              .    DLAQR3 above or from the eigenvalues
  453: *              .    of a trailing principal submatrix. ====
  454: *
  455:                IF( MOD( NDFL, KEXSH ).EQ.0 ) THEN
  456:                   KS = KBOT - NS + 1
  457:                   DO 30 I = KBOT, MAX( KS+1, KTOP+2 ), -2
  458:                      SS = ABS( H( I, I-1 ) ) + ABS( H( I-1, I-2 ) )
  459:                      AA = WILK1*SS + H( I, I )
  460:                      BB = SS
  461:                      CC = WILK2*SS
  462:                      DD = AA
  463:                      CALL DLANV2( AA, BB, CC, DD, WR( I-1 ), WI( I-1 ),
  464:      $                            WR( I ), WI( I ), CS, SN )
  465:    30             CONTINUE
  466:                   IF( KS.EQ.KTOP ) THEN
  467:                      WR( KS+1 ) = H( KS+1, KS+1 )
  468:                      WI( KS+1 ) = ZERO
  469:                      WR( KS ) = WR( KS+1 )
  470:                      WI( KS ) = WI( KS+1 )
  471:                   END IF
  472:                ELSE
  473: *
  474: *                 ==== Got NS/2 or fewer shifts? Use DLAQR4 or
  475: *                 .    DLAHQR on a trailing principal submatrix to
  476: *                 .    get more. (Since NS.LE.NSMAX.LE.(N+6)/9,
  477: *                 .    there is enough space below the subdiagonal
  478: *                 .    to fit an NS-by-NS scratch array.) ====
  479: *
  480:                   IF( KBOT-KS+1.LE.NS / 2 ) THEN
  481:                      KS = KBOT - NS + 1
  482:                      KT = N - NS + 1
  483:                      CALL DLACPY( 'A', NS, NS, H( KS, KS ), LDH,
  484:      $                            H( KT, 1 ), LDH )
  485:                      IF( NS.GT.NMIN ) THEN
  486:                         CALL DLAQR4( .false., .false., NS, 1, NS,
  487:      $                               H( KT, 1 ), LDH, WR( KS ),
  488:      $                               WI( KS ), 1, 1, ZDUM, 1, WORK,
  489:      $                               LWORK, INF )
  490:                      ELSE
  491:                         CALL DLAHQR( .false., .false., NS, 1, NS,
  492:      $                               H( KT, 1 ), LDH, WR( KS ),
  493:      $                               WI( KS ), 1, 1, ZDUM, 1, INF )
  494:                      END IF
  495:                      KS = KS + INF
  496: *
  497: *                    ==== In case of a rare QR failure use
  498: *                    .    eigenvalues of the trailing 2-by-2
  499: *                    .    principal submatrix.  ====
  500: *
  501:                      IF( KS.GE.KBOT ) THEN
  502:                         AA = H( KBOT-1, KBOT-1 )
  503:                         CC = H( KBOT, KBOT-1 )
  504:                         BB = H( KBOT-1, KBOT )
  505:                         DD = H( KBOT, KBOT )
  506:                         CALL DLANV2( AA, BB, CC, DD, WR( KBOT-1 ),
  507:      $                               WI( KBOT-1 ), WR( KBOT ),
  508:      $                               WI( KBOT ), CS, SN )
  509:                         KS = KBOT - 1
  510:                      END IF
  511:                   END IF
  512: *
  513:                   IF( KBOT-KS+1.GT.NS ) THEN
  514: *
  515: *                    ==== Sort the shifts (Helps a little)
  516: *                    .    Bubble sort keeps complex conjugate
  517: *                    .    pairs together. ====
  518: *
  519:                      SORTED = .false.
  520:                      DO 50 K = KBOT, KS + 1, -1
  521:                         IF( SORTED )
  522:      $                     GO TO 60
  523:                         SORTED = .true.
  524:                         DO 40 I = KS, K - 1
  525:                            IF( ABS( WR( I ) )+ABS( WI( I ) ).LT.
  526:      $                         ABS( WR( I+1 ) )+ABS( WI( I+1 ) ) ) THEN
  527:                               SORTED = .false.
  528: *
  529:                               SWAP = WR( I )
  530:                               WR( I ) = WR( I+1 )
  531:                               WR( I+1 ) = SWAP
  532: *
  533:                               SWAP = WI( I )
  534:                               WI( I ) = WI( I+1 )
  535:                               WI( I+1 ) = SWAP
  536:                            END IF
  537:    40                   CONTINUE
  538:    50                CONTINUE
  539:    60                CONTINUE
  540:                   END IF
  541: *
  542: *                 ==== Shuffle shifts into pairs of real shifts
  543: *                 .    and pairs of complex conjugate shifts
  544: *                 .    assuming complex conjugate shifts are
  545: *                 .    already adjacent to one another. (Yes,
  546: *                 .    they are.)  ====
  547: *
  548:                   DO 70 I = KBOT, KS + 2, -2
  549:                      IF( WI( I ).NE.-WI( I-1 ) ) THEN
  550: *
  551:                         SWAP = WR( I )
  552:                         WR( I ) = WR( I-1 )
  553:                         WR( I-1 ) = WR( I-2 )
  554:                         WR( I-2 ) = SWAP
  555: *
  556:                         SWAP = WI( I )
  557:                         WI( I ) = WI( I-1 )
  558:                         WI( I-1 ) = WI( I-2 )
  559:                         WI( I-2 ) = SWAP
  560:                      END IF
  561:    70             CONTINUE
  562:                END IF
  563: *
  564: *              ==== If there are only two shifts and both are
  565: *              .    real, then use only one.  ====
  566: *
  567:                IF( KBOT-KS+1.EQ.2 ) THEN
  568:                   IF( WI( KBOT ).EQ.ZERO ) THEN
  569:                      IF( ABS( WR( KBOT )-H( KBOT, KBOT ) ).LT.
  570:      $                   ABS( WR( KBOT-1 )-H( KBOT, KBOT ) ) ) THEN
  571:                         WR( KBOT-1 ) = WR( KBOT )
  572:                      ELSE
  573:                         WR( KBOT ) = WR( KBOT-1 )
  574:                      END IF
  575:                   END IF
  576:                END IF
  577: *
  578: *              ==== Use up to NS of the the smallest magnatiude
  579: *              .    shifts.  If there aren't NS shifts available,
  580: *              .    then use them all, possibly dropping one to
  581: *              .    make the number of shifts even. ====
  582: *
  583:                NS = MIN( NS, KBOT-KS+1 )
  584:                NS = NS - MOD( NS, 2 )
  585:                KS = KBOT - NS + 1
  586: *
  587: *              ==== Small-bulge multi-shift QR sweep:
  588: *              .    split workspace under the subdiagonal into
  589: *              .    - a KDU-by-KDU work array U in the lower
  590: *              .      left-hand-corner,
  591: *              .    - a KDU-by-at-least-KDU-but-more-is-better
  592: *              .      (KDU-by-NHo) horizontal work array WH along
  593: *              .      the bottom edge,
  594: *              .    - and an at-least-KDU-but-more-is-better-by-KDU
  595: *              .      (NVE-by-KDU) vertical work WV arrow along
  596: *              .      the left-hand-edge. ====
  597: *
  598:                KDU = 3*NS - 3
  599:                KU = N - KDU + 1
  600:                KWH = KDU + 1
  601:                NHO = ( N-KDU+1-4 ) - ( KDU+1 ) + 1
  602:                KWV = KDU + 4
  603:                NVE = N - KDU - KWV + 1
  604: *
  605: *              ==== Small-bulge multi-shift QR sweep ====
  606: *
  607:                CALL DLAQR5( WANTT, WANTZ, KACC22, N, KTOP, KBOT, NS,
  608:      $                      WR( KS ), WI( KS ), H, LDH, ILOZ, IHIZ, Z,
  609:      $                      LDZ, WORK, 3, H( KU, 1 ), LDH, NVE,
  610:      $                      H( KWV, 1 ), LDH, NHO, H( KU, KWH ), LDH )
  611:             END IF
  612: *
  613: *           ==== Note progress (or the lack of it). ====
  614: *
  615:             IF( LD.GT.0 ) THEN
  616:                NDFL = 1
  617:             ELSE
  618:                NDFL = NDFL + 1
  619:             END IF
  620: *
  621: *           ==== End of main loop ====
  622:    80    CONTINUE
  623: *
  624: *        ==== Iteration limit exceeded.  Set INFO to show where
  625: *        .    the problem occurred and exit. ====
  626: *
  627:          INFO = KBOT
  628:    90    CONTINUE
  629:       END IF
  630: *
  631: *     ==== Return the optimal value of LWORK. ====
  632: *
  633:       WORK( 1 ) = DBLE( LWKOPT )
  634: *
  635: *     ==== End of DLAQR0 ====
  636: *
  637:       END

CVSweb interface <joel.bertrand@systella.fr>