![]() |
LAPACK
3.4.1
LAPACK: Linear Algebra PACKage
|
00001 *> \brief \b DLAQR4 00002 * 00003 * =========== DOCUMENTATION =========== 00004 * 00005 * Online html documentation available at 00006 * http://www.netlib.org/lapack/explore-html/ 00007 * 00008 *> \htmlonly 00009 *> Download DLAQR4 + dependencies 00010 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaqr4.f"> 00011 *> [TGZ]</a> 00012 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaqr4.f"> 00013 *> [ZIP]</a> 00014 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaqr4.f"> 00015 *> [TXT]</a> 00016 *> \endhtmlonly 00017 * 00018 * Definition: 00019 * =========== 00020 * 00021 * SUBROUTINE DLAQR4( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, 00022 * ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO ) 00023 * 00024 * .. Scalar Arguments .. 00025 * INTEGER IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N 00026 * LOGICAL WANTT, WANTZ 00027 * .. 00028 * .. Array Arguments .. 00029 * DOUBLE PRECISION H( LDH, * ), WI( * ), WORK( * ), WR( * ), 00030 * $ Z( LDZ, * ) 00031 * .. 00032 * 00033 * 00034 *> \par Purpose: 00035 * ============= 00036 *> 00037 *> \verbatim 00038 *> 00039 *> DLAQR4 implements one level of recursion for DLAQR0. 00040 *> It is a complete implementation of the small bulge multi-shift 00041 *> QR algorithm. It may be called by DLAQR0 and, for large enough 00042 *> deflation window size, it may be called by DLAQR3. This 00043 *> subroutine is identical to DLAQR0 except that it calls DLAQR2 00044 *> instead of DLAQR3. 00045 *> 00046 *> DLAQR4 computes the eigenvalues of a Hessenberg matrix H 00047 *> and, optionally, the matrices T and Z from the Schur decomposition 00048 *> H = Z T Z**T, where T is an upper quasi-triangular matrix (the 00049 *> Schur form), and Z is the orthogonal matrix of Schur vectors. 00050 *> 00051 *> Optionally Z may be postmultiplied into an input orthogonal 00052 *> matrix Q so that this routine can give the Schur factorization 00053 *> of a matrix A which has been reduced to the Hessenberg form H 00054 *> by the orthogonal matrix Q: A = Q*H*Q**T = (QZ)*T*(QZ)**T. 00055 *> \endverbatim 00056 * 00057 * Arguments: 00058 * ========== 00059 * 00060 *> \param[in] WANTT 00061 *> \verbatim 00062 *> WANTT is LOGICAL 00063 *> = .TRUE. : the full Schur form T is required; 00064 *> = .FALSE.: only eigenvalues are required. 00065 *> \endverbatim 00066 *> 00067 *> \param[in] WANTZ 00068 *> \verbatim 00069 *> WANTZ is LOGICAL 00070 *> = .TRUE. : the matrix of Schur vectors Z is required; 00071 *> = .FALSE.: Schur vectors are not required. 00072 *> \endverbatim 00073 *> 00074 *> \param[in] N 00075 *> \verbatim 00076 *> N is INTEGER 00077 *> The order of the matrix H. N .GE. 0. 00078 *> \endverbatim 00079 *> 00080 *> \param[in] ILO 00081 *> \verbatim 00082 *> ILO is INTEGER 00083 *> \endverbatim 00084 *> 00085 *> \param[in] IHI 00086 *> \verbatim 00087 *> IHI is INTEGER 00088 *> It is assumed that H is already upper triangular in rows 00089 *> and columns 1:ILO-1 and IHI+1:N and, if ILO.GT.1, 00090 *> H(ILO,ILO-1) is zero. ILO and IHI are normally set by a 00091 *> previous call to DGEBAL, and then passed to DGEHRD when the 00092 *> matrix output by DGEBAL is reduced to Hessenberg form. 00093 *> Otherwise, ILO and IHI should be set to 1 and N, 00094 *> respectively. If N.GT.0, then 1.LE.ILO.LE.IHI.LE.N. 00095 *> If N = 0, then ILO = 1 and IHI = 0. 00096 *> \endverbatim 00097 *> 00098 *> \param[in,out] H 00099 *> \verbatim 00100 *> H is DOUBLE PRECISION array, dimension (LDH,N) 00101 *> On entry, the upper Hessenberg matrix H. 00102 *> On exit, if INFO = 0 and WANTT is .TRUE., then H contains 00103 *> the upper quasi-triangular matrix T from the Schur 00104 *> decomposition (the Schur form); 2-by-2 diagonal blocks 00105 *> (corresponding to complex conjugate pairs of eigenvalues) 00106 *> are returned in standard form, with H(i,i) = H(i+1,i+1) 00107 *> and H(i+1,i)*H(i,i+1).LT.0. If INFO = 0 and WANTT is 00108 *> .FALSE., then the contents of H are unspecified on exit. 00109 *> (The output value of H when INFO.GT.0 is given under the 00110 *> description of INFO below.) 00111 *> 00112 *> This subroutine may explicitly set H(i,j) = 0 for i.GT.j and 00113 *> j = 1, 2, ... ILO-1 or j = IHI+1, IHI+2, ... N. 00114 *> \endverbatim 00115 *> 00116 *> \param[in] LDH 00117 *> \verbatim 00118 *> LDH is INTEGER 00119 *> The leading dimension of the array H. LDH .GE. max(1,N). 00120 *> \endverbatim 00121 *> 00122 *> \param[out] WR 00123 *> \verbatim 00124 *> WR is DOUBLE PRECISION array, dimension (IHI) 00125 *> \endverbatim 00126 *> 00127 *> \param[out] WI 00128 *> \verbatim 00129 *> WI is DOUBLE PRECISION array, dimension (IHI) 00130 *> The real and imaginary parts, respectively, of the computed 00131 *> eigenvalues of H(ILO:IHI,ILO:IHI) are stored in WR(ILO:IHI) 00132 *> and WI(ILO:IHI). If two eigenvalues are computed as a 00133 *> complex conjugate pair, they are stored in consecutive 00134 *> elements of WR and WI, say the i-th and (i+1)th, with 00135 *> WI(i) .GT. 0 and WI(i+1) .LT. 0. If WANTT is .TRUE., then 00136 *> the eigenvalues are stored in the same order as on the 00137 *> diagonal of the Schur form returned in H, with 00138 *> WR(i) = H(i,i) and, if H(i:i+1,i:i+1) is a 2-by-2 diagonal 00139 *> block, WI(i) = sqrt(-H(i+1,i)*H(i,i+1)) and 00140 *> WI(i+1) = -WI(i). 00141 *> \endverbatim 00142 *> 00143 *> \param[in] ILOZ 00144 *> \verbatim 00145 *> ILOZ is INTEGER 00146 *> \endverbatim 00147 *> 00148 *> \param[in] IHIZ 00149 *> \verbatim 00150 *> IHIZ is INTEGER 00151 *> Specify the rows of Z to which transformations must be 00152 *> applied if WANTZ is .TRUE.. 00153 *> 1 .LE. ILOZ .LE. ILO; IHI .LE. IHIZ .LE. N. 00154 *> \endverbatim 00155 *> 00156 *> \param[in,out] Z 00157 *> \verbatim 00158 *> Z is DOUBLE PRECISION array, dimension (LDZ,IHI) 00159 *> If WANTZ is .FALSE., then Z is not referenced. 00160 *> If WANTZ is .TRUE., then Z(ILO:IHI,ILOZ:IHIZ) is 00161 *> replaced by Z(ILO:IHI,ILOZ:IHIZ)*U where U is the 00162 *> orthogonal Schur factor of H(ILO:IHI,ILO:IHI). 00163 *> (The output value of Z when INFO.GT.0 is given under 00164 *> the description of INFO below.) 00165 *> \endverbatim 00166 *> 00167 *> \param[in] LDZ 00168 *> \verbatim 00169 *> LDZ is INTEGER 00170 *> The leading dimension of the array Z. if WANTZ is .TRUE. 00171 *> then LDZ.GE.MAX(1,IHIZ). Otherwize, LDZ.GE.1. 00172 *> \endverbatim 00173 *> 00174 *> \param[out] WORK 00175 *> \verbatim 00176 *> WORK is DOUBLE PRECISION array, dimension LWORK 00177 *> On exit, if LWORK = -1, WORK(1) returns an estimate of 00178 *> the optimal value for LWORK. 00179 *> \endverbatim 00180 *> 00181 *> \param[in] LWORK 00182 *> \verbatim 00183 *> LWORK is INTEGER 00184 *> The dimension of the array WORK. LWORK .GE. max(1,N) 00185 *> is sufficient, but LWORK typically as large as 6*N may 00186 *> be required for optimal performance. A workspace query 00187 *> to determine the optimal workspace size is recommended. 00188 *> 00189 *> If LWORK = -1, then DLAQR4 does a workspace query. 00190 *> In this case, DLAQR4 checks the input parameters and 00191 *> estimates the optimal workspace size for the given 00192 *> values of N, ILO and IHI. The estimate is returned 00193 *> in WORK(1). No error message related to LWORK is 00194 *> issued by XERBLA. Neither H nor Z are accessed. 00195 *> \endverbatim 00196 *> 00197 *> \param[out] INFO 00198 *> \verbatim 00199 *> INFO is INTEGER 00200 *> = 0: successful exit 00201 *> .GT. 0: if INFO = i, DLAQR4 failed to compute all of 00202 *> the eigenvalues. Elements 1:ilo-1 and i+1:n of WR 00203 *> and WI contain those eigenvalues which have been 00204 *> successfully computed. (Failures are rare.) 00205 *> 00206 *> If INFO .GT. 0 and WANT is .FALSE., then on exit, 00207 *> the remaining unconverged eigenvalues are the eigen- 00208 *> values of the upper Hessenberg matrix rows and 00209 *> columns ILO through INFO of the final, output 00210 *> value of H. 00211 *> 00212 *> If INFO .GT. 0 and WANTT is .TRUE., then on exit 00213 *> 00214 *> (*) (initial value of H)*U = U*(final value of H) 00215 *> 00216 *> where U is a orthogonal matrix. The final 00217 *> value of H is upper Hessenberg and triangular in 00218 *> rows and columns INFO+1 through IHI. 00219 *> 00220 *> If INFO .GT. 0 and WANTZ is .TRUE., then on exit 00221 *> 00222 *> (final value of Z(ILO:IHI,ILOZ:IHIZ) 00223 *> = (initial value of Z(ILO:IHI,ILOZ:IHIZ)*U 00224 *> 00225 *> where U is the orthogonal matrix in (*) (regard- 00226 *> less of the value of WANTT.) 00227 *> 00228 *> If INFO .GT. 0 and WANTZ is .FALSE., then Z is not 00229 *> accessed. 00230 *> \endverbatim 00231 * 00232 * Authors: 00233 * ======== 00234 * 00235 *> \author Univ. of Tennessee 00236 *> \author Univ. of California Berkeley 00237 *> \author Univ. of Colorado Denver 00238 *> \author NAG Ltd. 00239 * 00240 *> \date November 2011 00241 * 00242 *> \ingroup doubleOTHERauxiliary 00243 * 00244 *> \par Contributors: 00245 * ================== 00246 *> 00247 *> Karen Braman and Ralph Byers, Department of Mathematics, 00248 *> University of Kansas, USA 00249 * 00250 *> \par References: 00251 * ================ 00252 *> 00253 *> K. Braman, R. Byers and R. Mathias, The Multi-Shift QR 00254 *> Algorithm Part I: Maintaining Well Focused Shifts, and Level 3 00255 *> Performance, SIAM Journal of Matrix Analysis, volume 23, pages 00256 *> 929--947, 2002. 00257 *> \n 00258 *> K. Braman, R. Byers and R. Mathias, The Multi-Shift QR 00259 *> Algorithm Part II: Aggressive Early Deflation, SIAM Journal 00260 *> of Matrix Analysis, volume 23, pages 948--973, 2002. 00261 *> 00262 * ===================================================================== 00263 SUBROUTINE DLAQR4( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, 00264 $ ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO ) 00265 * 00266 * -- LAPACK auxiliary routine (version 3.4.0) -- 00267 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00268 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00269 * November 2011 00270 * 00271 * .. Scalar Arguments .. 00272 INTEGER IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N 00273 LOGICAL WANTT, WANTZ 00274 * .. 00275 * .. Array Arguments .. 00276 DOUBLE PRECISION H( LDH, * ), WI( * ), WORK( * ), WR( * ), 00277 $ Z( LDZ, * ) 00278 * .. 00279 * 00280 * ================================================================ 00281 * .. Parameters .. 00282 * 00283 * ==== Matrices of order NTINY or smaller must be processed by 00284 * . DLAHQR because of insufficient subdiagonal scratch space. 00285 * . (This is a hard limit.) ==== 00286 INTEGER NTINY 00287 PARAMETER ( NTINY = 11 ) 00288 * 00289 * ==== Exceptional deflation windows: try to cure rare 00290 * . slow convergence by varying the size of the 00291 * . deflation window after KEXNW iterations. ==== 00292 INTEGER KEXNW 00293 PARAMETER ( KEXNW = 5 ) 00294 * 00295 * ==== Exceptional shifts: try to cure rare slow convergence 00296 * . with ad-hoc exceptional shifts every KEXSH iterations. 00297 * . ==== 00298 INTEGER KEXSH 00299 PARAMETER ( KEXSH = 6 ) 00300 * 00301 * ==== The constants WILK1 and WILK2 are used to form the 00302 * . exceptional shifts. ==== 00303 DOUBLE PRECISION WILK1, WILK2 00304 PARAMETER ( WILK1 = 0.75d0, WILK2 = -0.4375d0 ) 00305 DOUBLE PRECISION ZERO, ONE 00306 PARAMETER ( ZERO = 0.0d0, ONE = 1.0d0 ) 00307 * .. 00308 * .. Local Scalars .. 00309 DOUBLE PRECISION AA, BB, CC, CS, DD, SN, SS, SWAP 00310 INTEGER I, INF, IT, ITMAX, K, KACC22, KBOT, KDU, KS, 00311 $ KT, KTOP, KU, KV, KWH, KWTOP, KWV, LD, LS, 00312 $ LWKOPT, NDEC, NDFL, NH, NHO, NIBBLE, NMIN, NS, 00313 $ NSMAX, NSR, NVE, NW, NWMAX, NWR, NWUPBD 00314 LOGICAL SORTED 00315 CHARACTER JBCMPZ*2 00316 * .. 00317 * .. External Functions .. 00318 INTEGER ILAENV 00319 EXTERNAL ILAENV 00320 * .. 00321 * .. Local Arrays .. 00322 DOUBLE PRECISION ZDUM( 1, 1 ) 00323 * .. 00324 * .. External Subroutines .. 00325 EXTERNAL DLACPY, DLAHQR, DLANV2, DLAQR2, DLAQR5 00326 * .. 00327 * .. Intrinsic Functions .. 00328 INTRINSIC ABS, DBLE, INT, MAX, MIN, MOD 00329 * .. 00330 * .. Executable Statements .. 00331 INFO = 0 00332 * 00333 * ==== Quick return for N = 0: nothing to do. ==== 00334 * 00335 IF( N.EQ.0 ) THEN 00336 WORK( 1 ) = ONE 00337 RETURN 00338 END IF 00339 * 00340 IF( N.LE.NTINY ) THEN 00341 * 00342 * ==== Tiny matrices must use DLAHQR. ==== 00343 * 00344 LWKOPT = 1 00345 IF( LWORK.NE.-1 ) 00346 $ CALL DLAHQR( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, 00347 $ ILOZ, IHIZ, Z, LDZ, INFO ) 00348 ELSE 00349 * 00350 * ==== Use small bulge multi-shift QR with aggressive early 00351 * . deflation on larger-than-tiny matrices. ==== 00352 * 00353 * ==== Hope for the best. ==== 00354 * 00355 INFO = 0 00356 * 00357 * ==== Set up job flags for ILAENV. ==== 00358 * 00359 IF( WANTT ) THEN 00360 JBCMPZ( 1: 1 ) = 'S' 00361 ELSE 00362 JBCMPZ( 1: 1 ) = 'E' 00363 END IF 00364 IF( WANTZ ) THEN 00365 JBCMPZ( 2: 2 ) = 'V' 00366 ELSE 00367 JBCMPZ( 2: 2 ) = 'N' 00368 END IF 00369 * 00370 * ==== NWR = recommended deflation window size. At this 00371 * . point, N .GT. NTINY = 11, so there is enough 00372 * . subdiagonal workspace for NWR.GE.2 as required. 00373 * . (In fact, there is enough subdiagonal space for 00374 * . NWR.GE.3.) ==== 00375 * 00376 NWR = ILAENV( 13, 'DLAQR4', JBCMPZ, N, ILO, IHI, LWORK ) 00377 NWR = MAX( 2, NWR ) 00378 NWR = MIN( IHI-ILO+1, ( N-1 ) / 3, NWR ) 00379 * 00380 * ==== NSR = recommended number of simultaneous shifts. 00381 * . At this point N .GT. NTINY = 11, so there is at 00382 * . enough subdiagonal workspace for NSR to be even 00383 * . and greater than or equal to two as required. ==== 00384 * 00385 NSR = ILAENV( 15, 'DLAQR4', JBCMPZ, N, ILO, IHI, LWORK ) 00386 NSR = MIN( NSR, ( N+6 ) / 9, IHI-ILO ) 00387 NSR = MAX( 2, NSR-MOD( NSR, 2 ) ) 00388 * 00389 * ==== Estimate optimal workspace ==== 00390 * 00391 * ==== Workspace query call to DLAQR2 ==== 00392 * 00393 CALL DLAQR2( WANTT, WANTZ, N, ILO, IHI, NWR+1, H, LDH, ILOZ, 00394 $ IHIZ, Z, LDZ, LS, LD, WR, WI, H, LDH, N, H, LDH, 00395 $ N, H, LDH, WORK, -1 ) 00396 * 00397 * ==== Optimal workspace = MAX(DLAQR5, DLAQR2) ==== 00398 * 00399 LWKOPT = MAX( 3*NSR / 2, INT( WORK( 1 ) ) ) 00400 * 00401 * ==== Quick return in case of workspace query. ==== 00402 * 00403 IF( LWORK.EQ.-1 ) THEN 00404 WORK( 1 ) = DBLE( LWKOPT ) 00405 RETURN 00406 END IF 00407 * 00408 * ==== DLAHQR/DLAQR0 crossover point ==== 00409 * 00410 NMIN = ILAENV( 12, 'DLAQR4', JBCMPZ, N, ILO, IHI, LWORK ) 00411 NMIN = MAX( NTINY, NMIN ) 00412 * 00413 * ==== Nibble crossover point ==== 00414 * 00415 NIBBLE = ILAENV( 14, 'DLAQR4', JBCMPZ, N, ILO, IHI, LWORK ) 00416 NIBBLE = MAX( 0, NIBBLE ) 00417 * 00418 * ==== Accumulate reflections during ttswp? Use block 00419 * . 2-by-2 structure during matrix-matrix multiply? ==== 00420 * 00421 KACC22 = ILAENV( 16, 'DLAQR4', JBCMPZ, N, ILO, IHI, LWORK ) 00422 KACC22 = MAX( 0, KACC22 ) 00423 KACC22 = MIN( 2, KACC22 ) 00424 * 00425 * ==== NWMAX = the largest possible deflation window for 00426 * . which there is sufficient workspace. ==== 00427 * 00428 NWMAX = MIN( ( N-1 ) / 3, LWORK / 2 ) 00429 NW = NWMAX 00430 * 00431 * ==== NSMAX = the Largest number of simultaneous shifts 00432 * . for which there is sufficient workspace. ==== 00433 * 00434 NSMAX = MIN( ( N+6 ) / 9, 2*LWORK / 3 ) 00435 NSMAX = NSMAX - MOD( NSMAX, 2 ) 00436 * 00437 * ==== NDFL: an iteration count restarted at deflation. ==== 00438 * 00439 NDFL = 1 00440 * 00441 * ==== ITMAX = iteration limit ==== 00442 * 00443 ITMAX = MAX( 30, 2*KEXSH )*MAX( 10, ( IHI-ILO+1 ) ) 00444 * 00445 * ==== Last row and column in the active block ==== 00446 * 00447 KBOT = IHI 00448 * 00449 * ==== Main Loop ==== 00450 * 00451 DO 80 IT = 1, ITMAX 00452 * 00453 * ==== Done when KBOT falls below ILO ==== 00454 * 00455 IF( KBOT.LT.ILO ) 00456 $ GO TO 90 00457 * 00458 * ==== Locate active block ==== 00459 * 00460 DO 10 K = KBOT, ILO + 1, -1 00461 IF( H( K, K-1 ).EQ.ZERO ) 00462 $ GO TO 20 00463 10 CONTINUE 00464 K = ILO 00465 20 CONTINUE 00466 KTOP = K 00467 * 00468 * ==== Select deflation window size: 00469 * . Typical Case: 00470 * . If possible and advisable, nibble the entire 00471 * . active block. If not, use size MIN(NWR,NWMAX) 00472 * . or MIN(NWR+1,NWMAX) depending upon which has 00473 * . the smaller corresponding subdiagonal entry 00474 * . (a heuristic). 00475 * . 00476 * . Exceptional Case: 00477 * . If there have been no deflations in KEXNW or 00478 * . more iterations, then vary the deflation window 00479 * . size. At first, because, larger windows are, 00480 * . in general, more powerful than smaller ones, 00481 * . rapidly increase the window to the maximum possible. 00482 * . Then, gradually reduce the window size. ==== 00483 * 00484 NH = KBOT - KTOP + 1 00485 NWUPBD = MIN( NH, NWMAX ) 00486 IF( NDFL.LT.KEXNW ) THEN 00487 NW = MIN( NWUPBD, NWR ) 00488 ELSE 00489 NW = MIN( NWUPBD, 2*NW ) 00490 END IF 00491 IF( NW.LT.NWMAX ) THEN 00492 IF( NW.GE.NH-1 ) THEN 00493 NW = NH 00494 ELSE 00495 KWTOP = KBOT - NW + 1 00496 IF( ABS( H( KWTOP, KWTOP-1 ) ).GT. 00497 $ ABS( H( KWTOP-1, KWTOP-2 ) ) )NW = NW + 1 00498 END IF 00499 END IF 00500 IF( NDFL.LT.KEXNW ) THEN 00501 NDEC = -1 00502 ELSE IF( NDEC.GE.0 .OR. NW.GE.NWUPBD ) THEN 00503 NDEC = NDEC + 1 00504 IF( NW-NDEC.LT.2 ) 00505 $ NDEC = 0 00506 NW = NW - NDEC 00507 END IF 00508 * 00509 * ==== Aggressive early deflation: 00510 * . split workspace under the subdiagonal into 00511 * . - an nw-by-nw work array V in the lower 00512 * . left-hand-corner, 00513 * . - an NW-by-at-least-NW-but-more-is-better 00514 * . (NW-by-NHO) horizontal work array along 00515 * . the bottom edge, 00516 * . - an at-least-NW-but-more-is-better (NHV-by-NW) 00517 * . vertical work array along the left-hand-edge. 00518 * . ==== 00519 * 00520 KV = N - NW + 1 00521 KT = NW + 1 00522 NHO = ( N-NW-1 ) - KT + 1 00523 KWV = NW + 2 00524 NVE = ( N-NW ) - KWV + 1 00525 * 00526 * ==== Aggressive early deflation ==== 00527 * 00528 CALL DLAQR2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ, 00529 $ IHIZ, Z, LDZ, LS, LD, WR, WI, H( KV, 1 ), LDH, 00530 $ NHO, H( KV, KT ), LDH, NVE, H( KWV, 1 ), LDH, 00531 $ WORK, LWORK ) 00532 * 00533 * ==== Adjust KBOT accounting for new deflations. ==== 00534 * 00535 KBOT = KBOT - LD 00536 * 00537 * ==== KS points to the shifts. ==== 00538 * 00539 KS = KBOT - LS + 1 00540 * 00541 * ==== Skip an expensive QR sweep if there is a (partly 00542 * . heuristic) reason to expect that many eigenvalues 00543 * . will deflate without it. Here, the QR sweep is 00544 * . skipped if many eigenvalues have just been deflated 00545 * . or if the remaining active block is small. 00546 * 00547 IF( ( LD.EQ.0 ) .OR. ( ( 100*LD.LE.NW*NIBBLE ) .AND. ( KBOT- 00548 $ KTOP+1.GT.MIN( NMIN, NWMAX ) ) ) ) THEN 00549 * 00550 * ==== NS = nominal number of simultaneous shifts. 00551 * . This may be lowered (slightly) if DLAQR2 00552 * . did not provide that many shifts. ==== 00553 * 00554 NS = MIN( NSMAX, NSR, MAX( 2, KBOT-KTOP ) ) 00555 NS = NS - MOD( NS, 2 ) 00556 * 00557 * ==== If there have been no deflations 00558 * . in a multiple of KEXSH iterations, 00559 * . then try exceptional shifts. 00560 * . Otherwise use shifts provided by 00561 * . DLAQR2 above or from the eigenvalues 00562 * . of a trailing principal submatrix. ==== 00563 * 00564 IF( MOD( NDFL, KEXSH ).EQ.0 ) THEN 00565 KS = KBOT - NS + 1 00566 DO 30 I = KBOT, MAX( KS+1, KTOP+2 ), -2 00567 SS = ABS( H( I, I-1 ) ) + ABS( H( I-1, I-2 ) ) 00568 AA = WILK1*SS + H( I, I ) 00569 BB = SS 00570 CC = WILK2*SS 00571 DD = AA 00572 CALL DLANV2( AA, BB, CC, DD, WR( I-1 ), WI( I-1 ), 00573 $ WR( I ), WI( I ), CS, SN ) 00574 30 CONTINUE 00575 IF( KS.EQ.KTOP ) THEN 00576 WR( KS+1 ) = H( KS+1, KS+1 ) 00577 WI( KS+1 ) = ZERO 00578 WR( KS ) = WR( KS+1 ) 00579 WI( KS ) = WI( KS+1 ) 00580 END IF 00581 ELSE 00582 * 00583 * ==== Got NS/2 or fewer shifts? Use DLAHQR 00584 * . on a trailing principal submatrix to 00585 * . get more. (Since NS.LE.NSMAX.LE.(N+6)/9, 00586 * . there is enough space below the subdiagonal 00587 * . to fit an NS-by-NS scratch array.) ==== 00588 * 00589 IF( KBOT-KS+1.LE.NS / 2 ) THEN 00590 KS = KBOT - NS + 1 00591 KT = N - NS + 1 00592 CALL DLACPY( 'A', NS, NS, H( KS, KS ), LDH, 00593 $ H( KT, 1 ), LDH ) 00594 CALL DLAHQR( .false., .false., NS, 1, NS, 00595 $ H( KT, 1 ), LDH, WR( KS ), WI( KS ), 00596 $ 1, 1, ZDUM, 1, INF ) 00597 KS = KS + INF 00598 * 00599 * ==== In case of a rare QR failure use 00600 * . eigenvalues of the trailing 2-by-2 00601 * . principal submatrix. ==== 00602 * 00603 IF( KS.GE.KBOT ) THEN 00604 AA = H( KBOT-1, KBOT-1 ) 00605 CC = H( KBOT, KBOT-1 ) 00606 BB = H( KBOT-1, KBOT ) 00607 DD = H( KBOT, KBOT ) 00608 CALL DLANV2( AA, BB, CC, DD, WR( KBOT-1 ), 00609 $ WI( KBOT-1 ), WR( KBOT ), 00610 $ WI( KBOT ), CS, SN ) 00611 KS = KBOT - 1 00612 END IF 00613 END IF 00614 * 00615 IF( KBOT-KS+1.GT.NS ) THEN 00616 * 00617 * ==== Sort the shifts (Helps a little) 00618 * . Bubble sort keeps complex conjugate 00619 * . pairs together. ==== 00620 * 00621 SORTED = .false. 00622 DO 50 K = KBOT, KS + 1, -1 00623 IF( SORTED ) 00624 $ GO TO 60 00625 SORTED = .true. 00626 DO 40 I = KS, K - 1 00627 IF( ABS( WR( I ) )+ABS( WI( I ) ).LT. 00628 $ ABS( WR( I+1 ) )+ABS( WI( I+1 ) ) ) THEN 00629 SORTED = .false. 00630 * 00631 SWAP = WR( I ) 00632 WR( I ) = WR( I+1 ) 00633 WR( I+1 ) = SWAP 00634 * 00635 SWAP = WI( I ) 00636 WI( I ) = WI( I+1 ) 00637 WI( I+1 ) = SWAP 00638 END IF 00639 40 CONTINUE 00640 50 CONTINUE 00641 60 CONTINUE 00642 END IF 00643 * 00644 * ==== Shuffle shifts into pairs of real shifts 00645 * . and pairs of complex conjugate shifts 00646 * . assuming complex conjugate shifts are 00647 * . already adjacent to one another. (Yes, 00648 * . they are.) ==== 00649 * 00650 DO 70 I = KBOT, KS + 2, -2 00651 IF( WI( I ).NE.-WI( I-1 ) ) THEN 00652 * 00653 SWAP = WR( I ) 00654 WR( I ) = WR( I-1 ) 00655 WR( I-1 ) = WR( I-2 ) 00656 WR( I-2 ) = SWAP 00657 * 00658 SWAP = WI( I ) 00659 WI( I ) = WI( I-1 ) 00660 WI( I-1 ) = WI( I-2 ) 00661 WI( I-2 ) = SWAP 00662 END IF 00663 70 CONTINUE 00664 END IF 00665 * 00666 * ==== If there are only two shifts and both are 00667 * . real, then use only one. ==== 00668 * 00669 IF( KBOT-KS+1.EQ.2 ) THEN 00670 IF( WI( KBOT ).EQ.ZERO ) THEN 00671 IF( ABS( WR( KBOT )-H( KBOT, KBOT ) ).LT. 00672 $ ABS( WR( KBOT-1 )-H( KBOT, KBOT ) ) ) THEN 00673 WR( KBOT-1 ) = WR( KBOT ) 00674 ELSE 00675 WR( KBOT ) = WR( KBOT-1 ) 00676 END IF 00677 END IF 00678 END IF 00679 * 00680 * ==== Use up to NS of the the smallest magnatiude 00681 * . shifts. If there aren't NS shifts available, 00682 * . then use them all, possibly dropping one to 00683 * . make the number of shifts even. ==== 00684 * 00685 NS = MIN( NS, KBOT-KS+1 ) 00686 NS = NS - MOD( NS, 2 ) 00687 KS = KBOT - NS + 1 00688 * 00689 * ==== Small-bulge multi-shift QR sweep: 00690 * . split workspace under the subdiagonal into 00691 * . - a KDU-by-KDU work array U in the lower 00692 * . left-hand-corner, 00693 * . - a KDU-by-at-least-KDU-but-more-is-better 00694 * . (KDU-by-NHo) horizontal work array WH along 00695 * . the bottom edge, 00696 * . - and an at-least-KDU-but-more-is-better-by-KDU 00697 * . (NVE-by-KDU) vertical work WV arrow along 00698 * . the left-hand-edge. ==== 00699 * 00700 KDU = 3*NS - 3 00701 KU = N - KDU + 1 00702 KWH = KDU + 1 00703 NHO = ( N-KDU+1-4 ) - ( KDU+1 ) + 1 00704 KWV = KDU + 4 00705 NVE = N - KDU - KWV + 1 00706 * 00707 * ==== Small-bulge multi-shift QR sweep ==== 00708 * 00709 CALL DLAQR5( WANTT, WANTZ, KACC22, N, KTOP, KBOT, NS, 00710 $ WR( KS ), WI( KS ), H, LDH, ILOZ, IHIZ, Z, 00711 $ LDZ, WORK, 3, H( KU, 1 ), LDH, NVE, 00712 $ H( KWV, 1 ), LDH, NHO, H( KU, KWH ), LDH ) 00713 END IF 00714 * 00715 * ==== Note progress (or the lack of it). ==== 00716 * 00717 IF( LD.GT.0 ) THEN 00718 NDFL = 1 00719 ELSE 00720 NDFL = NDFL + 1 00721 END IF 00722 * 00723 * ==== End of main loop ==== 00724 80 CONTINUE 00725 * 00726 * ==== Iteration limit exceeded. Set INFO to show where 00727 * . the problem occurred and exit. ==== 00728 * 00729 INFO = KBOT 00730 90 CONTINUE 00731 END IF 00732 * 00733 * ==== Return the optimal value of LWORK. ==== 00734 * 00735 WORK( 1 ) = DBLE( LWKOPT ) 00736 * 00737 * ==== End of DLAQR4 ==== 00738 * 00739 END