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

◆ claqz0()

recursive subroutine claqz0 ( character, intent(in) wants,
character, intent(in) wantq,
character, intent(in) wantz,
integer, intent(in) n,
integer, intent(in) ilo,
integer, intent(in) ihi,
complex, dimension( lda, * ), intent(inout) a,
integer, intent(in) lda,
complex, dimension( ldb, * ), intent(inout) b,
integer, intent(in) ldb,
complex, dimension( * ), intent(inout) alpha,
complex, dimension( * ), intent(inout) beta,
complex, dimension( ldq, * ), intent(inout) q,
integer, intent(in) ldq,
complex, dimension( ldz, * ), intent(inout) z,
integer, intent(in) ldz,
complex, dimension( * ), intent(inout) work,
integer, intent(in) lwork,
real, dimension( * ), intent(out) rwork,
integer, intent(in) rec,
integer, intent(out) info )

CLAQZ0

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

Purpose:
!>
!> CLAQZ0 computes the eigenvalues of a matrix pair (H,T),
!> where H is an upper Hessenberg matrix and T is upper triangular,
!> using the double-shift QZ method.
!> Matrix pairs of this type are produced by the reduction to
!> generalized upper Hessenberg form of a matrix pair (A,B):
!>
!>    A = Q1*H*Z1**H,  B = Q1*T*Z1**H,
!>
!> as computed by CGGHRD.
!>
!> If JOB='S', then the Hessenberg-triangular pair (H,T) is
!> also reduced to generalized Schur form,
!>
!>    H = Q*S*Z**H,  T = Q*P*Z**H,
!>
!> where Q and Z are unitary matrices, P and S are an upper triangular
!> matrices.
!>
!> Optionally, the unitary matrix Q from the generalized Schur
!> factorization may be postmultiplied into an input matrix Q1, and the
!> unitary matrix Z may be postmultiplied into an input matrix Z1.
!> If Q1 and Z1 are the unitary matrices from CGGHRD that reduced
!> the matrix pair (A,B) to generalized upper Hessenberg form, then the
!> output matrices Q1*Q and Z1*Z are the unitary factors from the
!> generalized Schur factorization of (A,B):
!>
!>    A = (Q1*Q)*S*(Z1*Z)**H,  B = (Q1*Q)*P*(Z1*Z)**H.
!>
!> To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently,
!> of (A,B)) are computed as a pair of values (alpha,beta), where alpha is
!> complex and beta real.
!> If beta is nonzero, lambda = alpha / beta is an eigenvalue of the
!> generalized nonsymmetric eigenvalue problem (GNEP)
!>    A*x = lambda*B*x
!> and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the
!> alternate form of the GNEP
!>    mu*A*y = B*y.
!> Eigenvalues can be read directly from the generalized Schur
!> form:
!>   alpha = S(i,i), beta = P(i,i).
!>
!> Ref: C.B. Moler & G.W. Stewart, , SIAM J. Numer. Anal., 10(1973),
!>      pp. 241--256.
!>
!> Ref: B. Kagstrom, D. Kressner, , SIAM J. Numer.
!>      Anal., 29(2006), pp. 199--227.
!>
!> Ref: T. Steel, D. Camps, K. Meerbergen, R. Vandebril 
!> 
Parameters
[in]WANTS
!>          WANTS is CHARACTER*1
!>          = 'E': Compute eigenvalues only;
!>          = 'S': Compute eigenvalues and the Schur form.
!> 
[in]WANTQ
!>          WANTQ is CHARACTER*1
!>          = 'N': Left Schur vectors (Q) are not computed;
!>          = 'I': Q is initialized to the unit matrix and the matrix Q
!>                 of left Schur vectors of (A,B) is returned;
!>          = 'V': Q must contain an unitary matrix Q1 on entry and
!>                 the product Q1*Q is returned.
!> 
[in]WANTZ
!>          WANTZ is CHARACTER*1
!>          = 'N': Right Schur vectors (Z) are not computed;
!>          = 'I': Z is initialized to the unit matrix and the matrix Z
!>                 of right Schur vectors of (A,B) is returned;
!>          = 'V': Z must contain an unitary matrix Z1 on entry and
!>                 the product Z1*Z is returned.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices A, B, Q, and Z.  N >= 0.
!> 
[in]ILO
!>          ILO is INTEGER
!> 
[in]IHI
!>          IHI is INTEGER
!>          ILO and IHI mark the rows and columns of A which are in
!>          Hessenberg form.  It is assumed that A is already upper
!>          triangular in rows and columns 1:ILO-1 and IHI+1:N.
!>          If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA, N)
!>          On entry, the N-by-N upper Hessenberg matrix A.
!>          On exit, if JOB = 'S', A contains the upper triangular
!>          matrix S from the generalized Schur factorization.
!>          If JOB = 'E', the diagonal of A matches that of S, but
!>          the rest of A is unspecified.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max( 1, N ).
!> 
[in,out]B
!>          B is COMPLEX array, dimension (LDB, N)
!>          On entry, the N-by-N upper triangular matrix B.
!>          On exit, if JOB = 'S', B contains the upper triangular
!>          matrix P from the generalized Schur factorization.
!>          If JOB = 'E', the diagonal of B matches that of P, but
!>          the rest of B is unspecified.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max( 1, N ).
!> 
[out]ALPHA
!>          ALPHA is COMPLEX array, dimension (N)
!>          Each scalar alpha defining an eigenvalue
!>          of GNEP.
!> 
[out]BETA
!>          BETA is COMPLEX array, dimension (N)
!>          The scalars beta that define the eigenvalues of GNEP.
!>          Together, the quantities alpha = ALPHA(j) and
!>          beta = BETA(j) represent the j-th eigenvalue of the matrix
!>          pair (A,B), in one of the forms lambda = alpha/beta or
!>          mu = beta/alpha.  Since either lambda or mu may overflow,
!>          they should not, in general, be computed.
!> 
[in,out]Q
!>          Q is COMPLEX array, dimension (LDQ, N)
!>          On entry, if COMPQ = 'V', the unitary matrix Q1 used in
!>          the reduction of (A,B) to generalized Hessenberg form.
!>          On exit, if COMPQ = 'I', the unitary matrix of left Schur
!>          vectors of (A,B), and if COMPQ = 'V', the unitary matrix
!>          of left Schur vectors of (A,B).
!>          Not referenced if COMPQ = 'N'.
!> 
[in]LDQ
!>          LDQ is INTEGER
!>          The leading dimension of the array Q.  LDQ >= 1.
!>          If COMPQ='V' or 'I', then LDQ >= N.
!> 
[in,out]Z
!>          Z is COMPLEX array, dimension (LDZ, N)
!>          On entry, if COMPZ = 'V', the unitary matrix Z1 used in
!>          the reduction of (A,B) to generalized Hessenberg form.
!>          On exit, if COMPZ = 'I', the unitary matrix of
!>          right Schur vectors of (H,T), and if COMPZ = 'V', the
!>          unitary matrix of right Schur vectors of (A,B).
!>          Not referenced if COMPZ = 'N'.
!> 
[in]LDZ
!>          LDZ is INTEGER
!>          The leading dimension of the array Z.  LDZ >= 1.
!>          If COMPZ='V' or 'I', then LDZ >= N.
!> 
[out]WORK
!>          WORK is COMPLEX array, dimension (MAX(1,LWORK))
!>          On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.  LWORK >= max(1,N).
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]RWORK
!>          RWORK is REAL array, dimension (N)
!> 
[in]REC
!>          REC is INTEGER
!>             REC indicates the current recursion level. Should be set
!>             to 0 on first call.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          = 1,...,N: the QZ iteration did not converge.  (A,B) is not
!>                     in Schur form, but ALPHA(i) and
!>                     BETA(i), i=INFO+1,...,N should be correct.
!> 
Author
Thijs Steel, KU Leuven
Date
May 2020

