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

◆ stgsna()

subroutine stgsna ( character  job,
character  howmny,
logical, dimension( * )  select,
integer  n,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( ldb, * )  b,
integer  ldb,
real, dimension( ldvl, * )  vl,
integer  ldvl,
real, dimension( ldvr, * )  vr,
integer  ldvr,
real, dimension( * )  s,
real, dimension( * )  dif,
integer  mm,
integer  m,
real, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  info 
)

STGSNA

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

Purpose:
 STGSNA estimates reciprocal condition numbers for specified
 eigenvalues and/or eigenvectors of a matrix pair (A, B) in
 generalized real Schur canonical form (or of any matrix pair
 (Q*A*Z**T, Q*B*Z**T) with orthogonal matrices Q and Z, where
 Z**T denotes the transpose of Z.

 (A, B) must be in generalized real Schur form (as returned by SGGES),
 i.e. A is block upper triangular with 1-by-1 and 2-by-2 diagonal
 blocks. B is upper triangular.
Parameters
[in]JOB
          JOB is CHARACTER*1
          Specifies whether condition numbers are required for
          eigenvalues (S) or eigenvectors (DIF):
          = 'E': for eigenvalues only (S);
          = 'V': for eigenvectors only (DIF);
          = 'B': for both eigenvalues and eigenvectors (S and DIF).
[in]HOWMNY
          HOWMNY is CHARACTER*1
          = 'A': compute condition numbers for all eigenpairs;
          = 'S': compute condition numbers for selected eigenpairs
                 specified by the array SELECT.
[in]SELECT
          SELECT is LOGICAL array, dimension (N)
          If HOWMNY = 'S', SELECT specifies the eigenpairs for which
          condition numbers are required. To select condition numbers
          for the eigenpair corresponding to a real eigenvalue w(j),
          SELECT(j) must be set to .TRUE.. To select condition numbers
          corresponding to a complex conjugate pair of eigenvalues w(j)
          and w(j+1), either SELECT(j) or SELECT(j+1) or both, must be
          set to .TRUE..
          If HOWMNY = 'A', SELECT is not referenced.
[in]N
          N is INTEGER
          The order of the square matrix pair (A, B). N >= 0.
[in]A
          A is REAL array, dimension (LDA,N)
          The upper quasi-triangular matrix A in the pair (A,B).
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A. LDA >= max(1,N).
[in]B
          B is REAL array, dimension (LDB,N)
          The upper triangular matrix B in the pair (A,B).
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B. LDB >= max(1,N).
[in]VL
          VL is REAL array, dimension (LDVL,M)
          If JOB = 'E' or 'B', VL must contain left eigenvectors of
          (A, B), corresponding to the eigenpairs specified by HOWMNY
          and SELECT. The eigenvectors must be stored in consecutive
          columns of VL, as returned by STGEVC.
          If JOB = 'V', VL is not referenced.
[in]LDVL
          LDVL is INTEGER
          The leading dimension of the array VL. LDVL >= 1.
          If JOB = 'E' or 'B', LDVL >= N.
[in]VR
          VR is REAL array, dimension (LDVR,M)
          If JOB = 'E' or 'B', VR must contain right eigenvectors of
          (A, B), corresponding to the eigenpairs specified by HOWMNY
          and SELECT. The eigenvectors must be stored in consecutive
          columns ov VR, as returned by STGEVC.
          If JOB = 'V', VR is not referenced.
[in]LDVR
          LDVR is INTEGER
          The leading dimension of the array VR. LDVR >= 1.
          If JOB = 'E' or 'B', LDVR >= N.
[out]S
          S is REAL array, dimension (MM)
          If JOB = 'E' or 'B', the reciprocal condition numbers of the
          selected eigenvalues, stored in consecutive elements of the
          array. For a complex conjugate pair of eigenvalues two
          consecutive elements of S are set to the same value. Thus
          S(j), DIF(j), and the j-th columns of VL and VR all
          correspond to the same eigenpair (but not in general the
          j-th eigenpair, unless all eigenpairs are selected).
          If JOB = 'V', S is not referenced.
[out]DIF
          DIF is REAL array, dimension (MM)
          If JOB = 'V' or 'B', the estimated reciprocal condition
          numbers of the selected eigenvectors, stored in consecutive
          elements of the array. For a complex eigenvector two
          consecutive elements of DIF are set to the same value. If
          the eigenvalues cannot be reordered to compute DIF(j), DIF(j)
          is set to 0; this can only occur when the true value would be
          very small anyway.
          If JOB = 'E', DIF is not referenced.
[in]MM
          MM is INTEGER
          The number of elements in the arrays S and DIF. MM >= M.
[out]M
          M is INTEGER
          The number of elements of the arrays S and DIF used to store
          the specified condition numbers; for each selected real
          eigenvalue one element is used, and for each selected complex
          conjugate pair of eigenvalues, two elements are used.
          If HOWMNY = 'A', M is set to N.
[out]WORK
          WORK is REAL 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 JOB = 'V' or 'B' LWORK >= 2*N*(N+2)+16.

          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]IWORK
          IWORK is INTEGER array, dimension (N + 6)
          If JOB = 'E', IWORK is not referenced.
[out]INFO
          INFO is INTEGER
          =0: Successful exit
          <0: If INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  The reciprocal of the condition number of a generalized eigenvalue
  w = (a, b) is defined as

       S(w) = (|u**TAv|**2 + |u**TBv|**2)**(1/2) / (norm(u)*norm(v))

  where u and v are the left and right eigenvectors of (A, B)
  corresponding to w; |z| denotes the absolute value of the complex
  number, and norm(u) denotes the 2-norm of the vector u.
  The pair (a, b) corresponds to an eigenvalue w = a/b (= u**TAv/u**TBv)
  of the matrix pair (A, B). If both a and b equal zero, then (A B) is
  singular and S(I) = -1 is returned.

  An approximate error bound on the chordal distance between the i-th
  computed generalized eigenvalue w and the corresponding exact
  eigenvalue lambda is

       chord(w, lambda) <= EPS * norm(A, B) / S(I)

  where EPS is the machine precision.

  The reciprocal of the condition number DIF(i) of right eigenvector u
  and left eigenvector v corresponding to the generalized eigenvalue w
  is defined as follows:

  a) If the i-th eigenvalue w = (a,b) is real

     Suppose U and V are orthogonal transformations such that

              U**T*(A, B)*V  = (S, T) = ( a   *  ) ( b  *  )  1
                                        ( 0  S22 ),( 0 T22 )  n-1
                                          1  n-1     1 n-1

     Then the reciprocal condition number DIF(i) is

                Difl((a, b), (S22, T22)) = sigma-min( Zl ),

     where sigma-min(Zl) denotes the smallest singular value of the
     2(n-1)-by-2(n-1) matrix

         Zl = [ kron(a, In-1)  -kron(1, S22) ]
              [ kron(b, In-1)  -kron(1, T22) ] .

     Here In-1 is the identity matrix of size n-1. kron(X, Y) is the
     Kronecker product between the matrices X and Y.

     Note that if the default method for computing DIF(i) is wanted
     (see SLATDF), then the parameter DIFDRI (see below) should be
     changed from 3 to 4 (routine SLATDF(IJOB = 2 will be used)).
     See STGSYL for more details.

  b) If the i-th and (i+1)-th eigenvalues are complex conjugate pair,

     Suppose U and V are orthogonal transformations such that

              U**T*(A, B)*V = (S, T) = ( S11  *   ) ( T11  *  )  2
                                       ( 0    S22 ),( 0    T22) n-2
                                         2    n-2     2    n-2

     and (S11, T11) corresponds to the complex conjugate eigenvalue
     pair (w, conjg(w)). There exist unitary matrices U1 and V1 such
     that

       U1**T*S11*V1 = ( s11 s12 ) and U1**T*T11*V1 = ( t11 t12 )
                      (  0  s22 )                    (  0  t22 )

     where the generalized eigenvalues w = s11/t11 and
     conjg(w) = s22/t22.

     Then the reciprocal condition number DIF(i) is bounded by

         min( d1, max( 1, |real(s11)/real(s22)| )*d2 )

     where, d1 = Difl((s11, t11), (s22, t22)) = sigma-min(Z1), where
     Z1 is the complex 2-by-2 matrix

              Z1 =  [ s11  -s22 ]
                    [ t11  -t22 ],

     This is done by computing (using real arithmetic) the
     roots of the characteristical polynomial det(Z1**T * Z1 - lambda I),
     where Z1**T denotes the transpose of Z1 and det(X) denotes
     the determinant of X.

     and d2 is an upper bound on Difl((S11, T11), (S22, T22)), i.e. an
     upper bound on sigma-min(Z2), where Z2 is (2n-2)-by-(2n-2)

              Z2 = [ kron(S11**T, In-2)  -kron(I2, S22) ]
                   [ kron(T11**T, In-2)  -kron(I2, T22) ]

     Note that if the default method for computing DIF is wanted (see
     SLATDF), then the parameter DIFDRI (see below) should be changed
     from 3 to 4 (routine SLATDF(IJOB = 2 will be used)). See STGSYL
     for more details.

  For each eigenvalue/vector specified by SELECT, DIF stores a
  Frobenius norm-based estimate of Difl.

  An approximate error bound for the i-th computed eigenvector VL(i) or
  VR(i) is given by

             EPS * norm(A, B) / DIF(i).

  See ref. [2-3] for more details and further references.
