263 SUBROUTINE claqr3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
264 $ IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT,
265 $ NV, WV, LDWV, WORK, LWORK )
272 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
273 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
277 COMPLEX H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
278 $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
285 PARAMETER ( ZERO = ( 0.0e0, 0.0e0 ),
286 $ one = ( 1.0e0, 0.0e0 ) )
288 PARAMETER ( RZERO = 0.0e0, rone = 1.0e0 )
291 COMPLEX BETA, CDUM, S, TAU
292 REAL FOO, SAFMAX, SAFMIN, SMLNUM, ULP
293 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
294 $ knt, krow, kwtop, ltop, lwk1, lwk2, lwk3,
300 EXTERNAL slamch, ilaenv
307 INTRINSIC abs, aimag, cmplx, conjg, int, max, min, real
313 cabs1( cdum ) = abs( real( cdum ) ) + abs( aimag( cdum ) )
319 jw = min( nw, kbot-ktop+1 )
326 CALL cgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
327 lwk1 = int( work( 1 ) )
331 CALL cunmhr(
'R',
'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
333 lwk2 = int( work( 1 ) )
337 CALL claqr4( .true., .true., jw, 1, jw, t, ldt, sh, 1, jw, v,
338 $ ldv, work, -1, infqr )
339 lwk3 = int( work( 1 ) )
343 lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
348 IF( lwork.EQ.-1 )
THEN
349 work( 1 ) = cmplx( lwkopt, 0 )
366 safmin = slamch(
'SAFE MINIMUM' )
367 safmax = rone / safmin
368 ulp = slamch(
'PRECISION' )
369 smlnum = safmin*( real( n ) / ulp )
373 jw = min( nw, kbot-ktop+1 )
374 kwtop = kbot - jw + 1
375 IF( kwtop.EQ.ktop )
THEN
378 s = h( kwtop, kwtop-1 )
381 IF( kbot.EQ.kwtop )
THEN
385 sh( kwtop ) = h( kwtop, kwtop )
388 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
393 $ h( kwtop, kwtop-1 ) = zero
405 CALL clacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
406 CALL ccopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
408 CALL claset(
'A', jw, jw, zero, one, v, ldv )
409 nmin = ilaenv( 12,
'CLAQR3',
'SV', jw, 1, jw, lwork )
410 IF( jw.GT.nmin )
THEN
411 CALL claqr4( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
412 $ jw, v, ldv, work, lwork, infqr )
414 CALL clahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
415 $ jw, v, ldv, infqr )
422 DO 10 knt = infqr + 1, jw
426 foo = cabs1( t( ns, ns ) )
429 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
441 CALL ctrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
456 DO 30 i = infqr + 1, ns
459 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
464 $
CALL ctrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
470 DO 40 i = infqr + 1, jw
471 sh( kwtop+i-1 ) = t( i, i )
475 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
476 IF( ns.GT.1 .AND. s.NE.zero )
THEN
480 CALL ccopy( ns, v, ldv, work, 1 )
482 work( i ) = conjg( work( i ) )
485 CALL clarfg( ns, beta, work( 2 ), 1, tau )
488 CALL claset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
490 CALL clarf(
'L', ns, jw, work, 1, conjg( tau ), t, ldt,
492 CALL clarf(
'R', ns, ns, work, 1, tau, t, ldt,
494 CALL clarf(
'R', jw, ns, work, 1, tau, v, ldv,
497 CALL cgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
504 $ h( kwtop, kwtop-1 ) = s*conjg( v( 1, 1 ) )
505 CALL clacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
506 CALL ccopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
512 IF( ns.GT.1 .AND. s.NE.zero )
513 $
CALL cunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
514 $ work( jw+1 ), lwork-jw, info )
523 DO 60 krow = ltop, kwtop - 1, nv
524 kln = min( nv, kwtop-krow )
525 CALL cgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
526 $ ldh, v, ldv, zero, wv, ldwv )
527 CALL clacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
533 DO 70 kcol = kbot + 1, n, nh
534 kln = min( nh, n-kcol+1 )
535 CALL cgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
536 $ h( kwtop, kcol ), ldh, zero, t, ldt )
537 CALL clacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
545 DO 80 krow = iloz, ihiz, nv
546 kln = min( nv, ihiz-krow+1 )
547 CALL cgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
548 $ ldz, v, ldv, zero, wv, ldwv )
549 CALL clacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
569 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 claqr3(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)
CLAQR3 performs the unitary similarity transformation of a Hessenberg matrix to detect and deflate fu...
subroutine claqr4(wantt, wantz, n, ilo, ihi, h, ldh, w, iloz, ihiz, z, ldz, work, lwork, info)
CLAQR4 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur de...
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