Definition at line 278 of file claqz0.f.

283 IMPLICIT NONE
284
285* Arguments
286 CHARACTER, INTENT( IN ) :: WANTS, WANTQ, WANTZ
287 INTEGER, INTENT( IN ) :: N, ILO, IHI, LDA, LDB, LDQ, LDZ, LWORK,
288 $ REC
289 INTEGER, INTENT( OUT ) :: INFO
290 COMPLEX, INTENT( INOUT ) :: A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
291 $ Z( LDZ, * ), ALPHA( * ), BETA( * ), WORK( * )
292 REAL, INTENT( OUT ) :: RWORK( * )
293
294* Parameters
295 COMPLEX CZERO, CONE
296 parameter( czero = ( 0.0, 0.0 ), cone = ( 1.0, 0.0 ) )
297 REAL :: ZERO, ONE, HALF
298 parameter( zero = 0.0, one = 1.0, half = 0.5 )
299
300* Local scalars
301 REAL :: SMLNUM, ULP, SAFMIN, SAFMAX, C1, TEMPR, BNORM, BTOL
302 COMPLEX :: ESHIFT, S1, TEMP
303 INTEGER :: ISTART, ISTOP, IITER, MAXIT, ISTART2, K, LD, NSHIFTS,
304 $ NBLOCK, NW, NMIN, NIBBLE, N_UNDEFLATED, N_DEFLATED,
305 $ NS, SWEEP_INFO, SHIFTPOS, LWORKREQ, K2, ISTARTM,
306 $ ISTOPM, IWANTS, IWANTQ, IWANTZ, NORM_INFO, AED_INFO,
307 $ NWR, NBR, NSR, ITEMP1, ITEMP2, RCOST
308 LOGICAL :: ILSCHUR, ILQ, ILZ
309 CHARACTER :: JBCMPZ*3
310
311* External Functions
312 EXTERNAL :: xerbla, chgeqz, claqz2, claqz3, claset,
313 $ clartg, crot
314 REAL, EXTERNAL :: SLAMCH, CLANHS
315 LOGICAL, EXTERNAL :: LSAME
316 INTEGER, EXTERNAL :: ILAENV
317
318*
319* Decode wantS,wantQ,wantZ
320*
321 IF( lsame( wants, 'E' ) ) THEN
322 ilschur = .false.
323 iwants = 1
324 ELSE IF( lsame( wants, 'S' ) ) THEN
325 ilschur = .true.
326 iwants = 2
327 ELSE
328 iwants = 0
329 END IF
330
331 IF( lsame( wantq, 'N' ) ) THEN
332 ilq = .false.
333 iwantq = 1
334 ELSE IF( lsame( wantq, 'V' ) ) THEN
335 ilq = .true.
336 iwantq = 2
337 ELSE IF( lsame( wantq, 'I' ) ) THEN
338 ilq = .true.
339 iwantq = 3
340 ELSE
341 iwantq = 0
342 END IF
343
344 IF( lsame( wantz, 'N' ) ) THEN
345 ilz = .false.
346 iwantz = 1
347 ELSE IF( lsame( wantz, 'V' ) ) THEN
348 ilz = .true.
349 iwantz = 2
350 ELSE IF( lsame( wantz, 'I' ) ) THEN
351 ilz = .true.
352 iwantz = 3
353 ELSE
354 iwantz = 0
355 END IF
356*
357* Check Argument Values
358*
359 info = 0
360 IF( iwants.EQ.0 ) THEN
361 info = -1
362 ELSE IF( iwantq.EQ.0 ) THEN
363 info = -2
364 ELSE IF( iwantz.EQ.0 ) THEN
365 info = -3
366 ELSE IF( n.LT.0 ) THEN
367 info = -4
368 ELSE IF( ilo.LT.1 ) THEN
369 info = -5
370 ELSE IF( ihi.GT.n .OR. ihi.LT.ilo-1 ) THEN
371 info = -6
372 ELSE IF( lda.LT.n ) THEN
373 info = -8
374 ELSE IF( ldb.LT.n ) THEN
375 info = -10
376 ELSE IF( ldq.LT.1 .OR. ( ilq .AND. ldq.LT.n ) ) THEN
377 info = -15
378 ELSE IF( ldz.LT.1 .OR. ( ilz .AND. ldz.LT.n ) ) THEN
379 info = -17
380 END IF
381 IF( info.NE.0 ) THEN
382 CALL xerbla( 'CLAQZ0', -info )
383 RETURN
384 END IF
385
386*
387* Quick return if possible
388*
389 IF( n.LE.0 ) THEN
390 work( 1 ) = real( 1 )
391 RETURN
392 END IF
393
394*
395* Get the parameters
396*
397 jbcmpz( 1:1 ) = wants
398 jbcmpz( 2:2 ) = wantq
399 jbcmpz( 3:3 ) = wantz
400
401 nmin = ilaenv( 12, 'CLAQZ0', jbcmpz, n, ilo, ihi, lwork )
402
403 nwr = ilaenv( 13, 'CLAQZ0', jbcmpz, n, ilo, ihi, lwork )
404 nwr = max( 2, nwr )
405 nwr = min( ihi-ilo+1, ( n-1 ) / 3, nwr )
406
407 nibble = ilaenv( 14, 'CLAQZ0', jbcmpz, n, ilo, ihi, lwork )
408
409 nsr = ilaenv( 15, 'CLAQZ0', jbcmpz, n, ilo, ihi, lwork )
410 nsr = min( nsr, ( n+6 ) / 9, ihi-ilo )
411 nsr = max( 2, nsr-mod( nsr, 2 ) )
412
413 rcost = ilaenv( 17, 'CLAQZ0', jbcmpz, n, ilo, ihi, lwork )
414 itemp1 = int( real( nsr )/sqrt( 1+2*real( nsr )/
415 $ ( real( rcost )/100*real( n ) ) ) )
416 itemp1 = ( ( itemp1-1 )/4 )*4+4
417 nbr = nsr+itemp1
418
419 IF( n .LT. nmin .OR. rec .GE. 2 ) THEN
420 CALL chgeqz( wants, wantq, wantz, n, ilo, ihi, a, lda, b,
421 $ ldb,
422 $ alpha, beta, q, ldq, z, ldz, work, lwork, rwork,
423 $ info )
424 RETURN
425 END IF
426
427*
428* Find out required workspace
429*
430
431* Workspace query to CLAQZ2
432 nw = max( nwr, nmin )
433 CALL claqz2( ilschur, ilq, ilz, n, ilo, ihi, nw, a, lda, b,
434 $ ldb,
435 $ q, ldq, z, ldz, n_undeflated, n_deflated, alpha,
436 $ beta, work, nw, work, nw, work, -1, rwork, rec,
437 $ aed_info )
438 itemp1 = int( work( 1 ) )
439* Workspace query to CLAQZ3
440 CALL claqz3( ilschur, ilq, ilz, n, ilo, ihi, nsr, nbr, alpha,
441 $ beta, a, lda, b, ldb, q, ldq, z, ldz, work, nbr,
442 $ work, nbr, work, -1, sweep_info )
443 itemp2 = int( work( 1 ) )
444
445 lworkreq = max( itemp1+2*nw**2, itemp2+2*nbr**2 )
446 IF ( lwork .EQ.-1 ) THEN
447 work( 1 ) = real( lworkreq )
448 RETURN
449 ELSE IF ( lwork .LT. lworkreq ) THEN
450 info = -19
451 END IF
452 IF( info.NE.0 ) THEN
453 CALL xerbla( 'CLAQZ0', info )
454 RETURN
455 END IF
456*
457* Initialize Q and Z
458*
459 IF( iwantq.EQ.3 ) CALL claset( 'FULL', n, n, czero, cone, q,
460 $ ldq )
461 IF( iwantz.EQ.3 ) CALL claset( 'FULL', n, n, czero, cone, z,
462 $ ldz )
463
464* Get machine constants
465 safmin = slamch( 'SAFE MINIMUM' )
466 safmax = one/safmin
467 ulp = slamch( 'PRECISION' )
468 smlnum = safmin*( real( n )/ulp )
469
470 bnorm = clanhs( 'F', ihi-ilo+1, b( ilo, ilo ), ldb, rwork )
471 btol = max( safmin, ulp*bnorm )
472
473 istart = ilo
474 istop = ihi
475 maxit = 30*( ihi-ilo+1 )
476 ld = 0
477
478 DO iiter = 1, maxit
479 IF( iiter .GE. maxit ) THEN
480 info = istop+1
481 GOTO 80
482 END IF
483 IF ( istart+1 .GE. istop ) THEN
484 istop = istart
485 EXIT
486 END IF
487
488* Check deflations at the end
489 IF ( abs( a( istop, istop-1 ) ) .LE. max( smlnum,
490 $ ulp*( abs( a( istop, istop ) )+abs( a( istop-1,
491 $ istop-1 ) ) ) ) ) THEN
492 a( istop, istop-1 ) = czero
493 istop = istop-1
494 ld = 0
495 eshift = czero
496 END IF
497* Check deflations at the start
498 IF ( abs( a( istart+1, istart ) ) .LE. max( smlnum,
499 $ ulp*( abs( a( istart, istart ) )+abs( a( istart+1,
500 $ istart+1 ) ) ) ) ) THEN
501 a( istart+1, istart ) = czero
502 istart = istart+1
503 ld = 0
504 eshift = czero
505 END IF
506
507 IF ( istart+1 .GE. istop ) THEN
508 EXIT
509 END IF
510
511* Check interior deflations
512 istart2 = istart
513 DO k = istop, istart+1, -1
514 IF ( abs( a( k, k-1 ) ) .LE. max( smlnum, ulp*( abs( a( k,
515 $ k ) )+abs( a( k-1, k-1 ) ) ) ) ) THEN
516 a( k, k-1 ) = czero
517 istart2 = k
518 EXIT
519 END IF
520 END DO
521
522* Get range to apply rotations to
523 IF ( ilschur ) THEN
524 istartm = 1
525 istopm = n
526 ELSE
527 istartm = istart2
528 istopm = istop
529 END IF
530
531* Check infinite eigenvalues, this is done without blocking so might
532* slow down the method when many infinite eigenvalues are present
533 k = istop
534 DO WHILE ( k.GE.istart2 )
535
536 IF( abs( b( k, k ) ) .LT. btol ) THEN
537* A diagonal element of B is negligible, move it
538* to the top and deflate it
539
540 DO k2 = k, istart2+1, -1
541 CALL clartg( b( k2-1, k2 ), b( k2-1, k2-1 ), c1,
542 $ s1,
543 $ temp )
544 b( k2-1, k2 ) = temp
545 b( k2-1, k2-1 ) = czero
546
547 CALL crot( k2-2-istartm+1, b( istartm, k2 ), 1,
548 $ b( istartm, k2-1 ), 1, c1, s1 )
549 CALL crot( min( k2+1, istop )-istartm+1,
550 $ a( istartm,
551 $ k2 ), 1, a( istartm, k2-1 ), 1, c1, s1 )
552 IF ( ilz ) THEN
553 CALL crot( n, z( 1, k2 ), 1, z( 1, k2-1 ), 1,
554 $ c1,
555 $ s1 )
556 END IF
557
558 IF( k2.LT.istop ) THEN
559 CALL clartg( a( k2, k2-1 ), a( k2+1, k2-1 ), c1,
560 $ s1, temp )
561 a( k2, k2-1 ) = temp
562 a( k2+1, k2-1 ) = czero
563
564 CALL crot( istopm-k2+1, a( k2, k2 ), lda,
565 $ a( k2+1,
566 $ k2 ), lda, c1, s1 )
567 CALL crot( istopm-k2+1, b( k2, k2 ), ldb,
568 $ b( k2+1,
569 $ k2 ), ldb, c1, s1 )
570 IF( ilq ) THEN
571 CALL crot( n, q( 1, k2 ), 1, q( 1, k2+1 ), 1,
572 $ c1, conjg( s1 ) )
573 END IF
574 END IF
575
576 END DO
577
578 IF( istart2.LT.istop )THEN
579 CALL clartg( a( istart2, istart2 ), a( istart2+1,
580 $ istart2 ), c1, s1, temp )
581 a( istart2, istart2 ) = temp
582 a( istart2+1, istart2 ) = czero
583
584 CALL crot( istopm-( istart2+1 )+1, a( istart2,
585 $ istart2+1 ), lda, a( istart2+1,
586 $ istart2+1 ), lda, c1, s1 )
587 CALL crot( istopm-( istart2+1 )+1, b( istart2,
588 $ istart2+1 ), ldb, b( istart2+1,
589 $ istart2+1 ), ldb, c1, s1 )
590 IF( ilq ) THEN
591 CALL crot( n, q( 1, istart2 ), 1, q( 1,
592 $ istart2+1 ), 1, c1, conjg( s1 ) )
593 END IF
594 END IF
595
596 istart2 = istart2+1
597
598 END IF
599 k = k-1
600 END DO
601
602* istart2 now points to the top of the bottom right
603* unreduced Hessenberg block
604 IF ( istart2 .GE. istop ) THEN
605 istop = istart2-1
606 ld = 0
607 eshift = czero
608 cycle
609 END IF
610
611 nw = nwr
612 nshifts = nsr
613 nblock = nbr
614
615 IF ( istop-istart2+1 .LT. nmin ) THEN
616* Setting nw to the size of the subblock will make AED deflate
617* all the eigenvalues. This is slightly more efficient than just
618* using CHGEQZ because the off diagonal part gets updated via BLAS.
619 IF ( istop-istart+1 .LT. nmin ) THEN
620 nw = istop-istart+1
621 istart2 = istart
622 ELSE
623 nw = istop-istart2+1
624 END IF
625 END IF
626
627*
628* Time for AED
629*
630 CALL claqz2( ilschur, ilq, ilz, n, istart2, istop, nw, a,
631 $ lda,
632 $ b, ldb, q, ldq, z, ldz, n_undeflated, n_deflated,
633 $ alpha, beta, work, nw, work( nw**2+1 ), nw,
634 $ work( 2*nw**2+1 ), lwork-2*nw**2, rwork, rec,
635 $ aed_info )
636
637 IF ( n_deflated > 0 ) THEN
638 istop = istop-n_deflated
639 ld = 0
640 eshift = czero
641 END IF
642
643 IF ( 100*n_deflated > nibble*( n_deflated+n_undeflated ) .OR.
644 $ istop-istart2+1 .LT. nmin ) THEN
645* AED has uncovered many eigenvalues. Skip a QZ sweep and run
646* AED again.
647 cycle
648 END IF
649
650 ld = ld+1
651
652 ns = min( nshifts, istop-istart2 )
653 ns = min( ns, n_undeflated )
654 shiftpos = istop-n_undeflated+1
655
656 IF ( mod( ld, 6 ) .EQ. 0 ) THEN
657*
658* Exceptional shift. Chosen for no particularly good reason.
659*
660 IF( ( real( maxit )*safmin )*abs( a( istop,
661 $ istop-1 ) ).LT.abs( a( istop-1, istop-1 ) ) ) THEN
662 eshift = a( istop, istop-1 )/b( istop-1, istop-1 )
663 ELSE
664 eshift = eshift+cone/( safmin*real( maxit ) )
665 END IF
666 alpha( shiftpos ) = cone
667 beta( shiftpos ) = eshift
668 ns = 1
669 END IF
670
671*
672* Time for a QZ sweep
673*
674 CALL claqz3( ilschur, ilq, ilz, n, istart2, istop, ns,
675 $ nblock,
676 $ alpha( shiftpos ), beta( shiftpos ), a, lda, b,
677 $ ldb, q, ldq, z, ldz, work, nblock, work( nblock**
678 $ 2+1 ), nblock, work( 2*nblock**2+1 ),
679 $ lwork-2*nblock**2, sweep_info )
680
681 END DO
682
683*
684* Call CHGEQZ to normalize the eigenvalue blocks and set the eigenvalues
685* If all the eigenvalues have been found, CHGEQZ will not do any iterations
686* and only normalize the blocks. In case of a rare convergence failure,
687* the single shift might perform better.
688*
689 80 CALL chgeqz( wants, wantq, wantz, n, ilo, ihi, a, lda, b, ldb,
690 $ alpha, beta, q, ldq, z, ldz, work, lwork, rwork,
691 $ norm_info )
692
693 info = norm_info
694
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
CHGEQZ
Definition chgeqz.f:283
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clanhs(norm, n, a, lda, work)
CLANHS returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition clanhs.f:107
recursive subroutine claqz2(ilschur, ilq, ilz, n, ilo, ihi, nw, a, lda, b, ldb, q, ldq, z, ldz, ns, nd, alpha, beta, qc, ldqc, zc, ldzc, work, lwork, rwork, rec, info)
CLAQZ2
Definition claqz2.f:234
subroutine claqz3(ilschur, ilq, ilz, n, ilo, ihi, nshifts, nblock_desired, alpha, beta, a, lda, b, ldb, q, ldq, z, ldz, qc, ldqc, zc, ldzc, work, lwork, info)
CLAQZ3
Definition claqz3.f:205
subroutine clartg(f, g, c, s, r)
CLARTG generates a plane rotation with real cosine and complex sine.
Definition clartg.f90:116
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:104
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine crot(n, cx, incx, cy, incy, c, s)
CROT applies a plane rotation with real cosine and complex sine to a pair of complex vectors.
Definition crot.f:101
Here is the call graph for this function:
Here is the caller graph for this function: