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
314 INTRINSIC abs, dble, int, max, min, sqrt
320 jw = min( nw, kbot-ktop+1 )
327 CALL dgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
328 lwk1 = int( work( 1 ) )
332 CALL dormhr(
'R',
'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
334 lwk2 = int( work( 1 ) )
338 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr, si, 1, jw,
339 $ v, ldv, work, -1, infqr )
340 lwk3 = int( work( 1 ) )
344 lwkopt = max( jw+max( lwk1, lwk2 ), lwk3 )
349 IF( lwork.EQ.-1 )
THEN
350 work( 1 ) = dble( lwkopt )
367 safmin = dlamch(
'SAFE MINIMUM' )
368 safmax = one / safmin
369 ulp = dlamch(
'PRECISION' )
370 smlnum = safmin*( dble( 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 sr( kwtop ) = h( kwtop, kwtop )
390 IF( abs( s ).LE.max( smlnum, ulp*abs( h( kwtop, kwtop ) ) ) )
395 $ h( kwtop, kwtop-1 ) = zero
407 CALL dlacpy(
'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
408 CALL dcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
410 CALL dlaset(
'A', jw, jw, zero, one, v, ldv )
411 nmin = ilaenv( 12,
'DLAQR3',
'SV', jw, 1, jw, lwork )
412 IF( jw.GT.nmin )
THEN
413 CALL dlaqr4( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
414 $ si( kwtop ), 1, jw, v, ldv, work, lwork, infqr )
416 CALL dlahqr( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
417 $ si( kwtop ), 1, jw, v, ldv, infqr )
427 $ t( jw, jw-2 ) = zero
434 IF( ilst.LE.ns )
THEN
438 bulge = t( ns, ns-1 ).NE.zero
443 IF( .NOT. bulge )
THEN
447 foo = abs( t( ns, ns ) )
450 IF( abs( s*v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
THEN
461 CALL dtrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, work,
469 foo = abs( t( ns, ns ) ) + sqrt( abs( t( ns, ns-1 ) ) )*
470 $ sqrt( abs( t( ns-1, ns ) ) )
473 IF( max( abs( s*v( 1, ns ) ), abs( s*v( 1, ns-1 ) ) ).LE.
474 $ max( smlnum, ulp*foo ) )
THEN
486 CALL dtrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, work,
519 ELSE IF( t( i+1, i ).EQ.zero )
THEN
527 evi = abs( t( i, i ) )
529 evi = abs( t( i, i ) ) + sqrt( abs( t( i+1, i ) ) )*
530 $ sqrt( abs( t( i, i+1 ) ) )
534 evk = abs( t( k, k ) )
535 ELSE IF( t( k+1, k ).EQ.zero )
THEN
536 evk = abs( t( k, k ) )
538 evk = abs( t( k, k ) ) + sqrt( abs( t( k+1, k ) ) )*
539 $ sqrt( abs( t( k, k+1 ) ) )
542 IF( evi.GE.evk )
THEN
548 CALL dtrexc(
'V', jw, t, ldt, v, ldv, ifst, ilst, work,
558 ELSE IF( t( i+1, i ).EQ.zero )
THEN
573 IF( i.GE.infqr+1 )
THEN
574 IF( i.EQ.infqr+1 )
THEN
575 sr( kwtop+i-1 ) = t( i, i )
576 si( kwtop+i-1 ) = zero
578 ELSE IF( t( i, i-1 ).EQ.zero )
THEN
579 sr( kwtop+i-1 ) = t( i, i )
580 si( kwtop+i-1 ) = zero
587 CALL dlanv2( aa, bb, cc, dd, sr( kwtop+i-2 ),
588 $ si( kwtop+i-2 ), sr( kwtop+i-1 ),
589 $ si( kwtop+i-1 ), cs, sn )
595 IF( ns.LT.jw .OR. s.EQ.zero )
THEN
596 IF( ns.GT.1 .AND. s.NE.zero )
THEN
600 CALL dcopy( ns, v, ldv, work, 1 )
602 CALL dlarfg( ns, beta, work( 2 ), 1, tau )
605 CALL dlaset(
'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
607 CALL dlarf(
'L', ns, jw, work, 1, tau, t, ldt,
609 CALL dlarf(
'R', ns, ns, work, 1, tau, t, ldt,
611 CALL dlarf(
'R', jw, ns, work, 1, tau, v, ldv,
614 CALL dgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
621 $ h( kwtop, kwtop-1 ) = s*v( 1, 1 )
622 CALL dlacpy(
'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
623 CALL dcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
629 IF( ns.GT.1 .AND. s.NE.zero )
630 $
CALL dormhr(
'R',
'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
631 $ work( jw+1 ), lwork-jw, info )
640 DO 70 krow = ltop, kwtop - 1, nv
641 kln = min( nv, kwtop-krow )
642 CALL dgemm(
'N',
'N', kln, jw, jw, one, h( krow, kwtop ),
643 $ ldh, v, ldv, zero, wv, ldwv )
644 CALL dlacpy(
'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
650 DO 80 kcol = kbot + 1, n, nh
651 kln = min( nh, n-kcol+1 )
652 CALL dgemm(
'C',
'N', jw, kln, jw, one, v, ldv,
653 $ h( kwtop, kcol ), ldh, zero, t, ldt )
654 CALL dlacpy(
'A', jw, kln, t, ldt, h( kwtop, kcol ),
662 DO 90 krow = iloz, ihiz, nv
663 kln = min( nv, ihiz-krow+1 )
664 CALL dgemm(
'N',
'N', kln, jw, jw, one, z( krow, kwtop ),
665 $ ldz, v, ldv, zero, wv, ldwv )
666 CALL dlacpy(
'A', kln, jw, wv, ldwv, z( krow, kwtop ),
686 work( 1 ) = dble( lwkopt )
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
subroutine dgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
DGEHRD
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
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 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 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 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 dlarfg(n, alpha, x, incx, tau)
DLARFG generates an elementary reflector (Householder matrix).
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 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