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 CALL slabad( safmin, safmax )
369 ulp = slamch(
'PRECISION' )
370 smlnum = safmin*( real( n ) / ulp )
374 jw = min( nw, kbot-ktop+1 )
375 kwtop = kbot - jw + 1
376 IF( kwtop.EQ.ktop )
THEN
379 s = h( kwtop, kwtop-1 )
382 IF( kbot.EQ.kwtop )
THEN
386 sh( kwtop ) = h( kwtop, kwtop )
389 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
394 $ h( kwtop, kwtop-1 ) = zero
406 CALL clacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
407 CALL ccopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
409 CALL claset(
'A', jw, jw, zero, one, v, ldv )
410 nmin = ilaenv( 12,
'CLAQR3',
'SV', jw, 1, jw, lwork )
411 IF( jw.GT.nmin )
THEN
412 CALL claqr4( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
413 $ jw, v, ldv, work, lwork, infqr )
415 CALL clahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
416 $ jw, v, ldv, infqr )
423 DO 10 knt = infqr + 1, jw
427 foo = cabs1( t( ns, ns ) )
430 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
442 CALL ctrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
457 DO 30 i = infqr + 1, ns
460 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
465 $
CALL ctrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, info )
471 DO 40 i = infqr + 1, jw
472 sh( kwtop+i-1 ) = t( i, i )
476 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
477 IF( ns.GT.1 .AND. s.NE.zero )
THEN
481 CALL ccopy( ns, v, ldv, work, 1 )
483 work( i ) = conjg( work( i ) )
486 CALL clarfg( ns, beta, work( 2 ), 1, tau )
489 CALL claset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
491 CALL clarf(
'L', ns, jw, work, 1, conjg( tau ), t, ldt,
493 CALL clarf(
'R', ns, ns, work, 1, tau, t, ldt,
495 CALL clarf(
'R', jw, ns, work, 1, tau, v, ldv,
498 CALL cgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
505 $ h( kwtop, kwtop-1 ) = s*conjg( v( 1, 1 ) )
506 CALL clacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
507 CALL ccopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
513 IF( ns.GT.1 .AND. s.NE.zero )
514 $
CALL cunmhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
515 $ work( jw+1 ), lwork-jw, info )
524 DO 60 krow = ltop, kwtop - 1, nv
525 kln = min( nv, kwtop-krow )
526 CALL cgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
527 $ ldh, v, ldv, zero, wv, ldwv )
528 CALL clacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
534 DO 70 kcol = kbot + 1, n, nh
535 kln = min( nh, n-kcol+1 )
536 CALL cgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
537 $ h( kwtop, kcol ), ldh, zero, t, ldt )
538 CALL clacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
546 DO 80 krow = iloz, ihiz, nv
547 kln = min( nv, ihiz-krow+1 )
548 CALL cgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
549 $ ldz, v, ldv, zero, wv, ldwv )
550 CALL clacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
570 work( 1 ) = cmplx( lwkopt, 0 )
subroutine slabad(SMALL, LARGE)
SLABAD
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
subroutine cgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
CGEHRD
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 clarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
CLARF applies an elementary reflector to a general rectangular matrix.
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 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 clarfg(N, ALPHA, X, INCX, TAU)
CLARFG generates an elementary reflector (Householder matrix).
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 clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
subroutine cunmhr(SIDE, TRANS, M, N, ILO, IHI, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMHR
subroutine ctrexc(COMPQ, N, T, LDT, Q, LDQ, IFST, ILST, INFO)
CTREXC