LAPACK  3.9.1
LAPACK: Linear Algebra PACKage
zlaqr2.f
Go to the documentation of this file.
1 *> \brief \b ZLAQR2 performs the unitary similarity transformation of a Hessenberg matrix to detect and deflate fully converged eigenvalues from a trailing principal submatrix (aggressive early deflation).
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZLAQR2 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqr2.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqr2.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqr2.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZLAQR2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
22 * IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT,
23 * NV, WV, LDWV, WORK, LWORK )
24 *
25 * .. Scalar Arguments ..
26 * INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
27 * $ LDZ, LWORK, N, ND, NH, NS, NV, NW
28 * LOGICAL WANTT, WANTZ
29 * ..
30 * .. Array Arguments ..
31 * COMPLEX*16 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
32 * $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
33 * ..
34 *
35 *
36 *> \par Purpose:
37 * =============
38 *>
39 *> \verbatim
40 *>
41 *> ZLAQR2 is identical to ZLAQR3 except that it avoids
42 *> recursion by calling ZLAHQR instead of ZLAQR4.
43 *>
44 *> Aggressive early deflation:
45 *>
46 *> ZLAQR2 accepts as input an upper Hessenberg matrix
47 *> H and performs an unitary similarity transformation
48 *> designed to detect and deflate fully converged eigenvalues from
49 *> a trailing principal submatrix. On output H has been over-
50 *> written by a new Hessenberg matrix that is a perturbation of
51 *> an unitary similarity transformation of H. It is to be
52 *> hoped that the final version of H has many zero subdiagonal
53 *> entries.
54 *>
55 *> \endverbatim
56 *
57 * Arguments:
58 * ==========
59 *
60 *> \param[in] WANTT
61 *> \verbatim
62 *> WANTT is LOGICAL
63 *> If .TRUE., then the Hessenberg matrix H is fully updated
64 *> so that the triangular Schur factor may be
65 *> computed (in cooperation with the calling subroutine).
66 *> If .FALSE., then only enough of H is updated to preserve
67 *> the eigenvalues.
68 *> \endverbatim
69 *>
70 *> \param[in] WANTZ
71 *> \verbatim
72 *> WANTZ is LOGICAL
73 *> If .TRUE., then the unitary matrix Z is updated so
74 *> so that the unitary Schur factor may be computed
75 *> (in cooperation with the calling subroutine).
76 *> If .FALSE., then Z is not referenced.
77 *> \endverbatim
78 *>
79 *> \param[in] N
80 *> \verbatim
81 *> N is INTEGER
82 *> The order of the matrix H and (if WANTZ is .TRUE.) the
83 *> order of the unitary matrix Z.
84 *> \endverbatim
85 *>
86 *> \param[in] KTOP
87 *> \verbatim
88 *> KTOP is INTEGER
89 *> It is assumed that either KTOP = 1 or H(KTOP,KTOP-1)=0.
90 *> KBOT and KTOP together determine an isolated block
91 *> along the diagonal of the Hessenberg matrix.
92 *> \endverbatim
93 *>
94 *> \param[in] KBOT
95 *> \verbatim
96 *> KBOT is INTEGER
97 *> It is assumed without a check that either
98 *> KBOT = N or H(KBOT+1,KBOT)=0. KBOT and KTOP together
99 *> determine an isolated block along the diagonal of the
100 *> Hessenberg matrix.
101 *> \endverbatim
102 *>
103 *> \param[in] NW
104 *> \verbatim
105 *> NW is INTEGER
106 *> Deflation window size. 1 <= NW <= (KBOT-KTOP+1).
107 *> \endverbatim
108 *>
109 *> \param[in,out] H
110 *> \verbatim
111 *> H is COMPLEX*16 array, dimension (LDH,N)
112 *> On input the initial N-by-N section of H stores the
113 *> Hessenberg matrix undergoing aggressive early deflation.
114 *> On output H has been transformed by a unitary
115 *> similarity transformation, perturbed, and the returned
116 *> to Hessenberg form that (it is to be hoped) has some
117 *> zero subdiagonal entries.
118 *> \endverbatim
119 *>
120 *> \param[in] LDH
121 *> \verbatim
122 *> LDH is INTEGER
123 *> Leading dimension of H just as declared in the calling
124 *> subroutine. N <= LDH
125 *> \endverbatim
126 *>
127 *> \param[in] ILOZ
128 *> \verbatim
129 *> ILOZ is INTEGER
130 *> \endverbatim
131 *>
132 *> \param[in] IHIZ
133 *> \verbatim
134 *> IHIZ is INTEGER
135 *> Specify the rows of Z to which transformations must be
136 *> applied if WANTZ is .TRUE.. 1 <= ILOZ <= IHIZ <= N.
137 *> \endverbatim
138 *>
139 *> \param[in,out] Z
140 *> \verbatim
141 *> Z is COMPLEX*16 array, dimension (LDZ,N)
142 *> IF WANTZ is .TRUE., then on output, the unitary
143 *> similarity transformation mentioned above has been
144 *> accumulated into Z(ILOZ:IHIZ,ILOZ:IHIZ) from the right.
145 *> If WANTZ is .FALSE., then Z is unreferenced.
146 *> \endverbatim
147 *>
148 *> \param[in] LDZ
149 *> \verbatim
150 *> LDZ is INTEGER
151 *> The leading dimension of Z just as declared in the
152 *> calling subroutine. 1 <= LDZ.
153 *> \endverbatim
154 *>
155 *> \param[out] NS
156 *> \verbatim
157 *> NS is INTEGER
158 *> The number of unconverged (ie approximate) eigenvalues
159 *> returned in SR and SI that may be used as shifts by the
160 *> calling subroutine.
161 *> \endverbatim
162 *>
163 *> \param[out] ND
164 *> \verbatim
165 *> ND is INTEGER
166 *> The number of converged eigenvalues uncovered by this
167 *> subroutine.
168 *> \endverbatim
169 *>
170 *> \param[out] SH
171 *> \verbatim
172 *> SH is COMPLEX*16 array, dimension (KBOT)
173 *> On output, approximate eigenvalues that may
174 *> be used for shifts are stored in SH(KBOT-ND-NS+1)
175 *> through SR(KBOT-ND). Converged eigenvalues are
176 *> stored in SH(KBOT-ND+1) through SH(KBOT).
177 *> \endverbatim
178 *>
179 *> \param[out] V
180 *> \verbatim
181 *> V is COMPLEX*16 array, dimension (LDV,NW)
182 *> An NW-by-NW work array.
183 *> \endverbatim
184 *>
185 *> \param[in] LDV
186 *> \verbatim
187 *> LDV is INTEGER
188 *> The leading dimension of V just as declared in the
189 *> calling subroutine. NW <= LDV
190 *> \endverbatim
191 *>
192 *> \param[in] NH
193 *> \verbatim
194 *> NH is INTEGER
195 *> The number of columns of T. NH >= NW.
196 *> \endverbatim
197 *>
198 *> \param[out] T
199 *> \verbatim
200 *> T is COMPLEX*16 array, dimension (LDT,NW)
201 *> \endverbatim
202 *>
203 *> \param[in] LDT
204 *> \verbatim
205 *> LDT is INTEGER
206 *> The leading dimension of T just as declared in the
207 *> calling subroutine. NW <= LDT
208 *> \endverbatim
209 *>
210 *> \param[in] NV
211 *> \verbatim
212 *> NV is INTEGER
213 *> The number of rows of work array WV available for
214 *> workspace. NV >= NW.
215 *> \endverbatim
216 *>
217 *> \param[out] WV
218 *> \verbatim
219 *> WV is COMPLEX*16 array, dimension (LDWV,NW)
220 *> \endverbatim
221 *>
222 *> \param[in] LDWV
223 *> \verbatim
224 *> LDWV is INTEGER
225 *> The leading dimension of W just as declared in the
226 *> calling subroutine. NW <= LDV
227 *> \endverbatim
228 *>
229 *> \param[out] WORK
230 *> \verbatim
231 *> WORK is COMPLEX*16 array, dimension (LWORK)
232 *> On exit, WORK(1) is set to an estimate of the optimal value
233 *> of LWORK for the given values of N, NW, KTOP and KBOT.
234 *> \endverbatim
235 *>
236 *> \param[in] LWORK
237 *> \verbatim
238 *> LWORK is INTEGER
239 *> The dimension of the work array WORK. LWORK = 2*NW
240 *> suffices, but greater efficiency may result from larger
241 *> values of LWORK.
242 *>
243 *> If LWORK = -1, then a workspace query is assumed; ZLAQR2
244 *> only estimates the optimal workspace size for the given
245 *> values of N, NW, KTOP and KBOT. The estimate is returned
246 *> in WORK(1). No error message related to LWORK is issued
247 *> by XERBLA. Neither H nor Z are accessed.
248 *> \endverbatim
249 *
250 * Authors:
251 * ========
252 *
253 *> \author Univ. of Tennessee
254 *> \author Univ. of California Berkeley
255 *> \author Univ. of Colorado Denver
256 *> \author NAG Ltd.
257 *
258 *> \ingroup complex16OTHERauxiliary
259 *
260 *> \par Contributors:
261 * ==================
262 *>
263 *> Karen Braman and Ralph Byers, Department of Mathematics,
264 *> University of Kansas, USA
265 *>
266 * =====================================================================
267  SUBROUTINE zlaqr2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
268  $ IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT,
269  $ NV, WV, LDWV, WORK, LWORK )
270 *
271 * -- LAPACK auxiliary routine --
272 * -- LAPACK is a software package provided by Univ. of Tennessee, --
273 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
274 *
275 * .. Scalar Arguments ..
276  INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
277  $ LDZ, LWORK, N, ND, NH, NS, NV, NW
278  LOGICAL WANTT, WANTZ
279 * ..
280 * .. Array Arguments ..
281  COMPLEX*16 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
282  $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
283 * ..
284 *
285 * ================================================================
286 *
287 * .. Parameters ..
288  COMPLEX*16 ZERO, ONE
289  PARAMETER ( ZERO = ( 0.0d0, 0.0d0 ),
290  $ one = ( 1.0d0, 0.0d0 ) )
291  DOUBLE PRECISION RZERO, RONE
292  PARAMETER ( RZERO = 0.0d0, rone = 1.0d0 )
293 * ..
294 * .. Local Scalars ..
295  COMPLEX*16 BETA, CDUM, S, TAU
296  DOUBLE PRECISION FOO, SAFMAX, SAFMIN, SMLNUM, ULP
297  INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
298  $ knt, krow, kwtop, ltop, lwk1, lwk2, lwkopt
299 * ..
300 * .. External Functions ..
301  DOUBLE PRECISION DLAMCH
302  EXTERNAL DLAMCH
303 * ..
304 * .. External Subroutines ..
305  EXTERNAL dlabad, zcopy, zgehrd, zgemm, zlacpy, zlahqr,
307 * ..
308 * .. Intrinsic Functions ..
309  INTRINSIC abs, dble, dcmplx, dconjg, dimag, int, max, min
310 * ..
311 * .. Statement Functions ..
312  DOUBLE PRECISION CABS1
313 * ..
314 * .. Statement Function definitions ..
315  cabs1( cdum ) = abs( dble( cdum ) ) + abs( dimag( cdum ) )
316 * ..
317 * .. Executable Statements ..
318 *
319 * ==== Estimate optimal workspace. ====
320 *
321  jw = min( nw, kbot-ktop+1 )
322  IF( jw.LE.2 ) THEN
323  lwkopt = 1
324  ELSE
325 *
326 * ==== Workspace query call to ZGEHRD ====
327 *
328  CALL zgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
329  lwk1 = int( work( 1 ) )
330 *
331 * ==== Workspace query call to ZUNMHR ====
332 *
333  CALL zunmhr( 'R', 'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
334  $ work, -1, info )
335  lwk2 = int( work( 1 ) )
336 *
337 * ==== Optimal workspace ====
338 *
339  lwkopt = jw + max( lwk1, lwk2 )
340  END IF
341 *
342 * ==== Quick return in case of workspace query. ====
343 *
344  IF( lwork.EQ.-1 ) THEN
345  work( 1 ) = dcmplx( lwkopt, 0 )
346  RETURN
347  END IF
348 *
349 * ==== Nothing to do ...
350 * ... for an empty active block ... ====
351  ns = 0
352  nd = 0
353  work( 1 ) = one
354  IF( ktop.GT.kbot )
355  $ RETURN
356 * ... nor for an empty deflation window. ====
357  IF( nw.LT.1 )
358  $ RETURN
359 *
360 * ==== Machine constants ====
361 *
362  safmin = dlamch( 'SAFE MINIMUM' )
363  safmax = rone / safmin
364  CALL dlabad( safmin, safmax )
365  ulp = dlamch( 'PRECISION' )
366  smlnum = safmin*( dble( n ) / ulp )
367 *
368 * ==== Setup deflation window ====
369 *
370  jw = min( nw, kbot-ktop+1 )
371  kwtop = kbot - jw + 1
372  IF( kwtop.EQ.ktop ) THEN
373  s = zero
374  ELSE
375  s = h( kwtop, kwtop-1 )
376  END IF
377 *
378  IF( kbot.EQ.kwtop ) THEN
379 *
380 * ==== 1-by-1 deflation window: not much to do ====
381 *
382  sh( kwtop ) = h( kwtop, kwtop )
383  ns = 1
384  nd = 0
385  IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
386  $ kwtop ) ) ) ) THEN
387  ns = 0
388  nd = 1
389  IF( kwtop.GT.ktop )
390  $ h( kwtop, kwtop-1 ) = zero
391  END IF
392  work( 1 ) = one
393  RETURN
394  END IF
395 *
396 * ==== Convert to spike-triangular form. (In case of a
397 * . rare QR failure, this routine continues to do
398 * . aggressive early deflation using that part of
399 * . the deflation window that converged using INFQR
400 * . here and there to keep track.) ====
401 *
402  CALL zlacpy( 'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
403  CALL zcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
404 *
405  CALL zlaset( 'A', jw, jw, zero, one, v, ldv )
406  CALL zlahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
407  $ jw, v, ldv, infqr )
408 *
409 * ==== Deflation detection loop ====
410 *
411  ns = jw
412  ilst = infqr + 1
413  DO 10 knt = infqr + 1, jw
414 *
415 * ==== Small spike tip deflation test ====
416 *
417  foo = cabs1( t( ns, ns ) )
418  IF( foo.EQ.rzero )
419  $ foo = cabs1( s )
420  IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
421  $ THEN
422 *
423 * ==== One more converged eigenvalue ====
424 *
425  ns = ns - 1
426  ELSE
427 *
428 * ==== One undeflatable eigenvalue. Move it up out of the
429 * . way. (ZTREXC can not fail in this case.) ====
430 *
431  ifst = ns
432  CALL ztrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, info )
433  ilst = ilst + 1
434  END IF
435  10 CONTINUE
436 *
437 * ==== Return to Hessenberg form ====
438 *
439  IF( ns.EQ.0 )
440  $ s = zero
441 *
442  IF( ns.LT.jw ) THEN
443 *
444 * ==== sorting the diagonal of T improves accuracy for
445 * . graded matrices. ====
446 *
447  DO 30 i = infqr + 1, ns
448  ifst = i
449  DO 20 j = i + 1, ns
450  IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
451  $ ifst = j
452  20 CONTINUE
453  ilst = i
454  IF( ifst.NE.ilst )
455  $ CALL ztrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, info )
456  30 CONTINUE
457  END IF
458 *
459 * ==== Restore shift/eigenvalue array from T ====
460 *
461  DO 40 i = infqr + 1, jw
462  sh( kwtop+i-1 ) = t( i, i )
463  40 CONTINUE
464 *
465 *
466  IF( ns.LT.jw .OR. s.EQ.zero ) THEN
467  IF( ns.GT.1 .AND. s.NE.zero ) THEN
468 *
469 * ==== Reflect spike back into lower triangle ====
470 *
471  CALL zcopy( ns, v, ldv, work, 1 )
472  DO 50 i = 1, ns
473  work( i ) = dconjg( work( i ) )
474  50 CONTINUE
475  beta = work( 1 )
476  CALL zlarfg( ns, beta, work( 2 ), 1, tau )
477  work( 1 ) = one
478 *
479  CALL zlaset( 'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
480 *
481  CALL zlarf( 'L', ns, jw, work, 1, dconjg( tau ), t, ldt,
482  $ work( jw+1 ) )
483  CALL zlarf( 'R', ns, ns, work, 1, tau, t, ldt,
484  $ work( jw+1 ) )
485  CALL zlarf( 'R', jw, ns, work, 1, tau, v, ldv,
486  $ work( jw+1 ) )
487 *
488  CALL zgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
489  $ lwork-jw, info )
490  END IF
491 *
492 * ==== Copy updated reduced window into place ====
493 *
494  IF( kwtop.GT.1 )
495  $ h( kwtop, kwtop-1 ) = s*dconjg( v( 1, 1 ) )
496  CALL zlacpy( 'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
497  CALL zcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
498  $ ldh+1 )
499 *
500 * ==== Accumulate orthogonal matrix in order update
501 * . H and Z, if requested. ====
502 *
503  IF( ns.GT.1 .AND. s.NE.zero )
504  $ CALL zunmhr( 'R', 'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
505  $ work( jw+1 ), lwork-jw, info )
506 *
507 * ==== Update vertical slab in H ====
508 *
509  IF( wantt ) THEN
510  ltop = 1
511  ELSE
512  ltop = ktop
513  END IF
514  DO 60 krow = ltop, kwtop - 1, nv
515  kln = min( nv, kwtop-krow )
516  CALL zgemm( 'N', 'N', kln, jw, jw, one, h( krow, kwtop ),
517  $ ldh, v, ldv, zero, wv, ldwv )
518  CALL zlacpy( 'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
519  60 CONTINUE
520 *
521 * ==== Update horizontal slab in H ====
522 *
523  IF( wantt ) THEN
524  DO 70 kcol = kbot + 1, n, nh
525  kln = min( nh, n-kcol+1 )
526  CALL zgemm( 'C', 'N', jw, kln, jw, one, v, ldv,
527  $ h( kwtop, kcol ), ldh, zero, t, ldt )
528  CALL zlacpy( 'A', jw, kln, t, ldt, h( kwtop, kcol ),
529  $ ldh )
530  70 CONTINUE
531  END IF
532 *
533 * ==== Update vertical slab in Z ====
534 *
535  IF( wantz ) THEN
536  DO 80 krow = iloz, ihiz, nv
537  kln = min( nv, ihiz-krow+1 )
538  CALL zgemm( 'N', 'N', kln, jw, jw, one, z( krow, kwtop ),
539  $ ldz, v, ldv, zero, wv, ldwv )
540  CALL zlacpy( 'A', kln, jw, wv, ldwv, z( krow, kwtop ),
541  $ ldz )
542  80 CONTINUE
543  END IF
544  END IF
545 *
546 * ==== Return the number of deflations ... ====
547 *
548  nd = jw - ns
549 *
550 * ==== ... and the number of shifts. (Subtracting
551 * . INFQR from the spike length takes care
552 * . of the case of a rare QR failure while
553 * . calculating eigenvalues of the deflation
554 * . window.) ====
555 *
556  ns = ns - infqr
557 *
558 * ==== Return optimal workspace. ====
559 *
560  work( 1 ) = dcmplx( lwkopt, 0 )
561 *
562 * ==== End of ZLAQR2 ====
563 *
564  END
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:81
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
subroutine zgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
ZGEHRD
Definition: zgehrd.f:167
subroutine zlahqr(WANTT, WANTZ, N, ILO, IHI, H, LDH, W, ILOZ, IHIZ, Z, LDZ, INFO)
ZLAHQR computes the eigenvalues and Schur factorization of an upper Hessenberg matrix,...
Definition: zlahqr.f:195
subroutine zlaqr2(WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ, IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT, NV, WV, LDWV, WORK, LWORK)
ZLAQR2 performs the unitary similarity transformation of a Hessenberg matrix to detect and deflate fu...
Definition: zlaqr2.f:270
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:103
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition: zlaset.f:106
subroutine zlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition: zlarf.f:128
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:106
subroutine ztrexc(COMPQ, N, T, LDT, Q, LDQ, IFST, ILST, INFO)
ZTREXC
Definition: ztrexc.f:126
subroutine zunmhr(SIDE, TRANS, M, N, ILO, IHI, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMHR
Definition: zunmhr.f:178