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

◆ zlaqr3()

subroutine zlaqr3 ( 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 
)

ZLAQR3 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 ZLAQR3 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
    Aggressive early deflation:

    ZLAQR3 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; ZLAQR3
          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 264 of file zlaqr3.f.

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