265 SUBROUTINE zlaqr2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH,
267 $ IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT,
268 $ NV, WV, LDWV, WORK, LWORK )
275 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
276 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
280 COMPLEX*16 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
281 $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
288 PARAMETER ( ZERO = ( 0.0d0, 0.0d0 ),
289 $ one = ( 1.0d0, 0.0d0 ) )
290 DOUBLE PRECISION RZERO, RONE
291 parameter( rzero = 0.0d0, rone = 1.0d0 )
294 COMPLEX*16 CDUM, S, TAU
295 DOUBLE PRECISION FOO, SAFMAX, SAFMIN, SMLNUM, ULP
296 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
297 $ knt, krow, kwtop, ltop, lwk1, lwk2, lwkopt
300 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,
336 lwk2 = int( work( 1 ) )
340 lwkopt = jw + max( lwk1, lwk2 )
345 IF( lwork.EQ.-1 )
THEN
346 work( 1 ) = dcmplx( lwkopt, 0 )
363 safmin = dlamch(
'SAFE MINIMUM' )
364 safmax = rone / safmin
365 ulp = dlamch(
'PRECISION' )
366 smlnum = safmin*( dble( n ) / ulp )
370 jw = min( nw, kbot-ktop+1 )
371 kwtop = kbot - jw + 1
372 IF( kwtop.EQ.ktop )
THEN
375 s = h( kwtop, kwtop-1 )
378 IF( kbot.EQ.kwtop )
THEN
382 sh( kwtop ) = h( kwtop, kwtop )
385 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
390 $ h( kwtop, kwtop-1 ) = zero
402 CALL zlacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
403 CALL zcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ),
406 CALL zlaset(
'A', jw, jw, zero, one, v, ldv )
407 CALL zlahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
408 $ jw, v, ldv, infqr )
414 DO 10 knt = infqr + 1, jw
418 foo = cabs1( t( ns, ns ) )
421 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
433 CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
448 DO 30 i = infqr + 1, ns
451 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
456 $
CALL ztrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst,
463 DO 40 i = infqr + 1, jw
464 sh( kwtop+i-1 ) = t( i, i )
468 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
469 IF( ns.GT.1 .AND. s.NE.zero )
THEN
473 CALL zcopy( ns, v, ldv, work, 1 )
475 work( i ) = dconjg( work( i ) )
477 CALL zlarfg( ns, work( 1 ), work( 2 ), 1, tau )
479 CALL zlaset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ),
482 CALL zlarf1f(
'L', ns, jw, work, 1, conjg( tau ), t, ldt,
484 CALL zlarf1f(
'R', ns, ns, work, 1, tau, t, ldt,
486 CALL zlarf1f(
'R', jw, ns, work, 1, tau, v, ldv,
489 CALL zgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
496 $ h( kwtop, kwtop-1 ) = s*dconjg( v( 1, 1 ) )
497 CALL zlacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
498 CALL zcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
504 IF( ns.GT.1 .AND. s.NE.zero )
505 $
CALL zunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v,
507 $ work( jw+1 ), lwork-jw, info )
516 DO 60 krow = ltop, kwtop - 1, nv
517 kln = min( nv, kwtop-krow )
518 CALL zgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
519 $ ldh, v, ldv, zero, wv, ldwv )
520 CALL zlacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ),
527 DO 70 kcol = kbot + 1, n, nh
528 kln = min( nh, n-kcol+1 )
529 CALL zgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
530 $ h( kwtop, kcol ), ldh, zero, t, ldt )
531 CALL zlacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
539 DO 80 krow = iloz, ihiz, nv
540 kln = min( nv, ihiz-krow+1 )
541 CALL zgemm(
'N',
'N', kln, jw, jw, one, z( krow,
543 $ ldz, v, ldv, zero, wv, ldwv )
544 CALL zlacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
564 work( 1 ) = dcmplx( lwkopt, 0 )