LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ dlaqr2()

subroutine dlaqr2 ( logical  WANTT,
logical  WANTZ,
integer  N,
integer  KTOP,
integer  KBOT,
integer  NW,
double precision, dimension( ldh, * )  H,
integer  LDH,
integer  ILOZ,
integer  IHIZ,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
integer  NS,
integer  ND,
double precision, dimension( * )  SR,
double precision, dimension( * )  SI,
double precision, dimension( ldv, * )  V,
integer  LDV,
integer  NH,
double precision, dimension( ldt, * )  T,
integer  LDT,
integer  NV,
double precision, dimension( ldwv, * )  WV,
integer  LDWV,
double precision, dimension( * )  WORK,
integer  LWORK 
)

DLAQR2 performs the orthogonal similarity transformation of a Hessenberg matrix to detect and deflate fully converged eigenvalues from a trailing principal submatrix (aggressive early deflation).

Download DLAQR2 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
    DLAQR2 is identical to DLAQR3 except that it avoids
    recursion by calling DLAHQR instead of DLAQR4.

    Aggressive early deflation:

    This subroutine accepts as input an upper Hessenberg matrix
    H and performs an orthogonal similarity transformation
    designed to detect and deflate fully converged eigenvalues from
    a trailing principal submatrix.  On output H has been over-
    written by a new Hessenberg matrix that is a perturbation of
    an orthogonal similarity transformation of H.  It is to be
    hoped that the final version of H has many zero subdiagonal
    entries.
Parameters
[in]WANTT
          WANTT is LOGICAL
          If .TRUE., then the Hessenberg matrix H is fully updated
          so that the quasi-triangular Schur factor may be
          computed (in cooperation with the calling subroutine).
          If .FALSE., then only enough of H is updated to preserve
          the eigenvalues.
[in]WANTZ
          WANTZ is LOGICAL
          If .TRUE., then the orthogonal matrix Z is updated so
          so that the orthogonal Schur factor may be computed
          (in cooperation with the calling subroutine).
          If .FALSE., then Z is not referenced.
[in]N
          N is INTEGER
          The order of the matrix H and (if WANTZ is .TRUE.) the
          order of the orthogonal matrix Z.
[in]KTOP
          KTOP is INTEGER
          It is assumed that either KTOP = 1 or H(KTOP,KTOP-1)=0.
          KBOT and KTOP together determine an isolated block
          along the diagonal of the Hessenberg matrix.
[in]KBOT
          KBOT is INTEGER
          It is assumed without a check that either
          KBOT = N or H(KBOT+1,KBOT)=0.  KBOT and KTOP together
          determine an isolated block along the diagonal of the
          Hessenberg matrix.
[in]NW
          NW is INTEGER
          Deflation window size.  1 <= NW <= (KBOT-KTOP+1).
[in,out]H
          H is DOUBLE PRECISION array, dimension (LDH,N)
          On input the initial N-by-N section of H stores the
          Hessenberg matrix undergoing aggressive early deflation.
          On output H has been transformed by an orthogonal
          similarity transformation, perturbed, and the returned
          to Hessenberg form that (it is to be hoped) has some
          zero subdiagonal entries.
[in]LDH
          LDH is INTEGER
          Leading dimension of H just as declared in the calling
          subroutine.  N <= LDH
[in]ILOZ
          ILOZ is INTEGER
[in]IHIZ
          IHIZ is INTEGER
          Specify the rows of Z to which transformations must be
          applied if WANTZ is .TRUE.. 1 <= ILOZ <= IHIZ <= N.
[in,out]Z
          Z is DOUBLE PRECISION array, dimension (LDZ,N)
          IF WANTZ is .TRUE., then on output, the orthogonal
          similarity transformation mentioned above has been
          accumulated into Z(ILOZ:IHIZ,ILOZ:IHIZ) from the right.
          If WANTZ is .FALSE., then Z is unreferenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of Z just as declared in the
          calling subroutine.  1 <= LDZ.
[out]NS
          NS is INTEGER
          The number of unconverged (ie approximate) eigenvalues
          returned in SR and SI that may be used as shifts by the
          calling subroutine.
[out]ND
          ND is INTEGER
          The number of converged eigenvalues uncovered by this
          subroutine.
[out]SR
          SR is DOUBLE PRECISION array, dimension (KBOT)
