LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
ctgsna.f
Go to the documentation of this file.
1*> \brief \b CTGSNA
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CTGSNA + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ctgsna.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ctgsna.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ctgsna.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CTGSNA( JOB, HOWMNY, SELECT, N, A, LDA, B, LDB, VL,
20* LDVL, VR, LDVR, S, DIF, MM, M, WORK, LWORK,
21* IWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER HOWMNY, JOB
25* INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, M, MM, N
26* ..
27* .. Array Arguments ..
28* LOGICAL SELECT( * )
29* INTEGER IWORK( * )
30* REAL DIF( * ), S( * )
31* COMPLEX A( LDA, * ), B( LDB, * ), VL( LDVL, * ),
32* $ VR( LDVR, * ), WORK( * )
33* ..
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> CTGSNA estimates reciprocal condition numbers for specified
42*> eigenvalues and/or eigenvectors of a matrix pair (A, B).
43*>
44*> (A, B) must be in generalized Schur canonical form, that is, A and
45*> B are both upper triangular.
46*> \endverbatim
47*
48* Arguments:
49* ==========
50*
51*> \param[in] JOB
52*> \verbatim
53*> JOB is CHARACTER*1
54*> Specifies whether condition numbers are required for
55*> eigenvalues (S) or eigenvectors (DIF):
56*> = 'E': for eigenvalues only (S);
57*> = 'V': for eigenvectors only (DIF);
58*> = 'B': for both eigenvalues and eigenvectors (S and DIF).
59*> \endverbatim
60*>
61*> \param[in] HOWMNY
62*> \verbatim
63*> HOWMNY is CHARACTER*1
64*> = 'A': compute condition numbers for all eigenpairs;
65*> = 'S': compute condition numbers for selected eigenpairs
66*> specified by the array SELECT.
67*> \endverbatim
68*>
69*> \param[in] SELECT
70*> \verbatim
71*> SELECT is LOGICAL array, dimension (N)
72*> If HOWMNY = 'S', SELECT specifies the eigenpairs for which
73*> condition numbers are required. To select condition numbers
74*> for the corresponding j-th eigenvalue and/or eigenvector,
75*> SELECT(j) must be set to .TRUE..
76*> If HOWMNY = 'A', SELECT is not referenced.
77*> \endverbatim
78*>
79*> \param[in] N
80*> \verbatim
81*> N is INTEGER
82*> The order of the square matrix pair (A, B). N >= 0.
83*> \endverbatim
84*>
85*> \param[in] A
86*> \verbatim
87*> A is COMPLEX array, dimension (LDA,N)
88*> The upper triangular matrix A in the pair (A,B).
89*> \endverbatim
90*>
91*> \param[in] LDA
92*> \verbatim
93*> LDA is INTEGER
94*> The leading dimension of the array A. LDA >= max(1,N).
95*> \endverbatim
96*>
97*> \param[in] B
98*> \verbatim
99*> B is COMPLEX array, dimension (LDB,N)
100*> The upper triangular matrix B in the pair (A, B).
101*> \endverbatim
102*>
103*> \param[in] LDB
104*> \verbatim
105*> LDB is INTEGER
106*> The leading dimension of the array B. LDB >= max(1,N).
107*> \endverbatim
108*>
109*> \param[in] VL
110*> \verbatim
111*> VL is COMPLEX array, dimension (LDVL,M)
112*> IF JOB = 'E' or 'B', VL must contain left eigenvectors of
113*> (A, B), corresponding to the eigenpairs specified by HOWMNY
114*> and SELECT. The eigenvectors must be stored in consecutive
115*> columns of VL, as returned by CTGEVC.
116*> If JOB = 'V', VL is not referenced.
117*> \endverbatim
118*>
119*> \param[in] LDVL
120*> \verbatim
121*> LDVL is INTEGER
122*> The leading dimension of the array VL. LDVL >= 1; and
123*> If JOB = 'E' or 'B', LDVL >= N.
124*> \endverbatim
125*>
126*> \param[in] VR
127*> \verbatim
128*> VR is COMPLEX array, dimension (LDVR,M)
129*> IF JOB = 'E' or 'B', VR must contain right eigenvectors of
130*> (A, B), corresponding to the eigenpairs specified by HOWMNY
131*> and SELECT. The eigenvectors must be stored in consecutive
132*> columns of VR, as returned by CTGEVC.
133*> If JOB = 'V', VR is not referenced.
134*> \endverbatim
135*>
136*> \param[in] LDVR
137*> \verbatim
138*> LDVR is INTEGER
139*> The leading dimension of the array VR. LDVR >= 1;
140*> If JOB = 'E' or 'B', LDVR >= N.
141*> \endverbatim
142*>
143*> \param[out] S
144*> \verbatim
145*> S is REAL array, dimension (MM)
146*> If JOB = 'E' or 'B', the reciprocal condition numbers of the
147*> selected eigenvalues, stored in consecutive elements of the
148*> array.
149*> If JOB = 'V', S is not referenced.
150*> \endverbatim
151*>
152*> \param[out] DIF
153*> \verbatim
154*> DIF is REAL array, dimension (MM)
155*> If JOB = 'V' or 'B', the estimated reciprocal condition
156*> numbers of the selected eigenvectors, stored in consecutive
157*> elements of the array.
158*> If the eigenvalues cannot be reordered to compute DIF(j),
159*> DIF(j) is set to 0; this can only occur when the true value
160*> would be very small anyway.
161*> For each eigenvalue/vector specified by SELECT, DIF stores
162*> a Frobenius norm-based estimate of Difl.
163*> If JOB = 'E', DIF is not referenced.
164*> \endverbatim
165*>
166*> \param[in] MM
167*> \verbatim
168*> MM is INTEGER
169*> The number of elements in the arrays S and DIF. MM >= M.
170*> \endverbatim
171*>
172*> \param[out] M
173*> \verbatim
174*> M is INTEGER
175*> The number of elements of the arrays S and DIF used to store
176*> the specified condition numbers; for each selected eigenvalue
177*> one element is used. If HOWMNY = 'A', M is set to N.
178*> \endverbatim
179*>
180*> \param[out] WORK
181*> \verbatim
182*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
183*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
184*> \endverbatim
185*>
186*> \param[in] LWORK
187*> \verbatim
188*> LWORK is INTEGER
189*> The dimension of the array WORK. LWORK >= max(1,N).
190*> If JOB = 'V' or 'B', LWORK >= max(1,2*N*N).
191*> \endverbatim
192*>
193*> \param[out] IWORK
194*> \verbatim
195*> IWORK is INTEGER array, dimension (N+2)
196*> If JOB = 'E', IWORK is not referenced.
197*> \endverbatim
198*>
199*> \param[out] INFO
200*> \verbatim
201*> INFO is INTEGER
202*> = 0: Successful exit
203*> < 0: If INFO = -i, the i-th argument had an illegal value
204*> \endverbatim
205*
206* Authors:
207* ========
208*
209*> \author Univ. of Tennessee
210*> \author Univ. of California Berkeley
211*> \author Univ. of Colorado Denver
212*> \author NAG Ltd.
213*
214*> \ingroup tgsna
215*
216*> \par Further Details:
217* =====================
218*>
219*> \verbatim
220*>
221*> The reciprocal of the condition number of the i-th generalized
222*> eigenvalue w = (a, b) is defined as
223*>
224*> S(I) = (|v**HAu|**2 + |v**HBu|**2)**(1/2) / (norm(u)*norm(v))
225*>
226*> where u and v are the right and left eigenvectors of (A, B)
227*> corresponding to w; |z| denotes the absolute value of the complex
228*> number, and norm(u) denotes the 2-norm of the vector u. The pair
229*> (a, b) corresponds to an eigenvalue w = a/b (= v**HAu/v**HBu) of the
230*> matrix pair (A, B). If both a and b equal zero, then (A,B) is
231*> singular and S(I) = -1 is returned.
232*>
233*> An approximate error bound on the chordal distance between the i-th
234*> computed generalized eigenvalue w and the corresponding exact
235*> eigenvalue lambda is
236*>
237*> chord(w, lambda) <= EPS * norm(A, B) / S(I),
238*>
239*> where EPS is the machine precision.
240*>
241*> The reciprocal of the condition number of the right eigenvector u
242*> and left eigenvector v corresponding to the generalized eigenvalue w
243*> is defined as follows. Suppose
244*>
245*> (A, B) = ( a * ) ( b * ) 1
246*> ( 0 A22 ),( 0 B22 ) n-1
247*> 1 n-1 1 n-1
248*>
249*> Then the reciprocal condition number DIF(I) is
250*>
251*> Difl[(a, b), (A22, B22)] = sigma-min( Zl )
252*>
253*> where sigma-min(Zl) denotes the smallest singular value of
254*>
255*> Zl = [ kron(a, In-1) -kron(1, A22) ]
256*> [ kron(b, In-1) -kron(1, B22) ].
257*>
258*> Here In-1 is the identity matrix of size n-1 and X**H is the conjugate
259*> transpose of X. kron(X, Y) is the Kronecker product between the
260*> matrices X and Y.
261*>
262*> We approximate the smallest singular value of Zl with an upper
263*> bound. This is done by CLATDF.
264*>
265*> An approximate error bound for a computed eigenvector VL(i) or
266*> VR(i) is given by
267*>
268*> EPS * norm(A, B) / DIF(i).
269*>
270*> See ref. [2-3] for more details and further references.
271*> \endverbatim
272*
273*> \par Contributors:
274* ==================
275*>
276*> Bo Kagstrom and Peter Poromaa, Department of Computing Science,
277*> Umea University, S-901 87 Umea, Sweden.
278*
279*> \par References:
280* ================
281*>
282*> \verbatim
283*>
284*> [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
285*> Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
286*> M.S. Moonen et al (eds), Linear Algebra for Large Scale and
287*> Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.
288*>
289*> [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
290*> Eigenvalues of a Regular Matrix Pair (A, B) and Condition
291*> Estimation: Theory, Algorithms and Software, Report
292*> UMINF - 94.04, Department of Computing Science, Umea University,
293*> S-901 87 Umea, Sweden, 1994. Also as LAPACK Working Note 87.
294*> To appear in Numerical Algorithms, 1996.
295*>
296*> [3] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software
297*> for Solving the Generalized Sylvester Equation and Estimating the
298*> Separation between Regular Matrix Pairs, Report UMINF - 93.23,
299*> Department of Computing Science, Umea University, S-901 87 Umea,
300*> Sweden, December 1993, Revised April 1994, Also as LAPACK Working
301*> Note 75.
302*> To appear in ACM Trans. on Math. Software, Vol 22, No 1, 1996.
303*> \endverbatim
304*>
305* =====================================================================
306 SUBROUTINE ctgsna( JOB, HOWMNY, SELECT, N, A, LDA, B, LDB, VL,
307 $ LDVL, VR, LDVR, S, DIF, MM, M, WORK, LWORK,
308 $ IWORK, INFO )
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,
346 COMPLEX CDOTC
347 EXTERNAL LSAME, SCNRM2, SLAMCH,
348 $ slapy2, sroundup_lwork,
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*
520 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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 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 ctgsna(job, howmny, select, n, a, lda, b, ldb, vl, ldvl, vr, ldvr, s, dif, mm, m, work, lwork, iwork, info)
CTGSNA
Definition ctgsna.f:309
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