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