[out]SI
          SI is DOUBLE PRECISION array, dimension (KBOT)
          On output, the real and imaginary parts of approximate
          eigenvalues that may be used for shifts are stored in
          SR(KBOT-ND-NS+1) through SR(KBOT-ND) and
          SI(KBOT-ND-NS+1) through SI(KBOT-ND), respectively.
          The real and imaginary parts of converged eigenvalues
          are stored in SR(KBOT-ND+1) through SR(KBOT) and
          SI(KBOT-ND+1) through SI(KBOT), respectively.
[out]V
          V is DOUBLE PRECISION array, dimension (LDV,NW)
          An NW-by-NW work array.
[in]LDV
          LDV is INTEGER
          The leading dimension of V just as declared in the
          calling subroutine.  NW <= LDV
[in]NH
          NH is INTEGER
          The number of columns of T.  NH >= NW.
[out]T
          T is DOUBLE PRECISION array, dimension (LDT,NW)
[in]LDT
          LDT is INTEGER
          The leading dimension of T just as declared in the
          calling subroutine.  NW <= LDT
[in]NV
          NV is INTEGER
          The number of rows of work array WV available for
          workspace.  NV >= NW.
[out]WV
          WV is DOUBLE PRECISION array, dimension (LDWV,NW)
[in]LDWV
          LDWV is INTEGER
          The leading dimension of W just as declared in the
          calling subroutine.  NW <= LDV
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LWORK)
          On exit, WORK(1) is set to an estimate of the optimal value
          of LWORK for the given values of N, NW, KTOP and KBOT.
[in]LWORK
          LWORK is INTEGER
          The dimension of the work array WORK.  LWORK = 2*NW
          suffices, but greater efficiency may result from larger
          values of LWORK.

          If LWORK = -1, then a workspace query is assumed; DLAQR2
          only estimates the optimal workspace size for the given
          values of N, NW, KTOP and KBOT.  The estimate is returned
          in WORK(1).  No error message related to LWORK is issued
          by XERBLA.  Neither H nor Z are accessed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Karen Braman and Ralph Byers, Department of Mathematics, University of Kansas, USA

Definition at line 275 of file dlaqr2.f.

278 *
279 * -- LAPACK auxiliary routine --
280 * -- LAPACK is a software package provided by Univ. of Tennessee, --
281 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
282 *
283 * .. Scalar Arguments ..
284  INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
285  $ LDZ, LWORK, N, ND, NH, NS, NV, NW
286  LOGICAL WANTT, WANTZ
287 * ..
288 * .. Array Arguments ..
289  DOUBLE PRECISION H( LDH, * ), SI( * ), SR( * ), T( LDT, * ),
290  $ V( LDV, * ), WORK( * ), WV( LDWV, * ),
291  $ Z( LDZ, * )
292 * ..
293 *
294 * ================================================================
295 * .. Parameters ..
296  DOUBLE PRECISION ZERO, ONE
297  parameter( zero = 0.0d0, one = 1.0d0 )
298 * ..
299 * .. Local Scalars ..
300  DOUBLE PRECISION AA, BB, BETA, CC, CS, DD, EVI, EVK, FOO, S,
301  $ SAFMAX, SAFMIN, SMLNUM, SN, TAU, ULP
302  INTEGER I, IFST, ILST, INFO, INFQR, J, JW, K, KCOL,
303  $ KEND, KLN, KROW, KWTOP, LTOP, LWK1, LWK2,
304  $ LWKOPT
305  LOGICAL BULGE, SORTED
306 * ..
307 * .. External Functions ..
308  DOUBLE PRECISION DLAMCH
309  EXTERNAL dlamch
310 * ..
311 * .. External Subroutines ..
312  EXTERNAL dcopy, dgehrd, dgemm, dlabad, dlacpy, dlahqr,
314 * ..
315 * .. Intrinsic Functions ..
316  INTRINSIC abs, dble, int, max, min, sqrt
317 * ..
318 * .. Executable Statements ..
319 *
320 * ==== Estimate optimal workspace. ====
321 *
322  jw = min( nw, kbot-ktop+1 )
323  IF( jw.LE.2 ) THEN
324  lwkopt = 1
325  ELSE
326 *
327 * ==== Workspace query call to DGEHRD ====
328 *
329  CALL dgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
330  lwk1 = int( work( 1 ) )
331 *
332 * ==== Workspace query call to DORMHR ====
333 *
334  CALL dormhr( 'R', 'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
335  $ work, -1, info )
336  lwk2 = int( work( 1 ) )
337 *
338 * ==== Optimal workspace ====
339 *
340  lwkopt = jw + max( lwk1, lwk2 )
341  END IF
342 *
343 * ==== Quick return in case of workspace query. ====
344 *
345  IF( lwork.EQ.-1 ) THEN
346  work( 1 ) = dble( lwkopt )
347  RETURN
348  END IF
349 *
350 * ==== Nothing to do ...
351 * ... for an empty active block ... ====
352  ns = 0
353  nd = 0
354  work( 1 ) = one
355  IF( ktop.GT.kbot )
356  $ RETURN
357 * ... nor for an empty deflation window. ====
358  IF( nw.LT.1 )
359  $ RETURN
360 *
361 * ==== Machine constants ====
362 *
363  safmin = dlamch( 'SAFE MINIMUM' )
364  safmax = one / safmin
365  CALL dlabad( safmin, safmax )
366  ulp = dlamch( 'PRECISION' )
367  smlnum = safmin*( dble( n ) / ulp )
368 *
369 * ==== Setup deflation window ====
370 *
371  jw = min( nw, kbot-ktop+1 )
372  kwtop = kbot - jw + 1
373  IF( kwtop.EQ.ktop ) THEN
374  s = zero
375  ELSE
376  s = h( kwtop, kwtop-1 )
377  END IF
378 *
379  IF( kbot.EQ.kwtop ) THEN
380 *
381 * ==== 1-by-1 deflation window: not much to do ====
382 *
383  sr( kwtop ) = h( kwtop, kwtop )
384  si( kwtop ) = zero
385  ns = 1
386  nd = 0
387  IF( abs( s ).LE.max( smlnum, ulp*abs( h( kwtop, kwtop ) ) ) )
388  $ THEN
389  ns = 0
390  nd = 1
391  IF( kwtop.GT.ktop )
392  $ h( kwtop, kwtop-1 ) = zero
393  END IF
394  work( 1 ) = one
395  RETURN
396  END IF
397 *
398 * ==== Convert to spike-triangular form. (In case of a
399 * . rare QR failure, this routine continues to do
400 * . aggressive early deflation using that part of
401 * . the deflation window that converged using INFQR
402 * . here and there to keep track.) ====
403 *
404  CALL dlacpy( 'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
405  CALL dcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
406 *
407  CALL dlaset( 'A', jw, jw, zero, one, v, ldv )
408  CALL dlahqr( .true., .true., jw, 1, jw, t, ldt, sr( kwtop ),
409  $ si( kwtop ), 1, jw, v, ldv, infqr )
410 *
411 * ==== DTREXC needs a clean margin near the diagonal ====
412 *
413  DO 10 j = 1, jw - 3
414  t( j+2, j ) = zero
415  t( j+3, j ) = zero
416  10 CONTINUE
417  IF( jw.GT.2 )
418  $ t( jw, jw-2 ) = zero
419 *
420 * ==== Deflation detection loop ====
421 *
422  ns = jw
423  ilst = infqr + 1
424  20 CONTINUE
425  IF( ilst.LE.ns ) THEN
426  IF( ns.EQ.1 ) THEN
427  bulge = .false.
428  ELSE
429  bulge = t( ns, ns-1 ).NE.zero
430  END IF
431 *
432 * ==== Small spike tip test for deflation ====
433 *
434  IF( .NOT.bulge ) THEN
435 *
436 * ==== Real eigenvalue ====
437 *
438  foo = abs( t( ns, ns ) )
439  IF( foo.EQ.zero )
440  $ foo = abs( s )
441  IF( abs( s*v( 1, ns ) ).LE.max( smlnum, ulp*foo ) ) THEN
442 *
443 * ==== Deflatable ====
444 *
445  ns = ns - 1
446  ELSE
447 *
448 * ==== Undeflatable. Move it up out of the way.
449 * . (DTREXC can not fail in this case.) ====
450 *
451  ifst = ns
452  CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
453  $ info )
454  ilst = ilst + 1
455  END IF
456  ELSE
457 *
458 * ==== Complex conjugate pair ====
459 *
460  foo = abs( t( ns, ns ) ) + sqrt( abs( t( ns, ns-1 ) ) )*
461  $ sqrt( abs( t( ns-1, ns ) ) )
462  IF( foo.EQ.zero )
463  $ foo = abs( s )
464  IF( max( abs( s*v( 1, ns ) ), abs( s*v( 1, ns-1 ) ) ).LE.
465  $ max( smlnum, ulp*foo ) ) THEN
466 *
467 * ==== Deflatable ====
468 *
469  ns = ns - 2
470  ELSE
471 *
472 * ==== Undeflatable. Move them up out of the way.
473 * . Fortunately, DTREXC does the right thing with
474 * . ILST in case of a rare exchange failure. ====
475 *
476  ifst = ns
477  CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
478  $ info )
479  ilst = ilst + 2
480  END IF
481  END IF
482 *
483 * ==== End deflation detection loop ====
484 *
485  GO TO 20
486  END IF
487 *
488 * ==== Return to Hessenberg form ====
489 *
490  IF( ns.EQ.0 )
491  $ s = zero
492 *
493  IF( ns.LT.jw ) THEN
494 *
495 * ==== sorting diagonal blocks of T improves accuracy for
496 * . graded matrices. Bubble sort deals well with
497 * . exchange failures. ====
498 *
499  sorted = .false.
500  i = ns + 1
501  30 CONTINUE
502  IF( sorted )
503  $ GO TO 50
504  sorted = .true.
505 *
506  kend = i - 1
507  i = infqr + 1
508  IF( i.EQ.ns ) THEN
509  k = i + 1
510  ELSE IF( t( i+1, i ).EQ.zero ) THEN
511  k = i + 1
512  ELSE
513  k = i + 2
514  END IF
515  40 CONTINUE
516  IF( k.LE.kend ) THEN
517  IF( k.EQ.i+1 ) THEN
518  evi = abs( t( i, i ) )
519  ELSE
520  evi = abs( t( i, i ) ) + sqrt( abs( t( i+1, i ) ) )*
521  $ sqrt( abs( t( i, i+1 ) ) )
522  END IF
523 *
524  IF( k.EQ.kend ) THEN
525  evk = abs( t( k, k ) )
526  ELSE IF( t( k+1, k ).EQ.zero ) THEN
527  evk = abs( t( k, k ) )
528  ELSE
529  evk = abs( t( k, k ) ) + sqrt( abs( t( k+1, k ) ) )*
530  $ sqrt( abs( t( k, k+1 ) ) )
531  END IF
532 *
533  IF( evi.GE.evk ) THEN
534  i = k
535  ELSE
536  sorted = .false.
537  ifst = i
538  ilst = k
539  CALL dtrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, work,
540  $ info )
541  IF( info.EQ.0 ) THEN
542  i = ilst
543  ELSE
544  i = k
545  END IF
546  END IF
547  IF( i.EQ.kend ) THEN
548  k = i + 1
549  ELSE IF( t( i+1, i ).EQ.zero ) THEN
550  k = i + 1
551  ELSE
552  k = i + 2
553  END IF
554  GO TO 40
555  END IF
556  GO TO 30
557  50 CONTINUE
558  END IF
559 *
560 * ==== Restore shift/eigenvalue array from T ====
561 *
562  i = jw
563  60 CONTINUE
564  IF( i.GE.infqr+1 ) THEN
565  IF( i.EQ.infqr+1 ) THEN
566  sr( kwtop+i-1 ) = t( i, i )
567  si( kwtop+i-1 ) = zero
568  i = i - 1
569  ELSE IF( t( i, i-1 ).EQ.zero ) THEN
570  sr( kwtop+i-1 ) = t( i, i )
571  si( kwtop+i-1 ) = zero
572  i = i - 1
573  ELSE
574  aa = t( i-1, i-1 )
575  cc = t( i, i-1 )
576  bb = t( i-1, i )
577  dd = t( i, i )
578  CALL dlanv2( aa, bb, cc, dd, sr( kwtop+i-2 ),
579  $ si( kwtop+i-2 ), sr( kwtop+i-1 ),
580  $ si( kwtop+i-1 ), cs, sn )
581  i = i - 2
582  END IF
583  GO TO 60
584  END IF
585 *
586  IF( ns.LT.jw .OR. s.EQ.zero ) THEN
587  IF( ns.GT.1 .AND. s.NE.zero ) THEN
588 *
589 * ==== Reflect spike back into lower triangle ====
590 *
591  CALL dcopy( ns, v, ldv, work, 1 )
592  beta = work( 1 )
593  CALL dlarfg( ns, beta, work( 2 ), 1, tau )
594  work( 1 ) = one
595 *
596  CALL dlaset( 'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
597 *
598  CALL dlarf( 'L', ns, jw, work, 1, tau, t, ldt,
599  $ work( jw+1 ) )
600  CALL dlarf( 'R', ns, ns, work, 1, tau, t, ldt,
601  $ work( jw+1 ) )
602  CALL dlarf( 'R', jw, ns, work, 1, tau, v, ldv,
603  $ work( jw+1 ) )
604 *
605  CALL dgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
606  $ lwork-jw, info )
607  END IF
608 *
609 * ==== Copy updated reduced window into place ====
610 *
611  IF( kwtop.GT.1 )
612  $ h( kwtop, kwtop-1 ) = s*v( 1, 1 )
613  CALL dlacpy( 'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
614  CALL dcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
615  $ ldh+1 )
616 *
617 * ==== Accumulate orthogonal matrix in order update
618 * . H and Z, if requested. ====
619 *
620  IF( ns.GT.1 .AND. s.NE.zero )
621  $ CALL dormhr( 'R', 'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
622  $ work( jw+1 ), lwork-jw, info )
623 *
624 * ==== Update vertical slab in H ====
625 *
626  IF( wantt ) THEN
627  ltop = 1
628  ELSE
629  ltop = ktop
630  END IF
631  DO 70 krow = ltop, kwtop - 1, nv
632  kln = min( nv, kwtop-krow )
633  CALL dgemm( 'N', 'N', kln, jw, jw, one, h( krow, kwtop ),
634  $ ldh, v, ldv, zero, wv, ldwv )
635  CALL dlacpy( 'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
636  70 CONTINUE
637 *
638 * ==== Update horizontal slab in H ====
639 *
640  IF( wantt ) THEN
641  DO 80 kcol = kbot + 1, n, nh
642  kln = min( nh, n-kcol+1 )
643  CALL dgemm( 'C', 'N', jw, kln, jw, one, v, ldv,
644  $ h( kwtop, kcol ), ldh, zero, t, ldt )
645  CALL dlacpy( 'A', jw, kln, t, ldt, h( kwtop, kcol ),
646  $ ldh )
647  80 CONTINUE
648  END IF
649 *
650 * ==== Update vertical slab in Z ====
651 *
652  IF( wantz ) THEN
653  DO 90 krow = iloz, ihiz, nv
654  kln = min( nv, ihiz-krow+1 )
655  CALL dgemm( 'N', 'N', kln, jw, jw, one, z( krow, kwtop ),
656  $ ldz, v, ldv, zero, wv, ldwv )
657  CALL dlacpy( 'A', kln, jw, wv, ldwv, z( krow, kwtop ),
658  $ ldz )
659  90 CONTINUE
660  END IF
661  END IF
662 *
663 * ==== Return the number of deflations ... ====
664 *
665  nd = jw - ns
666 *
667 * ==== ... and the number of shifts. (Subtracting
668 * . INFQR from the spike length takes care
669 * . of the case of a rare QR failure while
670 * . calculating eigenvalues of the deflation
671 * . window.) ====
672 *
673  ns = ns - infqr
674 *
675 * ==== Return optimal workspace. ====
676 *
677  work( 1 ) = dble( lwkopt )
678 *
679 * ==== End of DLAQR2 ====
680 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:103
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.
Definition: dlaset.f:110
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:82
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:187
subroutine dgehrd(N, ILO, IHI, A, LDA, TAU, WORK, LWORK, INFO)
DGEHRD
Definition: dgehrd.f:167
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.
Definition: dlanv2.f:127
subroutine dlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
DLARF applies an elementary reflector to a general rectangular matrix.
Definition: dlarf.f:124
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,...
Definition: dlahqr.f:207
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
Definition: dlarfg.f:106
subroutine dtrexc(COMPQ, N, T, LDT, Q, LDQ, IFST, ILST, WORK, INFO)
DTREXC
Definition: dtrexc.f:148
subroutine dormhr(SIDE, TRANS, M, N, ILO, IHI, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMHR
Definition: dormhr.f:178
Here is the call graph for this function:
Here is the caller graph for this function: