LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zlaqr2()

subroutine zlaqr2 ( logical  wantt,
logical  wantz,
integer  n,
integer  ktop,
integer  kbot,
integer  nw,
complex*16, dimension( ldh, * )  h,
integer  ldh,
integer  iloz,
integer  ihiz,
complex*16, dimension( ldz, * )  z,
integer  ldz,
integer  ns,
integer  nd,
complex*16, dimension( * )  sh,
complex*16, dimension( ldv, * )  v,
integer  ldv,
integer  nh,
complex*16, dimension( ldt, * )  t,
integer  ldt,
integer  nv,
complex*16, dimension( ldwv, * )  wv,
integer  ldwv,
complex*16, dimension( * )  work,
integer  lwork 
)

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

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

Purpose:
    ZLAQR2 is identical to ZLAQR3 except that it avoids
    recursion by calling ZLAHQR instead of ZLAQR4.

    Aggressive early deflation:

    ZLAQR2 accepts as input an upper Hessenberg matrix
    H and performs an unitary 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 unitary 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 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 unitary matrix Z is updated so
          so that the unitary 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 unitary 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 COMPLEX*16 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 a unitary
          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 COMPLEX*16 array, dimension (LDZ,N)
          IF WANTZ is .TRUE., then on output, the unitary
          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]SH
          SH is COMPLEX*16 array, dimension (KBOT)
          On output, approximate eigenvalues that may
          be used for shifts are stored in SH(KBOT-ND-NS+1)
          through SR(KBOT-ND).  Converged eigenvalues are
          stored in SH(KBOT-ND+1) through SH(KBOT).
[out]V
          V is COMPLEX*16 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 COMPLEX*16 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 COMPLEX*16 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 COMPLEX*16 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; ZLAQR2
          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 267 of file zlaqr2.f.

270*
271* -- LAPACK auxiliary routine --
272* -- LAPACK is a software package provided by Univ. of Tennessee, --
273* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
274*
275* .. Scalar Arguments ..
276 INTEGER IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV,
277 $ LDZ, LWORK, N, ND, NH, NS, NV, NW
278 LOGICAL WANTT, WANTZ
279* ..
280* .. Array Arguments ..
281 COMPLEX*16 H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ),
282 $ WORK( * ), WV( LDWV, * ), Z( LDZ, * )
283* ..
284*
285* ================================================================
286*
287* .. Parameters ..
288 COMPLEX*16 ZERO, ONE
289 parameter( zero = ( 0.0d0, 0.0d0 ),
290 $ one = ( 1.0d0, 0.0d0 ) )
291 DOUBLE PRECISION RZERO, RONE
292 parameter( rzero = 0.0d0, rone = 1.0d0 )
293* ..
294* .. Local Scalars ..
295 COMPLEX*16 BETA, CDUM, S, TAU
296 DOUBLE PRECISION FOO, SAFMAX, SAFMIN, SMLNUM, ULP
297 INTEGER I, IFST, ILST, INFO, INFQR, J, JW, KCOL, KLN,
298 $ KNT, KROW, KWTOP, LTOP, LWK1, LWK2, LWKOPT
299* ..
300* .. External Functions ..
301 DOUBLE PRECISION DLAMCH
302 EXTERNAL dlamch
303* ..
304* .. External Subroutines ..
305 EXTERNAL zcopy, zgehrd, zgemm, zlacpy, zlahqr,
307* ..
308* .. Intrinsic Functions ..
309 INTRINSIC abs, dble, dcmplx, dconjg, dimag, int, max, min
310* ..
311* .. Statement Functions ..
312 DOUBLE PRECISION CABS1
313* ..
314* .. Statement Function definitions ..
315 cabs1( cdum ) = abs( dble( cdum ) ) + abs( dimag( cdum ) )
316* ..
317* .. Executable Statements ..
318*
319* ==== Estimate optimal workspace. ====
320*
321 jw = min( nw, kbot-ktop+1 )
322 IF( jw.LE.2 ) THEN
323 lwkopt = 1
324 ELSE
325*
326* ==== Workspace query call to ZGEHRD ====
327*
328 CALL zgehrd( jw, 1, jw-1, t, ldt, work, work, -1, info )
329 lwk1 = int( work( 1 ) )
330*
331* ==== Workspace query call to ZUNMHR ====
332*
333 CALL zunmhr( 'R', 'N', jw, jw, 1, jw-1, t, ldt, work, v, ldv,
334 $ work, -1, info )
335 lwk2 = int( work( 1 ) )
336*
337* ==== Optimal workspace ====
338*
339 lwkopt = jw + max( lwk1, lwk2 )
340 END IF
341*
342* ==== Quick return in case of workspace query. ====
343*
344 IF( lwork.EQ.-1 ) THEN
345 work( 1 ) = dcmplx( lwkopt, 0 )
346 RETURN
347 END IF
348*
349* ==== Nothing to do ...
350* ... for an empty active block ... ====
351 ns = 0
352 nd = 0
353 work( 1 ) = one
354 IF( ktop.GT.kbot )
355 $ RETURN
356* ... nor for an empty deflation window. ====
357 IF( nw.LT.1 )
358 $ RETURN
359*
360* ==== Machine constants ====
361*
362 safmin = dlamch( 'SAFE MINIMUM' )
363 safmax = rone / safmin
364 ulp = dlamch( 'PRECISION' )
365 smlnum = safmin*( dble( n ) / ulp )
366*
367* ==== Setup deflation window ====
368*
369 jw = min( nw, kbot-ktop+1 )
370 kwtop = kbot - jw + 1
371 IF( kwtop.EQ.ktop ) THEN
372 s = zero
373 ELSE
374 s = h( kwtop, kwtop-1 )
375 END IF
376*
377 IF( kbot.EQ.kwtop ) THEN
378*
379* ==== 1-by-1 deflation window: not much to do ====
380*
381 sh( kwtop ) = h( kwtop, kwtop )
382 ns = 1
383 nd = 0
384 IF( cabs1( s ).LE.max( smlnum, ulp*cabs1( h( kwtop,
385 $ kwtop ) ) ) ) THEN
386 ns = 0
387 nd = 1
388 IF( kwtop.GT.ktop )
389 $ h( kwtop, kwtop-1 ) = zero
390 END IF
391 work( 1 ) = one
392 RETURN
393 END IF
394*
395* ==== Convert to spike-triangular form. (In case of a
396* . rare QR failure, this routine continues to do
397* . aggressive early deflation using that part of
398* . the deflation window that converged using INFQR
399* . here and there to keep track.) ====
400*
401 CALL zlacpy( 'U', jw, jw, h( kwtop, kwtop ), ldh, t, ldt )
402 CALL zcopy( jw-1, h( kwtop+1, kwtop ), ldh+1, t( 2, 1 ), ldt+1 )
403*
404 CALL zlaset( 'A', jw, jw, zero, one, v, ldv )
405 CALL zlahqr( .true., .true., jw, 1, jw, t, ldt, sh( kwtop ), 1,
406 $ jw, v, ldv, infqr )
407*
408* ==== Deflation detection loop ====
409*
410 ns = jw
411 ilst = infqr + 1
412 DO 10 knt = infqr + 1, jw
413*
414* ==== Small spike tip deflation test ====
415*
416 foo = cabs1( t( ns, ns ) )
417 IF( foo.EQ.rzero )
418 $ foo = cabs1( s )
419 IF( cabs1( s )*cabs1( v( 1, ns ) ).LE.max( smlnum, ulp*foo ) )
420 $ THEN
421*
422* ==== One more converged eigenvalue ====
423*
424 ns = ns - 1
425 ELSE
426*
427* ==== One undeflatable eigenvalue. Move it up out of the
428* . way. (ZTREXC can not fail in this case.) ====
429*
430 ifst = ns
431 CALL ztrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, info )
432 ilst = ilst + 1
433 END IF
434 10 CONTINUE
435*
436* ==== Return to Hessenberg form ====
437*
438 IF( ns.EQ.0 )
439 $ s = zero
440*
441 IF( ns.LT.jw ) THEN
442*
443* ==== sorting the diagonal of T improves accuracy for
444* . graded matrices. ====
445*
446 DO 30 i = infqr + 1, ns
447 ifst = i
448 DO 20 j = i + 1, ns
449 IF( cabs1( t( j, j ) ).GT.cabs1( t( ifst, ifst ) ) )
450 $ ifst = j
451 20 CONTINUE
452 ilst = i
453 IF( ifst.NE.ilst )
454 $ CALL ztrexc( 'V', jw, t, ldt, v, ldv, ifst, ilst, info )
455 30 CONTINUE
456 END IF
457*
458* ==== Restore shift/eigenvalue array from T ====
459*
460 DO 40 i = infqr + 1, jw
461 sh( kwtop+i-1 ) = t( i, i )
462 40 CONTINUE
463*
464*
465 IF( ns.LT.jw .OR. s.EQ.zero ) THEN
466 IF( ns.GT.1 .AND. s.NE.zero ) THEN
467*
468* ==== Reflect spike back into lower triangle ====
469*
470 CALL zcopy( ns, v, ldv, work, 1 )
471 DO 50 i = 1, ns
472 work( i ) = dconjg( work( i ) )
473 50 CONTINUE
474 beta = work( 1 )
475 CALL zlarfg( ns, beta, work( 2 ), 1, tau )
476 work( 1 ) = one
477*
478 CALL zlaset( 'L', jw-2, jw-2, zero, zero, t( 3, 1 ), ldt )
479*
480 CALL zlarf( 'L', ns, jw, work, 1, dconjg( tau ), t, ldt,
481 $ work( jw+1 ) )
482 CALL zlarf( 'R', ns, ns, work, 1, tau, t, ldt,
483 $ work( jw+1 ) )
484 CALL zlarf( 'R', jw, ns, work, 1, tau, v, ldv,
485 $ work( jw+1 ) )
486*
487 CALL zgehrd( jw, 1, ns, t, ldt, work, work( jw+1 ),
488 $ lwork-jw, info )
489 END IF
490*
491* ==== Copy updated reduced window into place ====
492*
493 IF( kwtop.GT.1 )
494 $ h( kwtop, kwtop-1 ) = s*dconjg( v( 1, 1 ) )
495 CALL zlacpy( 'U', jw, jw, t, ldt, h( kwtop, kwtop ), ldh )
496 CALL zcopy( jw-1, t( 2, 1 ), ldt+1, h( kwtop+1, kwtop ),
497 $ ldh+1 )
498*
499* ==== Accumulate orthogonal matrix in order update
500* . H and Z, if requested. ====
501*
502 IF( ns.GT.1 .AND. s.NE.zero )
503 $ CALL zunmhr( 'R', 'N', jw, ns, 1, ns, t, ldt, work, v, ldv,
504 $ work( jw+1 ), lwork-jw, info )
505*
506* ==== Update vertical slab in H ====
507*
508 IF( wantt ) THEN
509 ltop = 1
510 ELSE
511 ltop = ktop
512 END IF
513 DO 60 krow = ltop, kwtop - 1, nv
514 kln = min( nv, kwtop-krow )
515 CALL zgemm( 'N', 'N', kln, jw, jw, one, h( krow, kwtop ),
516 $ ldh, v, ldv, zero, wv, ldwv )
517 CALL zlacpy( 'A', kln, jw, wv, ldwv, h( krow, kwtop ), ldh )
518 60 CONTINUE
519*
520* ==== Update horizontal slab in H ====
521*
522 IF( wantt ) THEN
523 DO 70 kcol = kbot + 1, n, nh
524 kln = min( nh, n-kcol+1 )
525 CALL zgemm( 'C', 'N', jw, kln, jw, one, v, ldv,
526 $ h( kwtop, kcol ), ldh, zero, t, ldt )
527 CALL zlacpy( 'A', jw, kln, t, ldt, h( kwtop, kcol ),
528 $ ldh )
529 70 CONTINUE
530 END IF
531*
532* ==== Update vertical slab in Z ====
533*
534 IF( wantz ) THEN
535 DO 80 krow = iloz, ihiz, nv
536 kln = min( nv, ihiz-krow+1 )
537 CALL zgemm( 'N', 'N', kln, jw, jw, one, z( krow, kwtop ),
538 $ ldz, v, ldv, zero, wv, ldwv )
539 CALL zlacpy( 'A', kln, jw, wv, ldwv, z( krow, kwtop ),
540 $ ldz )
541 80 CONTINUE
542 END IF
543 END IF
544*
545* ==== Return the number of deflations ... ====
546*
547 nd = jw - ns
548*
549* ==== ... and the number of shifts. (Subtracting
550* . INFQR from the spike length takes care
551* . of the case of a rare QR failure while
552* . calculating eigenvalues of the deflation
553* . window.) ====
554*
555 ns = ns - infqr
556*
557* ==== Return optimal workspace. ====
558*
559 work( 1 ) = dcmplx( lwkopt, 0 )
560*
561* ==== End of ZLAQR2 ====
562*
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZGEHRD
Definition zgehrd.f:167
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlahqr(wantt, wantz, n, ilo, ihi, h, ldh, w, iloz, ihiz, z, ldz, info)
ZLAHQR computes the eigenvalues and Schur factorization of an upper Hessenberg matrix,...
Definition zlahqr.f:195
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine zlarf(side, m, n, v, incv, tau, c, ldc, work)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition zlarf.f:128
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
Definition zlarfg.f:106
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine ztrexc(compq, n, t, ldt, q, ldq, ifst, ilst, info)
ZTREXC
Definition ztrexc.f:126
subroutine zunmhr(side, trans, m, n, ilo, ihi, a, lda, tau, c, ldc, work, lwork, info)
ZUNMHR
Definition zunmhr.f:178
Here is the call graph for this function:
Here is the caller graph for this function: