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

◆ ctgsna()

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

CTGSNA

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

Purpose:
!>
!> CTGSNA estimates reciprocal condition numbers for specified
!> eigenvalues and/or eigenvectors of a matrix pair (A, B).
!>
!> (A, B) must be in generalized Schur canonical form, that is, A and
!> B are both 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 corresponding j-th eigenvalue and/or eigenvector,
!>          SELECT(j) 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 COMPLEX array, dimension (LDA,N)
!>          The upper 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 COMPLEX 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 COMPLEX 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 CTGEVC.
!>          If JOB = 'V', VL is not referenced.
!> 
[in]LDVL
!>          LDVL is INTEGER
!>          The leading dimension of the array VL. LDVL >= 1; and
!>          If JOB = 'E' or 'B', LDVL >= N.
!> 
[in]VR
!>          VR is COMPLEX 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 of VR, as returned by CTGEVC.
!>          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.
!>          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.
!>          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.
!>          For each eigenvalue/vector specified by SELECT, DIF stores
!>          a Frobenius norm-based estimate of Difl.
!>          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 eigenvalue
!>          one element is used. If HOWMNY = 'A', M is set to 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 JOB = 'V' or 'B', LWORK >= max(1,2*N*N).
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (N+2)
!>          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 the i-th generalized
!>  eigenvalue w = (a, b) is defined as
!>
!>          S(I) = (|v**HAu|**2 + |v**HBu|**2)**(1/2) / (norm(u)*norm(v))
!>
!>  where u and v are the right and left 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 (= v**HAu/v**HBu) 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 of the right eigenvector u
!>  and left eigenvector v corresponding to the generalized eigenvalue w
!>  is defined as follows. Suppose
!>
!>                   (A, B) = ( a   *  ) ( b  *  )  1
!>                            ( 0  A22 ),( 0 B22 )  n-1
!>                              1  n-1     1 n-1
!>
!>  Then the reciprocal condition number DIF(I) is
!>
!>          Difl[(a, b), (A22, B22)]  = sigma-min( Zl )
!>
!>  where sigma-min(Zl) denotes the smallest singular value of
!>
!>         Zl = [ kron(a, In-1) -kron(1, A22) ]
!>              [ kron(b, In-1) -kron(1, B22) ].
!>
!>  Here In-1 is the identity matrix of size n-1 and X**H is the conjugate
!>  transpose of X. kron(X, Y) is the Kronecker product between the
!>  matrices X and Y.
!>
!>  We approximate the smallest singular value of Zl with an upper
!>  bound. This is done by CLATDF.
!>
!>  An approximate error bound for a 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 306 of file ctgsna.f.

309*
310* -- LAPACK computational routine --
311* -- LAPACK is a software package provided by Univ. of Tennessee, --
312* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
313*
314* .. Scalar Arguments ..
315 CHARACTER HOWMNY, JOB
316 INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, M, MM, N
317* ..
318* .. Array Arguments ..
319 LOGICAL SELECT( * )
320 INTEGER IWORK( * )
321 REAL DIF( * ), S( * )
322 COMPLEX A( LDA, * ), B( LDB, * ), VL( LDVL, * ),
323 $ VR( LDVR, * ), WORK( * )
324* ..
325*
326* =====================================================================
327*
328* .. Parameters ..
329 REAL ZERO, ONE
330 INTEGER IDIFJB
331 parameter( zero = 0.0e+0, one = 1.0e+0, idifjb = 3 )
332* ..
333* .. Local Scalars ..
334 LOGICAL LQUERY, SOMCON, WANTBH, WANTDF, WANTS
335 INTEGER I, IERR, IFST, ILST, K, KS, LWMIN, N1, N2
336 REAL BIGNUM, COND, EPS, LNRM, RNRM, SCALE, SMLNUM
337 COMPLEX YHAX, YHBX
338* ..
339* .. Local Arrays ..
340 COMPLEX DUMMY( 1 ), DUMMY1( 1 )
341* ..
342* .. External Functions ..
343 LOGICAL LSAME
344 REAL SCNRM2, SLAMCH, SLAPY2,
345 $ SROUNDUP_LWORK
346 COMPLEX CDOTC
347 EXTERNAL lsame, scnrm2, slamch,
349 $ cdotc
350* ..
351* .. External Subroutines ..
352 EXTERNAL cgemv, clacpy, ctgexc, ctgsyl,
353 $ xerbla
354* ..
355* .. Intrinsic Functions ..
356 INTRINSIC abs, cmplx, max
357* ..
358* .. Executable Statements ..
359*
360* Decode and test the input parameters
361*
362 wantbh = lsame( job, 'B' )
363 wants = lsame( job, 'E' ) .OR. wantbh
364 wantdf = lsame( job, 'V' ) .OR. wantbh
365*
366 somcon = lsame( howmny, 'S' )
367*
368 info = 0
369 lquery = ( lwork.EQ.-1 )
370*
371 IF( .NOT.wants .AND. .NOT.wantdf ) THEN
372 info = -1
373 ELSE IF( .NOT.lsame( howmny, 'A' ) .AND. .NOT.somcon ) THEN
374 info = -2
375 ELSE IF( n.LT.0 ) THEN
376 info = -4
377 ELSE IF( lda.LT.max( 1, n ) ) THEN
378 info = -6
379 ELSE IF( ldb.LT.max( 1, n ) ) THEN
380 info = -8
381 ELSE IF( wants .AND. ldvl.LT.n ) THEN
382 info = -10
383 ELSE IF( wants .AND. ldvr.LT.n ) THEN
384 info = -12
385 ELSE
386*
387* Set M to the number of eigenpairs for which condition numbers
388* are required, and test MM.
389*
390 IF( somcon ) THEN
391 m = 0
392 DO 10 k = 1, n
393 IF( SELECT( k ) )
394 $ m = m + 1
395 10 CONTINUE
396 ELSE
397 m = n
398 END IF
399*
400 IF( n.EQ.0 ) THEN
401 lwmin = 1
402 ELSE IF( lsame( job, 'V' ) .OR. lsame( job, 'B' ) ) THEN
403 lwmin = 2*n*n
404 ELSE
405 lwmin = n
406 END IF
407 work( 1 ) = sroundup_lwork(lwmin)
408*
409 IF( mm.LT.m ) THEN
410 info = -15
411 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
412 info = -18
413 END IF
414 END IF
415*
416 IF( info.NE.0 ) THEN
417 CALL xerbla( 'CTGSNA', -info )
418 RETURN
419 ELSE IF( lquery ) THEN
420 RETURN
421 END IF
422*
423* Quick return if possible
424*
425 IF( n.EQ.0 )
426 $ RETURN
427*
428* Get machine constants
429*
430 eps = slamch( 'P' )
431 smlnum = slamch( 'S' ) / eps
432 bignum = one / smlnum
433 ks = 0
434 DO 20 k = 1, n
435*
436* Determine whether condition numbers are required for the k-th
437* eigenpair.
438*
439 IF( somcon ) THEN
440 IF( .NOT.SELECT( k ) )
441 $ GO TO 20
442 END IF
443*
444 ks = ks + 1
445*
446 IF( wants ) THEN
447*
448* Compute the reciprocal condition number of the k-th
449* eigenvalue.
450*
451 rnrm = scnrm2( n, vr( 1, ks ), 1 )
452 lnrm = scnrm2( n, vl( 1, ks ), 1 )
453 CALL cgemv( 'N', n, n, cmplx( one, zero ), a, lda,
454 $ vr( 1, ks ), 1, cmplx( zero, zero ), work, 1 )
455 yhax = cdotc( n, work, 1, vl( 1, ks ), 1 )
456 CALL cgemv( 'N', n, n, cmplx( one, zero ), b, ldb,
457 $ vr( 1, ks ), 1, cmplx( zero, zero ), work, 1 )
458 yhbx = cdotc( n, work, 1, vl( 1, ks ), 1 )
459 cond = slapy2( abs( yhax ), abs( yhbx ) )
460 IF( cond.EQ.zero ) THEN
461 s( ks ) = -one
462 ELSE
463 s( ks ) = cond / ( rnrm*lnrm )
464 END IF
465 END IF
466*
467 IF( wantdf ) THEN
468 IF( n.EQ.1 ) THEN
469 dif( ks ) = slapy2( abs( a( 1, 1 ) ), abs( b( 1,
470 $ 1 ) ) )
471 ELSE
472*
473* Estimate the reciprocal condition number of the k-th
474* eigenvectors.
475*
476* Copy the matrix (A, B) to the array WORK and move the
477* (k,k)th pair to the (1,1) position.
478*
479 CALL clacpy( 'Full', n, n, a, lda, work, n )
480 CALL clacpy( 'Full', n, n, b, ldb, work( n*n+1 ), n )
481 ifst = k
482 ilst = 1
483*
484 CALL ctgexc( .false., .false., n, work, n,
485 $ work( n*n+1 ),
486 $ n, dummy, 1, dummy1, 1, ifst, ilst, ierr )
487*
488 IF( ierr.GT.0 ) THEN
489*
490* Ill-conditioned problem - swap rejected.
491*
492 dif( ks ) = zero
493 ELSE
494*
495* Reordering successful, solve generalized Sylvester
496* equation for R and L,
497* A22 * R - L * A11 = A12
498* B22 * R - L * B11 = B12,
499* and compute estimate of Difl[(A11,B11), (A22, B22)].
500*
501 n1 = 1
502 n2 = n - n1
503 i = n*n + 1
504 CALL ctgsyl( 'N', idifjb, n2, n1,
505 $ work( n*n1+n1+1 ),
506 $ n, work, n, work( n1+1 ), n,
507 $ work( n*n1+n1+i ), n, work( i ), n,
508 $ work( n1+i ), n, scale, dif( ks ), dummy,
509 $ 1, iwork, ierr )
510 END IF
511 END IF
512 END IF
513*
514 20 CONTINUE
515 work( 1 ) = sroundup_lwork(lwmin)
516 RETURN
517*
518* End of CTGSNA
519*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
complex function cdotc(n, cx, incx, cy, incy)
CDOTC
Definition cdotc.f:83
subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
CGEMV
Definition cgemv.f:160
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slapy2(x, y)
SLAPY2 returns sqrt(x2+y2).
Definition slapy2.f:61
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function scnrm2(n, x, incx)
SCNRM2
Definition scnrm2.f90:90
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine ctgexc(wantq, wantz, n, a, lda, b, ldb, q, ldq, z, ldz, ifst, ilst, info)
CTGEXC
Definition ctgexc.f:198
subroutine ctgsyl(trans, ijob, m, n, a, lda, b, ldb, c, ldc, d, ldd, e, lde, f, ldf, scale, dif, work, lwork, iwork, info)
CTGSYL
Definition ctgsyl.f:294
Here is the call graph for this function:
Here is the caller graph for this function: