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