LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zlaqr0()

subroutine zlaqr0 ( logical wantt,
logical wantz,
integer n,
integer ilo,
integer ihi,
complex*16, dimension( ldh, * ) h,
integer ldh,
complex*16, dimension( * ) w,
integer iloz,
integer ihiz,
complex*16, dimension( ldz, * ) z,
integer ldz,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZLAQR0 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur decomposition.

Download ZLAQR0 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!>    ZLAQR0 computes the eigenvalues of a Hessenberg matrix H
!>    and, optionally, the matrices T and Z from the Schur decomposition
!>    H = Z T Z**H, where T is an upper triangular matrix (the
!>    Schur form), and Z is the unitary matrix of Schur vectors.
!>
!>    Optionally Z may be postmultiplied into an input unitary
!>    matrix Q so that this routine can give the Schur factorization
!>    of a matrix A which has been reduced to the Hessenberg form H
!>    by the unitary matrix Q:  A = Q*H*Q**H = (QZ)*H*(QZ)**H.
!> 
Parameters
[in]WANTT
!>          WANTT is LOGICAL
!>          = .TRUE. : the full Schur form T is required;
!>          = .FALSE.: only eigenvalues are required.
!> 
[in]WANTZ
!>          WANTZ is LOGICAL
!>          = .TRUE. : the matrix of Schur vectors Z is required;
!>          = .FALSE.: Schur vectors are not required.
!> 
[in]N
!>          N is INTEGER
!>           The order of the matrix H.  N >= 0.
!> 
[in]ILO
!>          ILO is INTEGER
!> 
[in]IHI
!>          IHI is INTEGER
!>
!>           It is assumed that H is already upper triangular in rows
!>           and columns 1:ILO-1 and IHI+1:N and, if ILO > 1,
!>           H(ILO,ILO-1) is zero. ILO and IHI are normally set by a
!>           previous call to ZGEBAL, and then passed to ZGEHRD when the
!>           matrix output by ZGEBAL is reduced to Hessenberg form.
!>           Otherwise, ILO and IHI should be set to 1 and N,
!>           respectively.  If N > 0, then 1 <= ILO <= IHI <= N.
!>           If N = 0, then ILO = 1 and IHI = 0.
!> 
[in,out]H
!>          H is COMPLEX*16 array, dimension (LDH,N)
!>           On entry, the upper Hessenberg matrix H.
!>           On exit, if INFO = 0 and WANTT is .TRUE., then H
!>           contains the upper triangular matrix T from the Schur
!>           decomposition (the Schur form). If INFO = 0 and WANT is
!>           .FALSE., then the contents of H are unspecified on exit.
!>           (The output value of H when INFO > 0 is given under the
!>           description of INFO below.)
!>
!>           This subroutine may explicitly set H(i,j) = 0 for i > j and
!>           j = 1, 2, ... ILO-1 or j = IHI+1, IHI+2, ... N.
!> 
[in]LDH
!>          LDH is INTEGER
!>           The leading dimension of the array H. LDH >= max(1,N).
!> 
[out]W
!>          W is COMPLEX*16 array, dimension (N)
!>           The computed eigenvalues of H(ILO:IHI,ILO:IHI) are stored
!>           in W(ILO:IHI). If WANTT is .TRUE., then the eigenvalues are
!>           stored in the same order as on the diagonal of the Schur
!>           form returned in H, with W(i) = H(i,i).
!> 
[in]ILOZ
!>          ILOZ is INTEGER
!> 
[in]IHIZ
!>          IHIZ is INTEGER
!>           Specify the rows of Z to which transformations must be
!>           applied if WANTZ is .TRUE..
!>           1 <= ILOZ <= ILO; IHI <= IHIZ <= N.
!> 
[in,out]Z
!>          Z is COMPLEX*16 array, dimension (LDZ,IHI)
!>           If WANTZ is .FALSE., then Z is not referenced.
!>           If WANTZ is .TRUE., then Z(ILO:IHI,ILOZ:IHIZ) is
!>           replaced by Z(ILO:IHI,ILOZ:IHIZ)*U where U is the
!>           orthogonal Schur factor of H(ILO:IHI,ILO:IHI).
!>           (The output value of Z when INFO > 0 is given under
!>           the description of INFO below.)
!> 
[in]LDZ
!>          LDZ is INTEGER
!>           The leading dimension of the array Z.  if WANTZ is .TRUE.
!>           then LDZ >= MAX(1,IHIZ).  Otherwise, LDZ >= 1.
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension LWORK
!>           On exit, if LWORK = -1, WORK(1) returns an estimate of
!>           the optimal value for LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>           The dimension of the array WORK.  LWORK >= max(1,N)
!>           is sufficient, but LWORK typically as large as 6*N may
!>           be required for optimal performance.  A workspace query
!>           to determine the optimal workspace size is recommended.
!>
!>           If LWORK = -1, then ZLAQR0 does a workspace query.
!>           In this case, ZLAQR0 checks the input parameters and
!>           estimates the optimal workspace size for the given
!>           values of N, ILO and IHI.  The estimate is returned
!>           in WORK(1).  No error message related to LWORK is
!>           issued by XERBLA.  Neither H nor Z are accessed.
!> 
[out]INFO
!>          INFO is INTEGER
!>             = 0:  successful exit
!>             > 0:  if INFO = i, ZLAQR0 failed to compute all of
!>                the eigenvalues.  Elements 1:ilo-1 and i+1:n of WR
!>                and WI contain those eigenvalues which have been
!>                successfully computed.  (Failures are rare.)
!>
!>                If INFO > 0 and WANT is .FALSE., then on exit,
!>                the remaining unconverged eigenvalues are the eigen-
!>                values of the upper Hessenberg matrix rows and
!>                columns ILO through INFO of the final, output
!>                value of H.
!>
!>                If INFO > 0 and WANTT is .TRUE., then on exit
!>
!>           (*)  (initial value of H)*U  = U*(final value of H)
!>
!>                where U is a unitary matrix.  The final
!>                value of  H is upper Hessenberg and triangular in
!>                rows and columns INFO+1 through IHI.
!>
!>                If INFO > 0 and WANTZ is .TRUE., then on exit
!>
!>                  (final value of Z(ILO:IHI,ILOZ:IHIZ)
!>                   =  (initial value of Z(ILO:IHI,ILOZ:IHIZ)*U
!>
!>                where U is the unitary matrix in (*) (regard-
!>                less of the value of WANTT.)
!>
!>                If INFO > 0 and WANTZ is .FALSE., then Z is not
!>                accessed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Karen Braman and Ralph Byers, Department of Mathematics, University of Kansas, USA
References:
  K. Braman, R. Byers and R. Mathias, The Multi-Shift QR
  Algorithm Part I: Maintaining Well Focused Shifts, and Level 3
  Performance, SIAM Journal of Matrix Analysis, volume 23, pages
  929--947, 2002.

K. Braman, R. Byers and R. Mathias, The Multi-Shift QR Algorithm Part II: Aggressive Early Deflation, SIAM Journal of Matrix Analysis, volume 23, pages 948–973, 2002.

Definition at line 237 of file zlaqr0.f.

239*
240* -- LAPACK auxiliary routine --
241* -- LAPACK is a software package provided by Univ. of Tennessee, --
242* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
243*
244* .. Scalar Arguments ..
245 INTEGER IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N
246 LOGICAL WANTT, WANTZ
247* ..
248* .. Array Arguments ..
249 COMPLEX*16 H( LDH, * ), W( * ), WORK( * ), Z( LDZ, * )
250* ..
251*
252* ================================================================
253*
254* .. Parameters ..
255*
256* ==== Matrices of order NTINY or smaller must be processed by
257* . ZLAHQR because of insufficient subdiagonal scratch space.
258* . (This is a hard limit.) ====
259 INTEGER NTINY
260 parameter( ntiny = 15 )
261*
262* ==== Exceptional deflation windows: try to cure rare
263* . slow convergence by varying the size of the
264* . deflation window after KEXNW iterations. ====
265 INTEGER KEXNW
266 parameter( kexnw = 5 )
267*
268* ==== Exceptional shifts: try to cure rare slow convergence
269* . with ad-hoc exceptional shifts every KEXSH iterations.
270* . ====
271 INTEGER KEXSH
272 parameter( kexsh = 6 )
273*
274* ==== The constant WILK1 is used to form the exceptional
275* . shifts. ====
276 DOUBLE PRECISION WILK1
277 parameter( wilk1 = 0.75d0 )
278 COMPLEX*16 ZERO, ONE
279 parameter( zero = ( 0.0d0, 0.0d0 ),
280 $ one = ( 1.0d0, 0.0d0 ) )
281 DOUBLE PRECISION TWO
282 parameter( two = 2.0d0 )
283* ..
284* .. Local Scalars ..
285 COMPLEX*16 AA, BB, CC, CDUM, DD, DET, RTDISC, SWAP, TR2
286 DOUBLE PRECISION S
287 INTEGER I, INF, IT, ITMAX, K, KACC22, KBOT, KDU, KS,
288 $ KT, KTOP, KU, KV, KWH, KWTOP, KWV, LD, LS,
289 $ LWKOPT, NDEC, NDFL, NH, NHO, NIBBLE, NMIN, NS,
290 $ NSMAX, NSR, NVE, NW, NWMAX, NWR, NWUPBD
291 LOGICAL SORTED
292 CHARACTER JBCMPZ*2
293* ..
294* .. External Functions ..
295 INTEGER ILAENV
296 EXTERNAL ilaenv
297* ..
298* .. Local Arrays ..
299 COMPLEX*16 ZDUM( 1, 1 )
300* ..
301* .. External Subroutines ..
302 EXTERNAL zlacpy, zlahqr, zlaqr3, zlaqr4,
303 $ zlaqr5
304* ..
305* .. Intrinsic Functions ..
306 INTRINSIC abs, dble, dcmplx, dimag, int, max, min, mod,
307 $ sqrt
308* ..
309* .. Statement Functions ..
310 DOUBLE PRECISION CABS1
311* ..
312* .. Statement Function definitions ..
313 cabs1( cdum ) = abs( dble( cdum ) ) + abs( dimag( cdum ) )
314* ..
315* .. Executable Statements ..
316 info = 0
317*
318* ==== Quick return for N = 0: nothing to do. ====
319*
320 IF( n.EQ.0 ) THEN
321 work( 1 ) = one
322 RETURN
323 END IF
324*
325 IF( n.LE.ntiny ) THEN
326*
327* ==== Tiny matrices must use ZLAHQR. ====
328*
329 lwkopt = 1
330 IF( lwork.NE.-1 )
331 $ CALL zlahqr( wantt, wantz, n, ilo, ihi, h, ldh, w, iloz,
332 $ ihiz, z, ldz, info )
333 ELSE
334*
335* ==== Use small bulge multi-shift QR with aggressive early
336* . deflation on larger-than-tiny matrices. ====
337*
338* ==== Hope for the best. ====
339*
340 info = 0
341*
342* ==== Set up job flags for ILAENV. ====
343*
344 IF( wantt ) THEN
345 jbcmpz( 1: 1 ) = 'S'
346 ELSE
347 jbcmpz( 1: 1 ) = 'E'
348 END IF
349 IF( wantz ) THEN
350 jbcmpz( 2: 2 ) = 'V'
351 ELSE
352 jbcmpz( 2: 2 ) = 'N'
353 END IF
354*
355* ==== NWR = recommended deflation window size. At this
356* . point, N .GT. NTINY = 15, so there is enough
357* . subdiagonal workspace for NWR.GE.2 as required.
358* . (In fact, there is enough subdiagonal space for
359* . NWR.GE.4.) ====
360*
361 nwr = ilaenv( 13, 'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
362 nwr = max( 2, nwr )
363 nwr = min( ihi-ilo+1, ( n-1 ) / 3, nwr )
364*
365* ==== NSR = recommended number of simultaneous shifts.
366* . At this point N .GT. NTINY = 15, so there is at
367* . enough subdiagonal workspace for NSR to be even
368* . and greater than or equal to two as required. ====
369*
370 nsr = ilaenv( 15, 'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
371 nsr = min( nsr, ( n-3 ) / 6, ihi-ilo )
372 nsr = max( 2, nsr-mod( nsr, 2 ) )
373*
374* ==== Estimate optimal workspace ====
375*
376* ==== Workspace query call to ZLAQR3 ====
377*
378 CALL zlaqr3( wantt, wantz, n, ilo, ihi, nwr+1, h, ldh, iloz,
379 $ ihiz, z, ldz, ls, ld, w, h, ldh, n, h, ldh, n, h,
380 $ ldh, work, -1 )
381*
382* ==== Optimal workspace = MAX(ZLAQR5, ZLAQR3) ====
383*
384 lwkopt = max( 3*nsr / 2, int( work( 1 ) ) )
385*
386* ==== Quick return in case of workspace query. ====
387*
388 IF( lwork.EQ.-1 ) THEN
389 work( 1 ) = dcmplx( lwkopt, 0 )
390 RETURN
391 END IF
392*
393* ==== ZLAHQR/ZLAQR0 crossover point ====
394*
395 nmin = ilaenv( 12, 'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
396 nmin = max( ntiny, nmin )
397*
398* ==== Nibble crossover point ====
399*
400 nibble = ilaenv( 14, 'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
401 nibble = max( 0, nibble )
402*
403* ==== Accumulate reflections during ttswp? Use block
404* . 2-by-2 structure during matrix-matrix multiply? ====
405*
406 kacc22 = ilaenv( 16, 'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
407 kacc22 = max( 0, kacc22 )
408 kacc22 = min( 2, kacc22 )
409*
410* ==== NWMAX = the largest possible deflation window for
411* . which there is sufficient workspace. ====
412*
413 nwmax = min( ( n-1 ) / 3, lwork / 2 )
414 nw = nwmax
415*
416* ==== NSMAX = the Largest number of simultaneous shifts
417* . for which there is sufficient workspace. ====
418*
419 nsmax = min( ( n-3 ) / 6, 2*lwork / 3 )
420 nsmax = nsmax - mod( nsmax, 2 )
421*
422* ==== NDFL: an iteration count restarted at deflation. ====
423*
424 ndfl = 1
425*
426* ==== ITMAX = iteration limit ====
427*
428 itmax = max( 30, 2*kexsh )*max( 10, ( ihi-ilo+1 ) )
429*
430* ==== Last row and column in the active block ====
431*
432 kbot = ihi
433*
434* ==== Main Loop ====
435*
436 DO 70 it = 1, itmax
437*
438* ==== Done when KBOT falls below ILO ====
439*
440 IF( kbot.LT.ilo )
441 $ GO TO 80
442*
443* ==== Locate active block ====
444*
445 DO 10 k = kbot, ilo + 1, -1
446 IF( h( k, k-1 ).EQ.zero )
447 $ GO TO 20
448 10 CONTINUE
449 k = ilo
450 20 CONTINUE
451 ktop = k
452*
453* ==== Select deflation window size:
454* . Typical Case:
455* . If possible and advisable, nibble the entire
456* . active block. If not, use size MIN(NWR,NWMAX)
457* . or MIN(NWR+1,NWMAX) depending upon which has
458* . the smaller corresponding subdiagonal entry
459* . (a heuristic).
460* .
461* . Exceptional Case:
462* . If there have been no deflations in KEXNW or
463* . more iterations, then vary the deflation window
464* . size. At first, because, larger windows are,
465* . in general, more powerful than smaller ones,
466* . rapidly increase the window to the maximum possible.
467* . Then, gradually reduce the window size. ====
468*
469 nh = kbot - ktop + 1
470 nwupbd = min( nh, nwmax )
471 IF( ndfl.LT.kexnw ) THEN
472 nw = min( nwupbd, nwr )
473 ELSE
474 nw = min( nwupbd, 2*nw )
475 END IF
476 IF( nw.LT.nwmax ) THEN
477 IF( nw.GE.nh-1 ) THEN
478 nw = nh
479 ELSE
480 kwtop = kbot - nw + 1
481 IF( cabs1( h( kwtop, kwtop-1 ) ).GT.
482 $ cabs1( h( kwtop-1, kwtop-2 ) ) )nw = nw + 1
483 END IF
484 END IF
485 IF( ndfl.LT.kexnw ) THEN
486 ndec = -1
487 ELSE IF( ndec.GE.0 .OR. nw.GE.nwupbd ) THEN
488 ndec = ndec + 1
489 IF( nw-ndec.LT.2 )
490 $ ndec = 0
491 nw = nw - ndec
492 END IF
493*
494* ==== Aggressive early deflation:
495* . split workspace under the subdiagonal into
496* . - an nw-by-nw work array V in the lower
497* . left-hand-corner,
498* . - an NW-by-at-least-NW-but-more-is-better
499* . (NW-by-NHO) horizontal work array along
500* . the bottom edge,
501* . - an at-least-NW-but-more-is-better (NHV-by-NW)
502* . vertical work array along the left-hand-edge.
503* . ====
504*
505 kv = n - nw + 1
506 kt = nw + 1
507 nho = ( n-nw-1 ) - kt + 1
508 kwv = nw + 2
509 nve = ( n-nw ) - kwv + 1
510*
511* ==== Aggressive early deflation ====
512*
513 CALL zlaqr3( wantt, wantz, n, ktop, kbot, nw, h, ldh,
514 $ iloz,
515 $ ihiz, z, ldz, ls, ld, w, h( kv, 1 ), ldh, nho,
516 $ h( kv, kt ), ldh, nve, h( kwv, 1 ), ldh, work,
517 $ lwork )
518*
519* ==== Adjust KBOT accounting for new deflations. ====
520*
521 kbot = kbot - ld
522*
523* ==== KS points to the shifts. ====
524*
525 ks = kbot - ls + 1
526*
527* ==== Skip an expensive QR sweep if there is a (partly
528* . heuristic) reason to expect that many eigenvalues
529* . will deflate without it. Here, the QR sweep is
530* . skipped if many eigenvalues have just been deflated
531* . or if the remaining active block is small.
532*
533 IF( ( ld.EQ.0 ) .OR. ( ( 100*ld.LE.nw*nibble ) .AND. ( kbot-
534 $ ktop+1.GT.min( nmin, nwmax ) ) ) ) THEN
535*
536* ==== NS = nominal number of simultaneous shifts.
537* . This may be lowered (slightly) if ZLAQR3
538* . did not provide that many shifts. ====
539*
540 ns = min( nsmax, nsr, max( 2, kbot-ktop ) )
541 ns = ns - mod( ns, 2 )
542*
543* ==== If there have been no deflations
544* . in a multiple of KEXSH iterations,
545* . then try exceptional shifts.
546* . Otherwise use shifts provided by
547* . ZLAQR3 above or from the eigenvalues
548* . of a trailing principal submatrix. ====
549*
550 IF( mod( ndfl, kexsh ).EQ.0 ) THEN
551 ks = kbot - ns + 1
552 DO 30 i = kbot, ks + 1, -2
553 w( i ) = h( i, i ) + wilk1*cabs1( h( i, i-1 ) )
554 w( i-1 ) = w( i )
555 30 CONTINUE
556 ELSE
557*
558* ==== Got NS/2 or fewer shifts? Use ZLAQR4 or
559* . ZLAHQR on a trailing principal submatrix to
560* . get more. (Since NS.LE.NSMAX.LE.(N-3)/6,
561* . there is enough space below the subdiagonal
562* . to fit an NS-by-NS scratch array.) ====
563*
564 IF( kbot-ks+1.LE.ns / 2 ) THEN
565 ks = kbot - ns + 1
566 kt = n - ns + 1
567 CALL zlacpy( 'A', ns, ns, h( ks, ks ), ldh,
568 $ h( kt, 1 ), ldh )
569 IF( ns.GT.nmin ) THEN
570 CALL zlaqr4( .false., .false., ns, 1, ns,
571 $ h( kt, 1 ), ldh, w( ks ), 1, 1,
572 $ zdum, 1, work, lwork, inf )
573 ELSE
574 CALL zlahqr( .false., .false., ns, 1, ns,
575 $ h( kt, 1 ), ldh, w( ks ), 1, 1,
576 $ zdum, 1, inf )
577 END IF
578 ks = ks + inf
579*
580* ==== In case of a rare QR failure use
581* . eigenvalues of the trailing 2-by-2
582* . principal submatrix. Scale to avoid
583* . overflows, underflows and subnormals.
584* . (The scale factor S can not be zero,
585* . because H(KBOT,KBOT-1) is nonzero.) ====
586*
587 IF( ks.GE.kbot ) THEN
588 s = cabs1( h( kbot-1, kbot-1 ) ) +
589 $ cabs1( h( kbot, kbot-1 ) ) +
590 $ cabs1( h( kbot-1, kbot ) ) +
591 $ cabs1( h( kbot, kbot ) )
592 aa = h( kbot-1, kbot-1 ) / s
593 cc = h( kbot, kbot-1 ) / s
594 bb = h( kbot-1, kbot ) / s
595 dd = h( kbot, kbot ) / s
596 tr2 = ( aa+dd ) / two
597 det = ( aa-tr2 )*( dd-tr2 ) - bb*cc
598 rtdisc = sqrt( -det )
599 w( kbot-1 ) = ( tr2+rtdisc )*s
600 w( kbot ) = ( tr2-rtdisc )*s
601*
602 ks = kbot - 1
603 END IF
604 END IF
605*
606 IF( kbot-ks+1.GT.ns ) THEN
607*
608* ==== Sort the shifts (Helps a little) ====
609*
610 sorted = .false.
611 DO 50 k = kbot, ks + 1, -1
612 IF( sorted )
613 $ GO TO 60
614 sorted = .true.
615 DO 40 i = ks, k - 1
616 IF( cabs1( w( i ) ).LT.cabs1( w( i+1 ) ) )
617 $ THEN
618 sorted = .false.
619 swap = w( i )
620 w( i ) = w( i+1 )
621 w( i+1 ) = swap
622 END IF
623 40 CONTINUE
624 50 CONTINUE
625 60 CONTINUE
626 END IF
627 END IF
628*
629* ==== If there are only two shifts, then use
630* . only one. ====
631*
632 IF( kbot-ks+1.EQ.2 ) THEN
633 IF( cabs1( w( kbot )-h( kbot, kbot ) ).LT.
634 $ cabs1( w( kbot-1 )-h( kbot, kbot ) ) ) THEN
635 w( kbot-1 ) = w( kbot )
636 ELSE
637 w( kbot ) = w( kbot-1 )
638 END IF
639 END IF
640*
641* ==== Use up to NS of the the smallest magnitude
642* . shifts. If there aren't NS shifts available,
643* . then use them all, possibly dropping one to
644* . make the number of shifts even. ====
645*
646 ns = min( ns, kbot-ks+1 )
647 ns = ns - mod( ns, 2 )
648 ks = kbot - ns + 1
649*
650* ==== Small-bulge multi-shift QR sweep:
651* . split workspace under the subdiagonal into
652* . - a KDU-by-KDU work array U in the lower
653* . left-hand-corner,
654* . - a KDU-by-at-least-KDU-but-more-is-better
655* . (KDU-by-NHo) horizontal work array WH along
656* . the bottom edge,
657* . - and an at-least-KDU-but-more-is-better-by-KDU
658* . (NVE-by-KDU) vertical work WV arrow along
659* . the left-hand-edge. ====
660*
661 kdu = 2*ns
662 ku = n - kdu + 1
663 kwh = kdu + 1
664 nho = ( n-kdu+1-4 ) - ( kdu+1 ) + 1
665 kwv = kdu + 4
666 nve = n - kdu - kwv + 1
667*
668* ==== Small-bulge multi-shift QR sweep ====
669*
670 CALL zlaqr5( wantt, wantz, kacc22, n, ktop, kbot, ns,
671 $ w( ks ), h, ldh, iloz, ihiz, z, ldz, work,
672 $ 3, h( ku, 1 ), ldh, nve, h( kwv, 1 ), ldh,
673 $ nho, h( ku, kwh ), ldh )
674 END IF
675*
676* ==== Note progress (or the lack of it). ====
677*
678 IF( ld.GT.0 ) THEN
679 ndfl = 1
680 ELSE
681 ndfl = ndfl + 1
682 END IF
683*
684* ==== End of main loop ====
685 70 CONTINUE
686*
687* ==== Iteration limit exceeded. Set INFO to show where
688* . the problem occurred and exit. ====
689*
690 info = kbot
691 80 CONTINUE
692 END IF
693*
694* ==== Return the optimal value of LWORK. ====
695*
696 work( 1 ) = dcmplx( lwkopt, 0 )
697*
698* ==== End of ZLAQR0 ====
699*
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
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:193
subroutine zlaqr3(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)
ZLAQR3 performs the unitary similarity transformation of a Hessenberg matrix to detect and deflate fu...
Definition zlaqr3.f:266
subroutine zlaqr4(wantt, wantz, n, ilo, ihi, h, ldh, w, iloz, ihiz, z, ldz, work, lwork, info)
ZLAQR4 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur de...
Definition zlaqr4.f:245
subroutine zlaqr5(wantt, wantz, kacc22, n, ktop, kbot, nshfts, s, h, ldh, iloz, ihiz, z, ldz, v, ldv, u, ldu, nv, wv, ldwv, nh, wh, ldwh)
ZLAQR5 performs a single small-bulge multi-shift QR sweep.
Definition zlaqr5.f:256
Here is the call graph for this function:
Here is the caller graph for this function: