267 SUBROUTINE zlaqr2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
268 $ IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT,
269 $ NV, WV, LDWV, WORK, LWORK )
276 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
277 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
281 COMPLEX*16 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
282 $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
289 PARAMETER ( ZERO = ( 0.0d0, 0.0d0 ),
290 $ one = ( 1.0d0, 0.0d0 ) )
291 DOUBLE PRECISION RZERO, RONE
292 PARAMETER ( RZERO = 0.0d0, rone = 1.0d0 )
295 COMPLEX*16 BETA, CDUM, S, TAU
296 DOUBLE PRECISION FOO, SAFMAX, SAFMIN, SMLNUM, ULP
297 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
298 $ knt, krow, kwtop, ltop, lwk1, lwk2, lwkopt
301 DOUBLE PRECISION DLAMCH
309 INTRINSIC abs, dble, dcmplx, dconjg, dimag, int, max, min
312 DOUBLE PRECISION CABS1
315 cabs1( cdum ) = abs( dble( cdum ) ) + abs( dimag( cdum ) )
321 jw = min( nw, kbot-ktop+1 )
328 CALL zgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
329 lwk1 = int( work( 1 ) )
333 CALL zunmhr(
'R',
'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
335 lwk2 = int( work( 1 ) )
339 lwkopt = jw + max( lwk1, lwk2 )
344 IF( lwork.EQ.-1 )
THEN
345 work( 1 ) = dcmplx( lwkopt, 0 )
362 safmin = dlamch(
'SAFE MINIMUM' )
363 safmax = rone / safmin
364 ulp = dlamch(
'PRECISION' )
365 smlnum = safmin*( dble( n ) / ulp )
369 jw = min( nw, kbot-ktop+1 )
370 kwtop = kbot - jw + 1
371 IF( kwtop.EQ.ktop )
THEN
374 s = h( kwtop, kwtop-1 )
377 IF( kbot.EQ.kwtop )
THEN
381 sh( kwtop ) = h( kwtop, kwtop )
384 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
389 $ h( kwtop, kwtop-1 ) = zero
401 CALL zlacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
402 CALL zcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
404 CALL zlaset(
'A', jw, jw, zero, one, v, ldv )
405 CALL zlahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
406 $ jw, v, ldv, infqr )
412 DO 10 knt = infqr + 1, jw
416 foo = cabs1( t( ns, ns ) )
419 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
431 CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
446 DO 30 i = infqr + 1, ns
449 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
454 $
CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
460 DO 40 i = infqr + 1, jw
461 sh( kwtop+i-1 ) = t( i, i )
465 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
466 IF( ns.GT.1 .AND. s.NE.zero )
THEN
470 CALL zcopy( ns, v, ldv, work, 1 )
472 work( i ) = dconjg( work( i ) )
475 CALL zlarfg( ns, beta, work( 2 ), 1, tau )
478 CALL zlaset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
480 CALL zlarf(
'L', ns, jw, work, 1, dconjg( tau ), t, ldt,
482 CALL zlarf(
'R', ns, ns, work, 1, tau, t, ldt,
484 CALL zlarf(
'R', jw, ns, work, 1, tau, v, ldv,
487 CALL zgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
494 $ h( kwtop, kwtop-1 ) = s*dconjg( v( 1, 1 ) )
495 CALL zlacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
496 CALL zcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
502 IF( ns.GT.1 .AND. s.NE.zero )
503 $
CALL zunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
504 $ work( jw+1 ), lwork-jw, info )
513 DO 60 krow = ltop, kwtop - 1, nv
514 kln = min( nv, kwtop-krow )
515 CALL zgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
516 $ ldh, v, ldv, zero, wv, ldwv )
517 CALL zlacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
523 DO 70 kcol = kbot + 1, n, nh
524 kln = min( nh, n-kcol+1 )
525 CALL zgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
526 $ h( kwtop, kcol ), ldh, zero, t, ldt )
527 CALL zlacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
535 DO 80 krow = iloz, ihiz, nv
536 kln = min( nv, ihiz-krow+1 )
537 CALL zgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
538 $ ldz, v, ldv, zero, wv, ldwv )
539 CALL zlacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
559 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 zlaqr2(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)
ZLAQR2 performs the unitary similarity transformation of a Hessenberg matrix to detect and deflate fu...
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