LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dlaqr3.f
Go to the documentation of this file.
1*> \brief \b DLAQR3 performs the orthogonal 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 DLAQR3 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaqr3.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaqr3.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaqr3.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DLAQR3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
22* IHIZ, Z, LDZ, NS, ND, SR, SI, V, LDV, NH, T,
23* LDT, 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* DOUBLE PRECISION H( LDH, * ), SI( * ), SR( * ), T( LDT, * ),
32* $ V( LDV, * ), WORK( * ), WV( LDWV, * ),
33* $ Z( LDZ, * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> Aggressive early deflation:
43*>
44*> DLAQR3 accepts as input an upper Hessenberg matrix
45*> H and performs an orthogonal similarity transformation
46*> designed to detect and deflate fully converged eigenvalues from
47*> a trailing principal submatrix. On output H has been over-
48*> written by a new Hessenberg matrix that is a perturbation of
49*> an orthogonal similarity transformation of H. It is to be
50*> hoped that the final version of H has many zero subdiagonal
51*> entries.
52*> \endverbatim
53*
54* Arguments:
55* ==========
56*
57*> \param[in] WANTT
58*> \verbatim
59*> WANTT is LOGICAL
60*> If .TRUE., then the Hessenberg matrix H is fully updated
61*> so that the quasi-triangular Schur factor may be
62*> computed (in cooperation with the calling subroutine).
63*> If .FALSE., then only enough of H is updated to preserve
64*> the eigenvalues.
65*> \endverbatim
66*>
67*> \param[in] WANTZ
68*> \verbatim
69*> WANTZ is LOGICAL
70*> If .TRUE., then the orthogonal matrix Z is updated so
71*> so that the orthogonal Schur factor may be computed
72*> (in cooperation with the calling subroutine).
73*> If .FALSE., then Z is not referenced.
74*> \endverbatim
75*>
76*> \param[in] N
77*> \verbatim
78*> N is INTEGER
79*> The order of the matrix H and (if WANTZ is .TRUE.) the
80*> order of the orthogonal matrix Z.
81*> \endverbatim
82*>
83*> \param[in] KTOP
84*> \verbatim
85*> KTOP is INTEGER
86*> It is assumed that either KTOP = 1 or H(KTOP,KTOP-1)=0.
87*> KBOT and KTOP together determine an isolated block
88*> along the diagonal of the Hessenberg matrix.
89*> \endverbatim
90*>
91*> \param[in] KBOT
92*> \verbatim
93*> KBOT is INTEGER
94*> It is assumed without a check that either
95*> KBOT = N or H(KBOT+1,KBOT)=0. KBOT and KTOP together
96*> determine an isolated block along the diagonal of the
97*> Hessenberg matrix.
98*> \endverbatim
99*>
100*> \param[in] NW
101*> \verbatim
102*> NW is INTEGER
103*> Deflation window size. 1 <= NW <= (KBOT-KTOP+1).
104*> \endverbatim
105*>
106*> \param[in,out] H
107*> \verbatim
108*> H is DOUBLE PRECISION array, dimension (LDH,N)
109*> On input the initial N-by-N section of H stores the
110*> Hessenberg matrix undergoing aggressive early deflation.
111*> On output H has been transformed by an orthogonal
112*> similarity transformation, perturbed, and the returned
113*> to Hessenberg form that (it is to be hoped) has some
114*> zero subdiagonal entries.
115*> \endverbatim
116*>
117*> \param[in] LDH
118*> \verbatim
119*> LDH is INTEGER
120*> Leading dimension of H just as declared in the calling
121*> subroutine. N <= LDH
122*> \endverbatim
123*>
124*> \param[in] ILOZ
125*> \verbatim
126*> ILOZ is INTEGER
127*> \endverbatim
128*>
129*> \param[in] IHIZ
130*> \verbatim
131*> IHIZ is INTEGER
132*> Specify the rows of Z to which transformations must be
133*> applied if WANTZ is .TRUE.. 1 <= ILOZ <= IHIZ <= N.
134*> \endverbatim
135*>
136*> \param[in,out] Z
137*> \verbatim
138*> Z is DOUBLE PRECISION array, dimension (LDZ,N)
139*> IF WANTZ is .TRUE., then on output, the orthogonal
140*> similarity transformation mentioned above has been
141*> accumulated into Z(ILOZ:IHIZ,ILOZ:IHIZ) from the right.
142*> If WANTZ is .FALSE., then Z is unreferenced.
143*> \endverbatim
144*>
145*> \param[in] LDZ
146*> \verbatim
147*> LDZ is INTEGER
148*> The leading dimension of Z just as declared in the
149*> calling subroutine. 1 <= LDZ.
150*> \endverbatim
151*>
152*> \param[out] NS
153*> \verbatim
154*> NS is INTEGER
155*> The number of unconverged (ie approximate) eigenvalues
156*> returned in SR and SI that may be used as shifts by the
157*> calling subroutine.
158*> \endverbatim
159*>
160*> \param[out] ND
161*> \verbatim
162*> ND is INTEGER
163*> The number of converged eigenvalues uncovered by this
164*> subroutine.
165*> \endverbatim
166*>
167*> \param[out] SR
168*> \verbatim
169*> SR is DOUBLE PRECISION array, dimension (KBOT)
170*> \endverbatim
171*>
172*> \param[out] SI
173*> \verbatim
174*> SI is DOUBLE PRECISION array, dimension (KBOT)
175*> On output, the real and imaginary parts of approximate
176*> eigenvalues that may be used for shifts are stored in
177*> SR(KBOT-ND-NS+1) through SR(KBOT-ND) and
178*> SI(KBOT-ND-NS+1) through SI(KBOT-ND), respectively.
179*> The real and imaginary parts of converged eigenvalues
180*> are stored in SR(KBOT-ND+1) through SR(KBOT) and
181*> SI(KBOT-ND+1) through SI(KBOT), respectively.
182*> \endverbatim
183*>
184*> \param[out] V
185*> \verbatim
186*> V is DOUBLE PRECISION array, dimension (LDV,NW)
187*> An NW-by-NW work array.
188*> \endverbatim
189*>
190*> \param[in] LDV
191*> \verbatim
192*> LDV is INTEGER
193*> The leading dimension of V just as declared in the
194*> calling subroutine. NW <= LDV
195*> \endverbatim
196*>
197*> \param[in] NH
198*> \verbatim
199*> NH is INTEGER
200*> The number of columns of T. NH >= NW.
201*> \endverbatim
202*>
203*> \param[out] T
204*> \verbatim
205*> T is DOUBLE PRECISION array, dimension (LDT,NW)
206*> \endverbatim
207*>
208*> \param[in] LDT
209*> \verbatim
210*> LDT is INTEGER
211*> The leading dimension of T just as declared in the
212*> calling subroutine. NW <= LDT
213*> \endverbatim
214*>
215*> \param[in] NV
216*> \verbatim
217*> NV is INTEGER
218*> The number of rows of work array WV available for
219*> workspace. NV >= NW.
220*> \endverbatim
221*>
222*> \param[out] WV
223*> \verbatim
224*> WV is DOUBLE PRECISION array, dimension (LDWV,NW)
225*> \endverbatim
226*>
227*> \param[in] LDWV
228*> \verbatim
229*> LDWV is INTEGER
230*> The leading dimension of W just as declared in the
231*> calling subroutine. NW <= LDV
232*> \endverbatim
233*>
234*> \param[out] WORK
235*> \verbatim
236*> WORK is DOUBLE PRECISION array, dimension (LWORK)
237*> On exit, WORK(1) is set to an estimate of the optimal value
238*> of LWORK for the given values of N, NW, KTOP and KBOT.
239*> \endverbatim
240*>
241*> \param[in] LWORK
242*> \verbatim
243*> LWORK is INTEGER
244*> The dimension of the work array WORK. LWORK = 2*NW
245*> suffices, but greater efficiency may result from larger
246*> values of LWORK.
247*>
248*> If LWORK = -1, then a workspace query is assumed; DLAQR3
249*> only estimates the optimal workspace size for the given
250*> values of N, NW, KTOP and KBOT. The estimate is returned
251*> in WORK(1). No error message related to LWORK is issued
252*> by XERBLA. Neither H nor Z are accessed.
253*> \endverbatim
254*
255* Authors:
256* ========
257*
258*> \author Univ. of Tennessee
259*> \author Univ. of California Berkeley
260*> \author Univ. of Colorado Denver
261*> \author NAG Ltd.
262*
263*> \ingroup laqr3
264*
265*> \par Contributors:
266* ==================
267*>
268*> Karen Braman and Ralph Byers, Department of Mathematics,
269*> University of Kansas, USA
270*>
271* =====================================================================
272 SUBROUTINE dlaqr3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
273 $ IHIZ, Z, LDZ, NS, ND, SR, SI, V, LDV, NH, T,
274 $ LDT, NV, WV, LDWV, WORK, LWORK )
275*
276* -- LAPACK auxiliary routine --
277* -- LAPACK is a software package provided by Univ. of Tennessee, --
278* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
279*
280* .. Scalar Arguments ..
281 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
282 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
283 LOGICAL WANTT, WANTZ
284* ..
285* .. Array Arguments ..
286 DOUBLE PRECISION H( LDH, * ), SI( * ), SR( * ), T( LDT, * ),
287 $ V( LDV, * ), WORK( * ), WV( LDWV, * ),
288 $ z( ldz, * )
289* ..
290*
291* ================================================================
292* .. Parameters ..
293 DOUBLE PRECISION ZERO, ONE
294 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
295* ..
296* .. Local Scalars ..
297 DOUBLE PRECISION AA, BB, BETA, CC, CS, DD, EVI, EVK, FOO, S,
298 $ SAFMAX, SAFMIN, SMLNUM, SN, TAU, ULP
299 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, K, KCOL,
300 $ KEND, KLN, KROW, KWTOP, LTOP, LWK1, LWK2, LWK3,
301 $ lwkopt, nmin
302 LOGICAL BULGE, SORTED
303* ..
304* .. External Functions ..
305 DOUBLE PRECISION DLAMCH
306 INTEGER ILAENV
307 EXTERNAL dlamch, ilaenv
308* ..
309* .. External Subroutines ..
310 EXTERNAL dcopy, dgehrd, dgemm, dlacpy, dlahqr, dlanv2,
312* ..
313* .. Intrinsic Functions ..
314 INTRINSIC abs, dble, int, max, min, sqrt
315* ..
316* .. Executable Statements ..
317*
318* ==== Estimate optimal workspace. ====
319*
320 jw = min( nw, kbot-ktop+1 )
321 IF( jw.LE.2 ) THEN
322 lwkopt = 1
323 ELSE
324*
325* ==== Workspace query call to DGEHRD ====
326*
327 CALL dgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
328 lwk1 = int( work( 1 ) )
329*
330* ==== Workspace query call to DORMHR ====
331*
332 CALL dormhr( 'R', 'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
333 $ work, -1, info )
334 lwk2 = int( work( 1 ) )
335*
336* ==== Workspace query call to DLAQR4 ====
337*
338 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr, si, 1, jw,
339 $ v, ldv, work, -1, infqr )
340 lwk3 = int( work( 1 ) )
341*
342* ==== Optimal workspace ====
343*
344 lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
345 END IF
346*
347* ==== Quick return in case of workspace query. ====
348*
349 IF( lwork.EQ.-1 ) THEN
350 work( 1 ) = dble( lwkopt )
351 RETURN
352 END IF
353*
354* ==== Nothing to do ...
355* ... for an empty active block ... ====
356 ns = 0
357 nd = 0
358 work( 1 ) = one
359 IF( ktop.GT.kbot )
360 $ RETURN
361* ... nor for an empty deflation window. ====
362 IF( nw.LT.1 )
363 $ RETURN
364*
365* ==== Machine constants ====
366*
367 safmin = dlamch( 'SAFE MINIMUM' )
368 safmax = one / safmin
369 ulp = dlamch( 'PRECISION' )
370 smlnum = safmin*( dble( n ) / ulp )
371*
372* ==== Setup deflation window ====
373*
374 jw = min( nw, kbot-ktop+1 )
375 kwtop = kbot - jw + 1
376 IF( kwtop.EQ.ktop ) THEN
377 s = zero
378 ELSE
379 s = h( kwtop, kwtop-1 )
380 END IF
381*
382 IF( kbot.EQ.kwtop ) THEN
383*
384* ==== 1-by-1 deflation window: not much to do ====
385*
386 sr( kwtop ) = h( kwtop, kwtop )
387 si( kwtop ) = zero
388 ns = 1
389 nd = 0
390 IF( abs( s ).LE.max( smlnum, ulp*abs( h( kwtop, kwtop ) ) ) )
391 $ THEN
392 ns = 0
393 nd = 1
394 IF( kwtop.GT.ktop )
395 $ h( kwtop, kwtop-1 ) = zero
396 END IF
397 work( 1 ) = one
398 RETURN
399 END IF
400*
401* ==== Convert to spike-triangular form. (In case of a
402* . rare QR failure, this routine continues to do
403* . aggressive early deflation using that part of
404* . the deflation window that converged using INFQR
405* . here and there to keep track.) ====
406*
407 CALL dlacpy( 'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
408 CALL dcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
409*
410 CALL dlaset( 'A', jw, jw, zero, one, v, ldv )
411 nmin = ilaenv( 12, 'DLAQR3', 'SV', jw, 1, jw, lwork )
412 IF( jw.GT.nmin ) THEN
413 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
414 $ si( kwtop ), 1, jw, v, ldv, work, lwork, infqr )
415 ELSE
416 CALL dlahqr( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
417 $ si( kwtop ), 1, jw, v, ldv, infqr )
418 END IF
419*
420* ==== DTREXC needs a clean margin near the diagonal ====
421*
422 DO 10 j = 1, jw - 3
423 t( j+2, j ) = zero
424 t( j+3, j ) = zero
425 10 CONTINUE
426 IF( jw.GT.2 )
427 $ t( jw, jw-2 ) = zero
428*
429* ==== Deflation detection loop ====
430*
431 ns = jw
432 ilst = infqr + 1
433 20 CONTINUE
434 IF( ilst.LE.ns ) THEN
435 IF( ns.EQ.1 ) THEN
436 bulge = .false.
437 ELSE
438 bulge = t( ns, ns-1 ).NE.zero
439 END IF
440*
441* ==== Small spike tip test for deflation ====
442*
443 IF( .NOT. bulge ) THEN
444*
445* ==== Real eigenvalue ====
446*
447 foo = abs( t( ns, ns ) )
448 IF( foo.EQ.zero )
449 $ foo = abs( s )
450 IF( abs( s*v( 1, ns ) ).LE.max( smlnum, ulp*foo ) ) THEN
451*
452* ==== Deflatable ====
453*
454 ns = ns - 1
455 ELSE
456*
457* ==== Undeflatable. Move it up out of the way.
458* . (DTREXC can not fail in this case.) ====
459*
460 ifst = ns
461 CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
462 $ info )
463 ilst = ilst + 1
464 END IF
465 ELSE
466*
467* ==== Complex conjugate pair ====
468*
469 foo = abs( t( ns, ns ) ) + sqrt( abs( t( ns, ns-1 ) ) )*
470 $ sqrt( abs( t( ns-1, ns ) ) )
471 IF( foo.EQ.zero )
472 $ foo = abs( s )
473 IF( max( abs( s*v( 1, ns ) ), abs( s*v( 1, ns-1 ) ) ).LE.
474 $ max( smlnum, ulp*foo ) ) THEN
475*
476* ==== Deflatable ====
477*
478 ns = ns - 2
479 ELSE
480*
481* ==== Undeflatable. Move them up out of the way.
482* . Fortunately, DTREXC does the right thing with
483* . ILST in case of a rare exchange failure. ====
484*
485 ifst = ns
486 CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
487 $ info )
488 ilst = ilst + 2
489 END IF
490 END IF
491*
492* ==== End deflation detection loop ====
493*
494 GO TO 20
495 END IF
496*
497* ==== Return to Hessenberg form ====
498*
499 IF( ns.EQ.0 )
500 $ s = zero
501*
502 IF( ns.LT.jw ) THEN
503*
504* ==== sorting diagonal blocks of T improves accuracy for
505* . graded matrices. Bubble sort deals well with
506* . exchange failures. ====
507*
508 sorted = .false.
509 i = ns + 1
510 30 CONTINUE
511 IF( sorted )
512 $ GO TO 50
513 sorted = .true.
514*
515 kend = i - 1
516 i = infqr + 1
517 IF( i.EQ.ns ) THEN
518 k = i + 1
519 ELSE IF( t( i+1, i ).EQ.zero ) THEN
520 k = i + 1
521 ELSE
522 k = i + 2
523 END IF
524 40 CONTINUE
525 IF( k.LE.kend ) THEN
526 IF( k.EQ.i+1 ) THEN
527 evi = abs( t( i, i ) )
528 ELSE
529 evi = abs( t( i, i ) ) + sqrt( abs( t( i+1, i ) ) )*
530 $ sqrt( abs( t( i, i+1 ) ) )
531 END IF
532*
533 IF( k.EQ.kend ) THEN
534 evk = abs( t( k, k ) )
535 ELSE IF( t( k+1, k ).EQ.zero ) THEN
536 evk = abs( t( k, k ) )
537 ELSE
538 evk = abs( t( k, k ) ) + sqrt( abs( t( k+1, k ) ) )*
539 $ sqrt( abs( t( k, k+1 ) ) )
540 END IF
541*
542 IF( evi.GE.evk ) THEN
543 i = k
544 ELSE
545 sorted = .false.
546 ifst = i
547 ilst = k
548 CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
549 $ info )
550 IF( info.EQ.0 ) THEN
551 i = ilst
552 ELSE
553 i = k
554 END IF
555 END IF
556 IF( i.EQ.kend ) THEN
557 k = i + 1
558 ELSE IF( t( i+1, i ).EQ.zero ) THEN
559 k = i + 1
560 ELSE
561 k = i + 2
562 END IF
563 GO TO 40
564 END IF
565 GO TO 30
566 50 CONTINUE
567 END IF
568*
569* ==== Restore shift/eigenvalue array from T ====
570*
571 i = jw
572 60 CONTINUE
573 IF( i.GE.infqr+1 ) THEN
574 IF( i.EQ.infqr+1 ) THEN
575 sr( kwtop+i-1 ) = t( i, i )
576 si( kwtop+i-1 ) = zero
577 i = i - 1
578 ELSE IF( t( i, i-1 ).EQ.zero ) THEN
579 sr( kwtop+i-1 ) = t( i, i )
580 si( kwtop+i-1 ) = zero
581 i = i - 1
582 ELSE
583 aa = t( i-1, i-1 )
584 cc = t( i, i-1 )
585 bb = t( i-1, i )
586 dd = t( i, i )
587 CALL dlanv2( aa, bb, cc, dd, sr( kwtop+i-2 ),
588 $ si( kwtop+i-2 ), sr( kwtop+i-1 ),
589 $ si( kwtop+i-1 ), cs, sn )
590 i = i - 2
591 END IF
592 GO TO 60
593 END IF
594*
595 IF( ns.LT.jw .OR. s.EQ.zero ) THEN
596 IF( ns.GT.1 .AND. s.NE.zero ) THEN
597*
598* ==== Reflect spike back into lower triangle ====
599*
600 CALL dcopy( ns, v, ldv, work, 1 )
601 beta = work( 1 )
602 CALL dlarfg( ns, beta, work( 2 ), 1, tau )
603 work( 1 ) = one
604*
605 CALL dlaset( 'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
606*
607 CALL dlarf( 'L', ns, jw, work, 1, tau, t, ldt,
608 $ work( jw+1 ) )
609 CALL dlarf( 'R', ns, ns, work, 1, tau, t, ldt,
610 $ work( jw+1 ) )
611 CALL dlarf( 'R', jw, ns, work, 1, tau, v, ldv,
612 $ work( jw+1 ) )
613*
614 CALL dgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
615 $ lwork-jw, info )
616 END IF
617*
618* ==== Copy updated reduced window into place ====
619*
620 IF( kwtop.GT.1 )
621 $ h( kwtop, kwtop-1 ) = s*v( 1, 1 )
622 CALL dlacpy( 'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
623 CALL dcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
624 $ ldh+1 )
625*
626* ==== Accumulate orthogonal matrix in order update
627* . H and Z, if requested. ====
628*
629 IF( ns.GT.1 .AND. s.NE.zero )
630 $ CALL dormhr( 'R', 'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
631 $ work( jw+1 ), lwork-jw, info )
632*
633* ==== Update vertical slab in H ====
634*
635 IF( wantt ) THEN
636 ltop = 1
637 ELSE
638 ltop = ktop
639 END IF
640 DO 70 krow = ltop, kwtop - 1, nv
641 kln = min( nv, kwtop-krow )
642 CALL dgemm( 'N', 'N', kln, jw, jw, one, h( krow, kwtop ),
643 $ ldh, v, ldv, zero, wv, ldwv )
644 CALL dlacpy( 'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
645 70 CONTINUE
646*
647* ==== Update horizontal slab in H ====
648*
649 IF( wantt ) THEN
650 DO 80 kcol = kbot + 1, n, nh
651 kln = min( nh, n-kcol+1 )
652 CALL dgemm( 'C', 'N', jw, kln, jw, one, v, ldv,
653 $ h( kwtop, kcol ), ldh, zero, t, ldt )
654 CALL dlacpy( 'A', jw, kln, t, ldt, h( kwtop, kcol ),
655 $ ldh )
656 80 CONTINUE
657 END IF
658*
659* ==== Update vertical slab in Z ====
660*
661 IF( wantz ) THEN
662 DO 90 krow = iloz, ihiz, nv
663 kln = min( nv, ihiz-krow+1 )
664 CALL dgemm( 'N', 'N', kln, jw, jw, one, z( krow, kwtop ),
665 $ ldz, v, ldv, zero, wv, ldwv )
666 CALL dlacpy( 'A', kln, jw, wv, ldwv, z( krow, kwtop ),
667 $ ldz )
668 90 CONTINUE
669 END IF
670 END IF
671*
672* ==== Return the number of deflations ... ====
673*
674 nd = jw - ns
675*
676* ==== ... and the number of shifts. (Subtracting
677* . INFQR from the spike length takes care
678* . of the case of a rare QR failure while
679* . calculating eigenvalues of the deflation
680* . window.) ====
681*
682 ns = ns - infqr
683*
684* ==== Return optimal workspace. ====
685*
686 work( 1 ) = dble( lwkopt )
687*
688* ==== End of DLAQR3 ====
689*
690 END
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
DGEHRD
Definition dgehrd.f:167
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
subroutine dlahqr(wantt, wantz, n, ilo, ihi, h, ldh, wr, wi, iloz, ihiz, z, ldz, info)
DLAHQR computes the eigenvalues and Schur factorization of an upper Hessenberg matrix,...
Definition dlahqr.f:207
subroutine dlanv2(a, b, c, d, rt1r, rt1i, rt2r, rt2i, cs, sn)
DLANV2 computes the Schur factorization of a real 2-by-2 nonsymmetric matrix in standard form.
Definition dlanv2.f:127
subroutine dlaqr3(wantt, wantz, n, ktop, kbot, nw, h, ldh, iloz, ihiz, z, ldz, ns, nd, sr, si, v, ldv, nh, t, ldt, nv, wv, ldwv, work, lwork)
DLAQR3 performs the orthogonal similarity transformation of a Hessenberg matrix to detect and deflate...
Definition dlaqr3.f:275
subroutine dlaqr4(wantt, wantz, n, ilo, ihi, h, ldh, wr, wi, iloz, ihiz, z, ldz, work, lwork, info)
DLAQR4 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur de...
Definition dlaqr4.f:263
subroutine dlarf(side, m, n, v, incv, tau, c, ldc, work)
DLARF applies an elementary reflector to a general rectangular matrix.
Definition dlarf.f:124
subroutine dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
Definition dlarfg.f:106
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:110
subroutine dtrexc(compq, n, t, ldt, q, ldq, ifst, ilst, work, info)
DTREXC
Definition dtrexc.f:148
subroutine dormhr(side, trans, m, n, ilo, ihi, a, lda, tau, c, ldc, work, lwork, info)
DORMHR
Definition dormhr.f:178