Contributors:
Bo Kagstrom and Peter Poromaa, Department of Computing Science, Umea University, S-901 87 Umea, Sweden.
References:
  [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
      Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
      M.S. Moonen et al (eds), Linear Algebra for Large Scale and
      Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.

  [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
      Eigenvalues of a Regular Matrix Pair (A, B) and Condition
      Estimation: Theory, Algorithms and Software,
      Report UMINF - 94.04, Department of Computing Science, Umea
      University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
      Note 87. To appear in Numerical Algorithms, 1996.

  [3] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software
      for Solving the Generalized Sylvester Equation and Estimating the
      Separation between Regular Matrix Pairs, Report UMINF - 93.23,
      Department of Computing Science, Umea University, S-901 87 Umea,
      Sweden, December 1993, Revised April 1994, Also as LAPACK Working
      Note 75.  To appear in ACM Trans. on Math. Software, Vol 22,
      No 1, 1996.

Definition at line 378 of file stgsna.f.

381*
382* -- LAPACK computational routine --
383* -- LAPACK is a software package provided by Univ. of Tennessee, --
384* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
385*
386* .. Scalar Arguments ..
387 CHARACTER HOWMNY, JOB
388 INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, M, MM, N
389* ..
390* .. Array Arguments ..
391 LOGICAL SELECT( * )
392 INTEGER IWORK( * )
393 REAL A( LDA, * ), B( LDB, * ), DIF( * ), S( * ),
394 $ VL( LDVL, * ), VR( LDVR, * ), WORK( * )
395* ..
396*
397* =====================================================================
398*
399* .. Parameters ..
400 INTEGER DIFDRI
401 parameter( difdri = 3 )
402 REAL ZERO, ONE, TWO, FOUR
403 parameter( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0,
404 $ four = 4.0e+0 )
405* ..
406* .. Local Scalars ..
407 LOGICAL LQUERY, PAIR, SOMCON, WANTBH, WANTDF, WANTS
408 INTEGER I, IERR, IFST, ILST, IZ, K, KS, LWMIN, N1, N2
409 REAL ALPHAI, ALPHAR, ALPRQT, BETA, C1, C2, COND,
410 $ EPS, LNRM, RNRM, ROOT1, ROOT2, SCALE, SMLNUM,
411 $ TMPII, TMPIR, TMPRI, TMPRR, UHAV, UHAVI, UHBV,
412 $ UHBVI
413* ..
414* .. Local Arrays ..
415 REAL DUMMY( 1 ), DUMMY1( 1 )
416* ..
417* .. External Functions ..
418 LOGICAL LSAME
419 REAL SDOT, SLAMCH, SLAPY2, SNRM2, SROUNDUP_LWORK
420 EXTERNAL lsame, sdot, slamch, slapy2, snrm2,
422* ..
423* .. External Subroutines ..
424 EXTERNAL sgemv, slacpy, slag2, stgexc, stgsyl, xerbla
425* ..
426* .. Intrinsic Functions ..
427 INTRINSIC max, min, sqrt
428* ..
429* .. Executable Statements ..
430*
431* Decode and test the input parameters
432*
433 wantbh = lsame( job, 'B' )
434 wants = lsame( job, 'E' ) .OR. wantbh
435 wantdf = lsame( job, 'V' ) .OR. wantbh
436*
437 somcon = lsame( howmny, 'S' )
438*
439 info = 0
440 lquery = ( lwork.EQ.-1 )
441*
442 IF( .NOT.wants .AND. .NOT.wantdf ) THEN
443 info = -1
444 ELSE IF( .NOT.lsame( howmny, 'A' ) .AND. .NOT.somcon ) THEN
445 info = -2
446 ELSE IF( n.LT.0 ) THEN
447 info = -4
448 ELSE IF( lda.LT.max( 1, n ) ) THEN
449 info = -6
450 ELSE IF( ldb.LT.max( 1, n ) ) THEN
451 info = -8
452 ELSE IF( wants .AND. ldvl.LT.n ) THEN
453 info = -10
454 ELSE IF( wants .AND. ldvr.LT.n ) THEN
455 info = -12
456 ELSE
457*
458* Set M to the number of eigenpairs for which condition numbers
459* are required, and test MM.
460*
461 IF( somcon ) THEN
462 m = 0
463 pair = .false.
464 DO 10 k = 1, n
465 IF( pair ) THEN
466 pair = .false.
467 ELSE
468 IF( k.LT.n ) THEN
469 IF( a( k+1, k ).EQ.zero ) THEN
470 IF( SELECT( k ) )
471 $ m = m + 1
472 ELSE
473 pair = .true.
474 IF( SELECT( k ) .OR. SELECT( k+1 ) )
475 $ m = m + 2
476 END IF
477 ELSE
478 IF( SELECT( n ) )
479 $ m = m + 1
480 END IF
481 END IF
482 10 CONTINUE
483 ELSE
484 m = n
485 END IF
486*
487 IF( n.EQ.0 ) THEN
488 lwmin = 1
489 ELSE IF( lsame( job, 'V' ) .OR. lsame( job, 'B' ) ) THEN
490 lwmin = 2*n*( n + 2 ) + 16
491 ELSE
492 lwmin = n
493 END IF
494 work( 1 ) = sroundup_lwork(lwmin)
495*
496 IF( mm.LT.m ) THEN
497 info = -15
498 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
499 info = -18
500 END IF
501 END IF
502*
503 IF( info.NE.0 ) THEN
504 CALL xerbla( 'STGSNA', -info )
505 RETURN
506 ELSE IF( lquery ) THEN
507 RETURN
508 END IF
509*
510* Quick return if possible
511*
512 IF( n.EQ.0 )
513 $ RETURN
514*
515* Get machine constants
516*
517 eps = slamch( 'P' )
518 smlnum = slamch( 'S' ) / eps
519 ks = 0
520 pair = .false.
521*
522 DO 20 k = 1, n
523*
524* Determine whether A(k,k) begins a 1-by-1 or 2-by-2 block.
525*
526 IF( pair ) THEN
527 pair = .false.
528 GO TO 20
529 ELSE
530 IF( k.LT.n )
531 $ pair = a( k+1, k ).NE.zero
532 END IF
533*
534* Determine whether condition numbers are required for the k-th
535* eigenpair.
536*
537 IF( somcon ) THEN
538 IF( pair ) THEN
539 IF( .NOT.SELECT( k ) .AND. .NOT.SELECT( k+1 ) )
540 $ GO TO 20
541 ELSE
542 IF( .NOT.SELECT( k ) )
543 $ GO TO 20
544 END IF
545 END IF
546*
547 ks = ks + 1
548*
549 IF( wants ) THEN
550*
551* Compute the reciprocal condition number of the k-th
552* eigenvalue.
553*
554 IF( pair ) THEN
555*
556* Complex eigenvalue pair.
557*
558 rnrm = slapy2( snrm2( n, vr( 1, ks ), 1 ),
559 $ snrm2( n, vr( 1, ks+1 ), 1 ) )
560 lnrm = slapy2( snrm2( n, vl( 1, ks ), 1 ),
561 $ snrm2( n, vl( 1, ks+1 ), 1 ) )
562 CALL sgemv( 'N', n, n, one, a, lda, vr( 1, ks ), 1, zero,
563 $ work, 1 )
564 tmprr = sdot( n, work, 1, vl( 1, ks ), 1 )
565 tmpri = sdot( n, work, 1, vl( 1, ks+1 ), 1 )
566 CALL sgemv( 'N', n, n, one, a, lda, vr( 1, ks+1 ), 1,
567 $ zero, work, 1 )
568 tmpii = sdot( n, work, 1, vl( 1, ks+1 ), 1 )
569 tmpir = sdot( n, work, 1, vl( 1, ks ), 1 )
570 uhav = tmprr + tmpii
571 uhavi = tmpir - tmpri
572 CALL sgemv( 'N', n, n, one, b, ldb, vr( 1, ks ), 1, zero,
573 $ work, 1 )
574 tmprr = sdot( n, work, 1, vl( 1, ks ), 1 )
575 tmpri = sdot( n, work, 1, vl( 1, ks+1 ), 1 )
576 CALL sgemv( 'N', n, n, one, b, ldb, vr( 1, ks+1 ), 1,
577 $ zero, work, 1 )
578 tmpii = sdot( n, work, 1, vl( 1, ks+1 ), 1 )
579 tmpir = sdot( n, work, 1, vl( 1, ks ), 1 )
580 uhbv = tmprr + tmpii
581 uhbvi = tmpir - tmpri
582 uhav = slapy2( uhav, uhavi )
583 uhbv = slapy2( uhbv, uhbvi )
584 cond = slapy2( uhav, uhbv )
585 s( ks ) = cond / ( rnrm*lnrm )
586 s( ks+1 ) = s( ks )
587*
588 ELSE
589*
590* Real eigenvalue.
591*
592 rnrm = snrm2( n, vr( 1, ks ), 1 )
593 lnrm = snrm2( n, vl( 1, ks ), 1 )
594 CALL sgemv( 'N', n, n, one, a, lda, vr( 1, ks ), 1, zero,
595 $ work, 1 )
596 uhav = sdot( n, work, 1, vl( 1, ks ), 1 )
597 CALL sgemv( 'N', n, n, one, b, ldb, vr( 1, ks ), 1, zero,
598 $ work, 1 )
599 uhbv = sdot( n, work, 1, vl( 1, ks ), 1 )
600 cond = slapy2( uhav, uhbv )
601 IF( cond.EQ.zero ) THEN
602 s( ks ) = -one
603 ELSE
604 s( ks ) = cond / ( rnrm*lnrm )
605 END IF
606 END IF
607 END IF
608*
609 IF( wantdf ) THEN
610 IF( n.EQ.1 ) THEN
611 dif( ks ) = slapy2( a( 1, 1 ), b( 1, 1 ) )
612 GO TO 20
613 END IF
614*
615* Estimate the reciprocal condition number of the k-th
616* eigenvectors.
617 IF( pair ) THEN
618*
619* Copy the 2-by 2 pencil beginning at (A(k,k), B(k, k)).
620* Compute the eigenvalue(s) at position K.
621*
622 work( 1 ) = a( k, k )
623 work( 2 ) = a( k+1, k )
624 work( 3 ) = a( k, k+1 )
625 work( 4 ) = a( k+1, k+1 )
626 work( 5 ) = b( k, k )
627 work( 6 ) = b( k+1, k )
628 work( 7 ) = b( k, k+1 )
629 work( 8 ) = b( k+1, k+1 )
630 CALL slag2( work, 2, work( 5 ), 2, smlnum*eps, beta,
631 $ dummy1( 1 ), alphar, dummy( 1 ), alphai )
632 alprqt = one
633 c1 = two*( alphar*alphar+alphai*alphai+beta*beta )
634 c2 = four*beta*beta*alphai*alphai
635 root1 = c1 + sqrt( c1*c1-4.0*c2 )
636 root1 = root1 / two
637 root2 = c2 / root1
638 cond = min( sqrt( root1 ), sqrt( root2 ) )
639 END IF
640*
641* Copy the matrix (A, B) to the array WORK and swap the
642* diagonal block beginning at A(k,k) to the (1,1) position.
643*
644 CALL slacpy( 'Full', n, n, a, lda, work, n )
645 CALL slacpy( 'Full', n, n, b, ldb, work( n*n+1 ), n )
646 ifst = k
647 ilst = 1
648*
649 CALL stgexc( .false., .false., n, work, n, work( n*n+1 ), n,
650 $ dummy, 1, dummy1, 1, ifst, ilst,
651 $ work( n*n*2+1 ), lwork-2*n*n, ierr )
652*
653 IF( ierr.GT.0 ) THEN
654*
655* Ill-conditioned problem - swap rejected.
656*
657 dif( ks ) = zero
658 ELSE
659*
660* Reordering successful, solve generalized Sylvester
661* equation for R and L,
662* A22 * R - L * A11 = A12
663* B22 * R - L * B11 = B12,
664* and compute estimate of Difl((A11,B11), (A22, B22)).
665*
666 n1 = 1
667 IF( work( 2 ).NE.zero )
668 $ n1 = 2
669 n2 = n - n1
670 IF( n2.EQ.0 ) THEN
671 dif( ks ) = cond
672 ELSE
673 i = n*n + 1
674 iz = 2*n*n + 1
675 CALL stgsyl( 'N', difdri, n2, n1, work( n*n1+n1+1 ),
676 $ n, work, n, work( n1+1 ), n,
677 $ work( n*n1+n1+i ), n, work( i ), n,
678 $ work( n1+i ), n, scale, dif( ks ),
679 $ work( iz+1 ), lwork-2*n*n, iwork, ierr )
680*
681 IF( pair )
682 $ dif( ks ) = min( max( one, alprqt )*dif( ks ),
683 $ cond )
684 END IF
685 END IF
686 IF( pair )
687 $ dif( ks+1 ) = dif( ks )
688 END IF
689 IF( pair )
690 $ ks = ks + 1
691*
692 20 CONTINUE
693 work( 1 ) = sroundup_lwork(lwmin)
694 RETURN
695*
696* End of STGSNA
697*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
real function sdot(n, sx, incx, sy, incy)
SDOT
Definition sdot.f:82
subroutine sgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
SGEMV
Definition sgemv.f:158
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slag2(a, lda, b, ldb, safmin, scale1, scale2, wr1, wr2, wi)
SLAG2 computes the eigenvalues of a 2-by-2 generalized eigenvalue problem, with scaling as necessary ...
Definition slag2.f:156
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slapy2(x, y)
SLAPY2 returns sqrt(x2+y2).
Definition slapy2.f:63
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function snrm2(n, x, incx)
SNRM2
Definition snrm2.f90:89
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine stgexc(wantq, wantz, n, a, lda, b, ldb, q, ldq, z, ldz, ifst, ilst, work, lwork, info)
STGEXC
Definition stgexc.f:220
subroutine stgsyl(trans, ijob, m, n, a, lda, b, ldb, c, ldc, d, ldd, e, lde, f, ldf, scale, dif, work, lwork, iwork, info)
STGSYL
Definition stgsyl.f:299
Here is the call graph for this function:
Here is the caller graph for this function: