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

◆ zlaqz0()

recursive subroutine zlaqz0 ( 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*16, dimension( lda, * ), intent(inout)  A,
integer, intent(in)  LDA,
complex*16, dimension( ldb, * ), intent(inout)  B,
integer, intent(in)  LDB,
complex*16, dimension( * ), intent(inout)  ALPHA,
complex*16, dimension( * ), intent(inout)  BETA,
complex*16, dimension( ldq, * ), intent(inout)  Q,
integer, intent(in)  LDQ,
complex*16, dimension( ldz, * ), intent(inout)  Z,
integer, intent(in)  LDZ,
complex*16, dimension( * ), intent(inout)  WORK,
integer, intent(in)  LWORK,
double precision, dimension( * ), intent(out)  RWORK,
integer, intent(in)  REC,
integer, intent(out)  INFO 
)

ZLAQZ0

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

Purpose:
 ZLAQZ0 computes the eigenvalues of a real 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 real matrix pair (A,B):

    A = Q1*H*Z1**H,  B = Q1*T*Z1**H,

 as computed by ZGGHRD.

 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 ZGGHRD 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, "An Algorithm for Generalized Matrix
      Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973),
      pp. 241--256.

 Ref: B. Kagstrom, D. Kressner, "Multishift Variants of the QZ
      Algorithm with Aggressive Early Deflation", SIAM J. Numer.
      Anal., 29(2006), pp. 199--227.

 Ref: T. Steel, D. Camps, K. Meerbergen, R. Vandebril "A multishift,
      multipole rational QZ method with agressive early deflation"
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*16 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 blocks of A match those 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*16 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 blocks of B match those 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*16 array, dimension (N)
          Each scalar alpha defining an eigenvalue
          of GNEP.
[out]BETA
          BETA is COMPLEX*16 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*16 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*16 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*16 array, dimension (MAX(1,LWORK))
          On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[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.
[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 280 of file zlaqz0.f.

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