266 SUBROUTINE claqr2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
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 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
281 $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
288 PARAMETER ( ZERO = ( 0.0e0, 0.0e0 ),
289 $ one = ( 1.0e0, 0.0e0 ) )
291 PARAMETER ( RZERO = 0.0e0, rone = 1.0e0 )
294 COMPLEX BETA, CDUM, S, TAU
295 REAL FOO, SAFMAX, SAFMIN, SMLNUM, ULP
296 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
297 $ knt, krow, kwtop, ltop, lwk1, lwk2, lwkopt
308 INTRINSIC abs, aimag, cmplx, conjg, int, max, min, real
314 cabs1( cdum ) = abs( real( cdum ) ) + abs( aimag( cdum ) )
320 jw = min( nw, kbot-ktop+1 )
327 CALL cgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
328 lwk1 = int( work( 1 ) )
332 CALL cunmhr(
'R',
'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
334 lwk2 = int( work( 1 ) )
338 lwkopt = jw + max( lwk1, lwk2 )
343 IF( lwork.EQ.-1 )
THEN
344 work( 1 ) = cmplx( lwkopt, 0 )
361 safmin = slamch(
'SAFE MINIMUM' )
362 safmax = rone / safmin
363 ulp = slamch(
'PRECISION' )
364 smlnum = safmin*( real( n ) / ulp )
368 jw = min( nw, kbot-ktop+1 )
369 kwtop = kbot - jw + 1
370 IF( kwtop.EQ.ktop )
THEN
373 s = h( kwtop, kwtop-1 )
376 IF( kbot.EQ.kwtop )
THEN
380 sh( kwtop ) = h( kwtop, kwtop )
383 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
388 $ h( kwtop, kwtop-1 ) = zero
400 CALL clacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
401 CALL ccopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
403 CALL claset(
'A', jw, jw, zero, one, v, ldv )
404 CALL clahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
405 $ jw, v, ldv, infqr )
411 DO 10 knt = infqr + 1, jw
415 foo = cabs1( t( ns, ns ) )
418 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
430 CALL ctrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
445 DO 30 i = infqr + 1, ns
448 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
453 $
CALL ctrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
459 DO 40 i = infqr + 1, jw
460 sh( kwtop+i-1 ) = t( i, i )
464 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
465 IF( ns.GT.1 .AND. s.NE.zero )
THEN
469 CALL ccopy( ns, v, ldv, work, 1 )
471 work( i ) = conjg( work( i ) )
474 CALL clarfg( ns, beta, work( 2 ), 1, tau )
477 CALL claset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
479 CALL clarf(
'L', ns, jw, work, 1, conjg( tau ), t, ldt,
481 CALL clarf(
'R', ns, ns, work, 1, tau, t, ldt,
483 CALL clarf(
'R', jw, ns, work, 1, tau, v, ldv,
486 CALL cgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
493 $ h( kwtop, kwtop-1 ) = s*conjg( v( 1, 1 ) )
494 CALL clacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
495 CALL ccopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
501 IF( ns.GT.1 .AND. s.NE.zero )
502 $
CALL cunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
503 $ work( jw+1 ), lwork-jw, info )
512 DO 60 krow = ltop, kwtop - 1, nv
513 kln = min( nv, kwtop-krow )
514 CALL cgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
515 $ ldh, v, ldv, zero, wv, ldwv )
516 CALL clacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
522 DO 70 kcol = kbot + 1, n, nh
523 kln = min( nh, n-kcol+1 )
524 CALL cgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
525 $ h( kwtop, kcol ), ldh, zero, t, ldt )
526 CALL clacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
534 DO 80 krow = iloz, ihiz, nv
535 kln = min( nv, ihiz-krow+1 )
536 CALL cgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
537 $ ldz, v, ldv, zero, wv, ldwv )
538 CALL clacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
558 work( 1 ) = cmplx( lwkopt, 0 )
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
subroutine cgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
CGEHRD
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
subroutine clahqr(wantt, wantz, n, ilo, ihi, h, ldh, w, iloz, ihiz, z, ldz, info)
CLAHQR computes the eigenvalues and Schur factorization of an upper Hessenberg matrix,...
subroutine claqr2(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)
CLAQR2 performs the unitary similarity transformation of a Hessenberg matrix to detect and deflate fu...
subroutine clarf(side, m, n, v, incv, tau, c, ldc, work)
CLARF applies an elementary reflector to a general rectangular matrix.
subroutine clarfg(n, alpha, x, incx, tau)
CLARFG generates an elementary reflector (Householder matrix).
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
subroutine ctrexc(compq, n, t, ldt, q, ldq, ifst, ilst, info)
CTREXC
subroutine cunmhr(side, trans, m, n, ilo, ihi, a, lda, tau, c, ldc, work, lwork, info)
CUNMHR