264 SUBROUTINE zlaqr3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
265 $ IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT,
266 $ NV, WV, LDWV, WORK, LWORK )
273 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
274 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
278 COMPLEX*16 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
279 $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
286 PARAMETER ( ZERO = ( 0.0d0, 0.0d0 ),
287 $ one = ( 1.0d0, 0.0d0 ) )
288 DOUBLE PRECISION RZERO, RONE
289 PARAMETER ( RZERO = 0.0d0, rone = 1.0d0 )
292 COMPLEX*16 BETA, CDUM, S, TAU
293 DOUBLE PRECISION FOO, SAFMAX, SAFMIN, SMLNUM, ULP
294 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
295 $ knt, krow, kwtop, ltop, lwk1, lwk2, lwk3,
299 DOUBLE PRECISION DLAMCH
301 EXTERNAL dlamch, ilaenv
308 INTRINSIC abs, dble, dcmplx, dconjg, dimag, int, max, min
311 DOUBLE PRECISION CABS1
314 cabs1( cdum ) = abs( dble( cdum ) ) + abs( dimag( cdum ) )
320 jw = min( nw, kbot-ktop+1 )
327 CALL zgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
328 lwk1 = int( work( 1 ) )
332 CALL zunmhr(
'R',
'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
334 lwk2 = int( work( 1 ) )
338 CALL zlaqr4( .true., .true., jw, 1, jw, t, ldt, sh, 1, jw, v,
339 $ ldv, work, -1, infqr )
340 lwk3 = int( work( 1 ) )
344 lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
349 IF( lwork.EQ.-1 )
THEN
350 work( 1 ) = dcmplx( lwkopt, 0 )
367 safmin = dlamch(
'SAFE MINIMUM' )
368 safmax = rone / safmin
369 ulp = dlamch(
'PRECISION' )
370 smlnum = safmin*( dble( n ) / ulp )
374 jw = min( nw, kbot-ktop+1 )
375 kwtop = kbot - jw + 1
376 IF( kwtop.EQ.ktop )
THEN
379 s = h( kwtop, kwtop-1 )
382 IF( kbot.EQ.kwtop )
THEN
386 sh( kwtop ) = h( kwtop, kwtop )
389 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
394 $ h( kwtop, kwtop-1 ) = zero
406 CALL zlacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
407 CALL zcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
409 CALL zlaset(
'A', jw, jw, zero, one, v, ldv )
410 nmin = ilaenv( 12,
'ZLAQR3',
'SV', jw, 1, jw, lwork )
411 IF( jw.GT.nmin )
THEN
412 CALL zlaqr4( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
413 $ jw, v, ldv, work, lwork, infqr )
415 CALL zlahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
416 $ jw, v, ldv, infqr )
423 DO 10 knt = infqr + 1, jw
427 foo = cabs1( t( ns, ns ) )
430 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
442 CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
457 DO 30 i = infqr + 1, ns
460 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
465 $
CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
471 DO 40 i = infqr + 1, jw
472 sh( kwtop+i-1 ) = t( i, i )
476 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
477 IF( ns.GT.1 .AND. s.NE.zero )
THEN
481 CALL zcopy( ns, v, ldv, work, 1 )
483 work( i ) = dconjg( work( i ) )
486 CALL zlarfg( ns, beta, work( 2 ), 1, tau )
489 CALL zlaset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
491 CALL zlarf(
'L', ns, jw, work, 1, dconjg( tau ), t, ldt,
493 CALL zlarf(
'R', ns, ns, work, 1, tau, t, ldt,
495 CALL zlarf(
'R', jw, ns, work, 1, tau, v, ldv,
498 CALL zgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
505 $ h( kwtop, kwtop-1 ) = s*dconjg( v( 1, 1 ) )
506 CALL zlacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
507 CALL zcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
513 IF( ns.GT.1 .AND. s.NE.zero )
514 $
CALL zunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
515 $ work( jw+1 ), lwork-jw, info )
524 DO 60 krow = ltop, kwtop - 1, nv
525 kln = min( nv, kwtop-krow )
526 CALL zgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
527 $ ldh, v, ldv, zero, wv, ldwv )
528 CALL zlacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
534 DO 70 kcol = kbot + 1, n, nh
535 kln = min( nh, n-kcol+1 )
536 CALL zgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
537 $ h( kwtop, kcol ), ldh, zero, t, ldt )
538 CALL zlacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
546 DO 80 krow = iloz, ihiz, nv
547 kln = min( nv, ihiz-krow+1 )
548 CALL zgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
549 $ ldz, v, ldv, zero, wv, ldwv )
550 CALL zlacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
570 work( 1 ) = dcmplx( lwkopt, 0 )
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
subroutine zgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZGEHRD
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
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 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 zlarf(side, m, n, v, incv, tau, c, ldc, work)
ZLARF applies an elementary reflector to a general rectangular matrix.
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
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.
subroutine ztrexc(compq, n, t, ldt, q, ldq, ifst, ilst, info)
ZTREXC
subroutine zunmhr(side, trans, m, n, ilo, ihi, a, lda, tau, c, ldc, work, lwork, info)
ZUNMHR