272 SUBROUTINE dlaqr3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ,
273 $ IHIZ, Z, LDZ, NS, ND, SR, SI, V, LDV, NH, T,
274 $ LDT, NV, WV, LDWV, WORK, LWORK )
281 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
282 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
286 DOUBLE PRECISION H( LDH, * ), SI( * ), SR( * ), T( LDT, * ),
287 $ V( LDV, * ), WORK( * ), WV( LDWV, * ),
293 DOUBLE PRECISION ZERO, ONE
294 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
297 DOUBLE PRECISION AA, BB, BETA, CC, CS, DD, EVI, EVK, FOO, S,
298 $ SAFMAX, SAFMIN, SMLNUM, SN, TAU, ULP
299 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, K, KCOL,
300 $ KEND, KLN, KROW, KWTOP, LTOP, LWK1, LWK2, LWK3,
302 LOGICAL BULGE, SORTED
305 DOUBLE PRECISION DLAMCH
307 EXTERNAL dlamch, ilaenv
315 INTRINSIC abs, dble, int, max, min, sqrt
321 jw = min( nw, kbot-ktop+1 )
328 CALL dgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
329 lwk1 = int( work( 1 ) )
333 CALL dormhr(
'R',
'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
335 lwk2 = int( work( 1 ) )
339 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr, si, 1, jw,
340 $ v, ldv, work, -1, infqr )
341 lwk3 = int( work( 1 ) )
345 lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
350 IF( lwork.EQ.-1 )
THEN
351 work( 1 ) = dble( lwkopt )
368 safmin = dlamch(
'SAFE MINIMUM' )
369 safmax = one / safmin
370 CALL dlabad( safmin, safmax )
371 ulp = dlamch(
'PRECISION' )
372 smlnum = safmin*( dble( n ) / ulp )
376 jw = min( nw, kbot-ktop+1 )
377 kwtop = kbot - jw + 1
378 IF( kwtop.EQ.ktop )
THEN
381 s = h( kwtop, kwtop-1 )
384 IF( kbot.EQ.kwtop )
THEN
388 sr( kwtop ) = h( kwtop, kwtop )
392 IF( abs( s ).LE.max( smlnum, ulp*abs( h( kwtop, kwtop ) ) ) )
397 $ h( kwtop, kwtop-1 ) = zero
409 CALL dlacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
410 CALL dcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
412 CALL dlaset(
'A', jw, jw, zero, one, v, ldv )
413 nmin = ilaenv( 12,
'DLAQR3',
'SV', jw, 1, jw, lwork )
414 IF( jw.GT.nmin )
THEN
415 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
416 $ si( kwtop ), 1, jw, v, ldv, work, lwork, infqr )
418 CALL dlahqr( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
419 $ si( kwtop ), 1, jw, v, ldv, infqr )
429 $ t( jw, jw-2 ) = zero
436 IF( ilst.LE.ns )
THEN
440 bulge = t( ns, ns-1 ).NE.zero
445 IF( .NOT. bulge )
THEN
449 foo = abs( t( ns, ns ) )
452 IF( abs( s*v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
THEN
463 CALL dtrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, work,
471 foo = abs( t( ns, ns ) ) + sqrt( abs( t( ns, ns-1 ) ) )*
472 $ sqrt( abs( t( ns-1, ns ) ) )
475 IF( max( abs( s*v( 1, ns ) ), abs( s*v( 1, ns-1 ) ) ).LE.
476 $ max( smlnum, ulp*foo ) )
THEN
488 CALL dtrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, work,
521 ELSE IF( t( i+1, i ).EQ.zero )
THEN
529 evi = abs( t( i, i ) )
531 evi = abs( t( i, i ) ) + sqrt( abs( t( i+1, i ) ) )*
532 $ sqrt( abs( t( i, i+1 ) ) )
536 evk = abs( t( k, k ) )
537 ELSE IF( t( k+1, k ).EQ.zero )
THEN
538 evk = abs( t( k, k ) )
540 evk = abs( t( k, k ) ) + sqrt( abs( t( k+1, k ) ) )*
541 $ sqrt( abs( t( k, k+1 ) ) )
544 IF( evi.GE.evk )
THEN
550 CALL dtrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, work,
560 ELSE IF( t( i+1, i ).EQ.zero )
THEN
575 IF( i.GE.infqr+1 )
THEN
576 IF( i.EQ.infqr+1 )
THEN
577 sr( kwtop+i-1 ) = t( i, i )
578 si( kwtop+i-1 ) = zero
580 ELSE IF( t( i, i-1 ).EQ.zero )
THEN
581 sr( kwtop+i-1 ) = t( i, i )
582 si( kwtop+i-1 ) = zero
589 CALL dlanv2( aa, bb, cc, dd, sr( kwtop+i-2 ),
590 $ si( kwtop+i-2 ), sr( kwtop+i-1 ),
591 $ si( kwtop+i-1 ), cs, sn )
597 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
598 IF( ns.GT.1 .AND. s.NE.zero )
THEN
602 CALL dcopy( ns, v, ldv, work, 1 )
604 CALL dlarfg( ns, beta, work( 2 ), 1, tau )
607 CALL dlaset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
609 CALL dlarf(
'L', ns, jw, work, 1, tau, t, ldt,
611 CALL dlarf(
'R', ns, ns, work, 1, tau, t, ldt,
613 CALL dlarf(
'R', jw, ns, work, 1, tau, v, ldv,
616 CALL dgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
623 $ h( kwtop, kwtop-1 ) = s*v( 1, 1 )
624 CALL dlacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
625 CALL dcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
631 IF( ns.GT.1 .AND. s.NE.zero )
632 $
CALL dormhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
633 $ work( jw+1 ), lwork-jw, info )
642 DO 70 krow = ltop, kwtop - 1, nv
643 kln = min( nv, kwtop-krow )
644 CALL dgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
645 $ ldh, v, ldv, zero, wv, ldwv )
646 CALL dlacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
652 DO 80 kcol = kbot + 1, n, nh
653 kln = min( nh, n-kcol+1 )
654 CALL dgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
655 $ h( kwtop, kcol ), ldh, zero, t, ldt )
656 CALL dlacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
664 DO 90 krow = iloz, ihiz, nv
665 kln = min( nv, ihiz-krow+1 )
666 CALL dgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
667 $ ldz, v, ldv, zero, wv, ldwv )
668 CALL dlacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
688 work( 1 ) = dble( lwkopt )
subroutine dlabad(SMALL, LARGE)
DLABAD
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DGEHRD
subroutine dlanv2(A, B, C, D, RT1R, RT1I, RT2R, RT2I, CS, SN)
DLANV2 computes the Schur factorization of a real 2-by-2 nonsymmetric matrix in standard form.
subroutine dlaqr4(WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO)
DLAQR4 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur de...
subroutine dlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
DLARF applies an elementary reflector to a general rectangular matrix.
subroutine dlaqr3(WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ, IHIZ, Z, LDZ, NS, ND, SR, SI, V, LDV, NH, T, LDT, NV, WV, LDWV, WORK, LWORK)
DLAQR3 performs the orthogonal similarity transformation of a Hessenberg matrix to detect and deflate...
subroutine dlahqr(WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, ILOZ, IHIZ, Z, LDZ, INFO)
DLAHQR computes the eigenvalues and Schur factorization of an upper Hessenberg matrix,...
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
subroutine dtrexc(COMPQ, N, T, LDT, Q, LDQ, IFST, ILST, WORK, INFO)
DTREXC
subroutine dormhr(SIDE, TRANS, M, N, ILO, IHI, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMHR