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 CALL dlabad( safmin, safmax )
370 ulp = dlamch(
'PRECISION' )
371 smlnum = safmin*( dble( n ) / ulp )
375 jw = min( nw, kbot-ktop+1 )
376 kwtop = kbot - jw + 1
377 IF( kwtop.EQ.ktop )
THEN
380 s = h( kwtop, kwtop-1 )
383 IF( kbot.EQ.kwtop )
THEN
387 sh( kwtop ) = h( kwtop, kwtop )
390 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
395 $ h( kwtop, kwtop-1 ) = zero
407 CALL zlacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
408 CALL zcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
410 CALL zlaset(
'A', jw, jw, zero, one, v, ldv )
411 nmin = ilaenv( 12,
'ZLAQR3',
'SV', jw, 1, jw, lwork )
412 IF( jw.GT.nmin )
THEN
413 CALL zlaqr4( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
414 $ jw, v, ldv, work, lwork, infqr )
416 CALL zlahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
417 $ jw, v, ldv, infqr )
424 DO 10 knt = infqr + 1, jw
428 foo = cabs1( t( ns, ns ) )
431 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
443 CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
458 DO 30 i = infqr + 1, ns
461 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
466 $
CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
472 DO 40 i = infqr + 1, jw
473 sh( kwtop+i-1 ) = t( i, i )
477 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
478 IF( ns.GT.1 .AND. s.NE.zero )
THEN
482 CALL zcopy( ns, v, ldv, work, 1 )
484 work( i ) = dconjg( work( i ) )
487 CALL zlarfg( ns, beta, work( 2 ), 1, tau )
490 CALL zlaset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
492 CALL zlarf(
'L', ns, jw, work, 1, dconjg( tau ), t, ldt,
494 CALL zlarf(
'R', ns, ns, work, 1, tau, t, ldt,
496 CALL zlarf(
'R', jw, ns, work, 1, tau, v, ldv,
499 CALL zgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
506 $ h( kwtop, kwtop-1 ) = s*dconjg( v( 1, 1 ) )
507 CALL zlacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
508 CALL zcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
514 IF( ns.GT.1 .AND. s.NE.zero )
515 $
CALL zunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
516 $ work( jw+1 ), lwork-jw, info )
525 DO 60 krow = ltop, kwtop - 1, nv
526 kln = min( nv, kwtop-krow )
527 CALL zgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
528 $ ldh, v, ldv, zero, wv, ldwv )
529 CALL zlacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
535 DO 70 kcol = kbot + 1, n, nh
536 kln = min( nh, n-kcol+1 )
537 CALL zgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
538 $ h( kwtop, kcol ), ldh, zero, t, ldt )
539 CALL zlacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
547 DO 80 krow = iloz, ihiz, nv
548 kln = min( nv, ihiz-krow+1 )
549 CALL zgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
550 $ ldz, v, ldv, zero, wv, ldwv )
551 CALL zlacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
571 work( 1 ) = dcmplx( lwkopt, 0 )
subroutine dlabad(SMALL, LARGE)
DLABAD
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
subroutine zgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
ZGEHRD
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 zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
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 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 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 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 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