239 SUBROUTINE zlaqr0( WANTT, WANTZ, N, ILO, IHI, H, LDH, W, ILOZ,
240 $ IHIZ, Z, LDZ, WORK, LWORK, INFO )
247 INTEGER IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N
251 COMPLEX*16 H( LDH, * ), W( * ), WORK( * ), Z( LDZ, * )
262 parameter( ntiny = 15 )
268 parameter( kexnw = 5 )
274 parameter( kexsh = 6 )
278 DOUBLE PRECISION WILK1
279 parameter( wilk1 = 0.75d0 )
281 parameter( zero = ( 0.0d0, 0.0d0 ),
282 $ one = ( 1.0d0, 0.0d0 ) )
284 parameter( two = 2.0d0 )
287 COMPLEX*16 AA, BB, CC, CDUM, DD, DET, RTDISC, SWAP, TR2
289 INTEGER I, INF, IT, ITMAX, K, KACC22, KBOT, KDU, KS,
290 $ kt, ktop, ku, kv, kwh, kwtop, kwv, ld, ls,
291 $ lwkopt, ndec, ndfl, nh, nho, nibble, nmin, ns,
292 $ nsmax, nsr, nve, nw, nwmax, nwr, nwupbd
301 COMPLEX*16 ZDUM( 1, 1 )
307 INTRINSIC abs, dble, dcmplx, dimag, int, max, min, mod,
311 DOUBLE PRECISION CABS1
314 cabs1( cdum ) = abs( dble( cdum ) ) + abs( dimag( cdum ) )
326 IF( n.LE.ntiny )
THEN
332 $
CALL zlahqr( wantt, wantz, n, ilo, ihi, h, ldh, w, iloz,
333 $ ihiz, z, ldz, info )
362 nwr = ilaenv( 13,
'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
364 nwr = min( ihi-ilo+1, ( n-1 ) / 3, nwr )
371 nsr = ilaenv( 15,
'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
372 nsr = min( nsr, ( n-3 ) / 6, ihi-ilo )
373 nsr = max( 2, nsr-mod( nsr, 2 ) )
379 CALL zlaqr3( wantt, wantz, n, ilo, ihi, nwr+1, h, ldh, iloz,
380 $ ihiz, z, ldz, ls, ld, w, h, ldh, n, h, ldh, n, h,
385 lwkopt = max( 3*nsr / 2, int( work( 1 ) ) )
389 IF( lwork.EQ.-1 )
THEN
390 work( 1 ) = dcmplx( lwkopt, 0 )
396 nmin = ilaenv( 12,
'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
397 nmin = max( ntiny, nmin )
401 nibble = ilaenv( 14,
'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
402 nibble = max( 0, nibble )
407 kacc22 = ilaenv( 16,
'ZLAQR0', jbcmpz, n, ilo, ihi, lwork )
408 kacc22 = max( 0, kacc22 )
409 kacc22 = min( 2, kacc22 )
414 nwmax = min( ( n-1 ) / 3, lwork / 2 )
420 nsmax = min( ( n-3 ) / 6, 2*lwork / 3 )
421 nsmax = nsmax - mod( nsmax, 2 )
429 itmax = max( 30, 2*kexsh )*max( 10, ( ihi-ilo+1 ) )
446 DO 10 k = kbot, ilo + 1, -1
447 IF( h( k, k-1 ).EQ.zero )
471 nwupbd = min( nh, nwmax )
472 IF( ndfl.LT.kexnw )
THEN
473 nw = min( nwupbd, nwr )
475 nw = min( nwupbd, 2*nw )
477 IF( nw.LT.nwmax )
THEN
478 IF( nw.GE.nh-1 )
THEN
481 kwtop = kbot - nw + 1
482 IF( cabs1( h( kwtop, kwtop-1 ) ).GT.
483 $ cabs1( h( kwtop-1, kwtop-2 ) ) )nw = nw + 1
486 IF( ndfl.LT.kexnw )
THEN
488 ELSE IF( ndec.GE.0 .OR. nw.GE.nwupbd )
THEN
508 nho = ( n-nw-1 ) - kt + 1
510 nve = ( n-nw ) - kwv + 1
514 CALL zlaqr3( wantt, wantz, n, ktop, kbot, nw, h, ldh, iloz,
515 $ ihiz, z, ldz, ls, ld, w, h( kv, 1 ), ldh, nho,
516 $ h( kv, kt ), ldh, nve, h( kwv, 1 ), ldh, work,
533 IF( ( ld.EQ.0 ) .OR. ( ( 100*ld.LE.nw*nibble ) .AND. ( kbot-
534 $ ktop+1.GT.min( nmin, nwmax ) ) ) )
THEN
540 ns = min( nsmax, nsr, max( 2, kbot-ktop ) )
541 ns = ns - mod( ns, 2 )
550 IF( mod( ndfl, kexsh ).EQ.0 )
THEN
552 DO 30 i = kbot, ks + 1, -2
553 w( i ) = h( i, i ) + wilk1*cabs1( h( i, i-1 ) )
564 IF( kbot-ks+1.LE.ns / 2 )
THEN
567 CALL zlacpy(
'A', ns, ns, h( ks, ks ), 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 )
574 CALL zlahqr( .false., .false., ns, 1, ns,
575 $ h( kt, 1 ), ldh, w( ks ), 1, 1,
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
606 IF( kbot-ks+1.GT.ns )
THEN
611 DO 50 k = kbot, ks + 1, -1
616 IF( cabs1( w( i ) ).LT.cabs1( w( i+1 ) ) )
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 )
637 w( kbot ) = w( kbot-1 )
646 ns = min( ns, kbot-ks+1 )
647 ns = ns - mod( ns, 2 )
664 nho = ( n-kdu+1-4 ) - ( kdu+1 ) + 1
666 nve = n - kdu - kwv + 1
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 )
696 work( 1 ) = dcmplx( lwkopt, 0 )
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
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,...
subroutine zlaqr0(wantt, wantz, n, ilo, ihi, h, ldh, w, iloz, ihiz, z, ldz, work, lwork, info)
ZLAQR0 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur de...
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...
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...
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.