LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ztgsna.f
Go to the documentation of this file.
1*> \brief \b ZTGSNA
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZTGSNA + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ztgsna.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ztgsna.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ztgsna.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZTGSNA( 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* DOUBLE PRECISION DIF( * ), S( * )
31* COMPLEX*16 A( LDA, * ), B( LDB, * ), VL( LDVL, * ),
32* $ VR( LDVR, * ), WORK( * )
33* ..
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> ZTGSNA 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*16 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*16 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*16 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 ZTGEVC.
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*16 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 ZTGEVC.
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 DOUBLE PRECISION 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 DOUBLE PRECISION 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*16 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 ZLATDF.
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 ztgsna( 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 DOUBLE PRECISION DIF( * ), S( * )
322 COMPLEX*16 A( LDA, * ), B( LDB, * ), VL( LDVL, * ),
323 $ vr( ldvr, * ), work( * )
324* ..
325*
326* =====================================================================
327*
328* .. Parameters ..
329 DOUBLE PRECISION ZERO, ONE
330 INTEGER IDIFJB
331 parameter( zero = 0.0d+0, one = 1.0d+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 DOUBLE PRECISION BIGNUM, COND, EPS, LNRM, RNRM, SCALE, SMLNUM
337 COMPLEX*16 YHAX, YHBX
338* ..
339* .. Local Arrays ..
340 COMPLEX*16 DUMMY( 1 ), DUMMY1( 1 )
341* ..
342* .. External Functions ..
343 LOGICAL LSAME
344 DOUBLE PRECISION DLAMCH, DLAPY2, DZNRM2
345 COMPLEX*16 ZDOTC
346 EXTERNAL lsame, dlamch, dlapy2, dznrm2,
347 $ zdotc
348* ..
349* .. External Subroutines ..
350 EXTERNAL xerbla, zgemv, zlacpy, ztgexc,
351 $ ztgsyl
352* ..
353* .. Intrinsic Functions ..
354 INTRINSIC abs, dcmplx, max
355* ..
356* .. Executable Statements ..
357*
358* Decode and test the input parameters
359*
360 wantbh = lsame( job, 'B' )
361 wants = lsame( job, 'E' ) .OR. wantbh
362 wantdf = lsame( job, 'V' ) .OR. wantbh
363*
364 somcon = lsame( howmny, 'S' )
365*
366 info = 0
367 lquery = ( lwork.EQ.-1 )
368*
369 IF( .NOT.wants .AND. .NOT.wantdf ) THEN
370 info = -1
371 ELSE IF( .NOT.lsame( howmny, 'A' ) .AND. .NOT.somcon ) THEN
372 info = -2
373 ELSE IF( n.LT.0 ) THEN
374 info = -4
375 ELSE IF( lda.LT.max( 1, n ) ) THEN
376 info = -6
377 ELSE IF( ldb.LT.max( 1, n ) ) THEN
378 info = -8
379 ELSE IF( wants .AND. ldvl.LT.n ) THEN
380 info = -10
381 ELSE IF( wants .AND. ldvr.LT.n ) THEN
382 info = -12
383 ELSE
384*
385* Set M to the number of eigenpairs for which condition numbers
386* are required, and test MM.
387*
388 IF( somcon ) THEN
389 m = 0
390 DO 10 k = 1, n
391 IF( SELECT( k ) )
392 $ m = m + 1
393 10 CONTINUE
394 ELSE
395 m = n
396 END IF
397*
398 IF( n.EQ.0 ) THEN
399 lwmin = 1
400 ELSE IF( lsame( job, 'V' ) .OR. lsame( job, 'B' ) ) THEN
401 lwmin = 2*n*n
402 ELSE
403 lwmin = n
404 END IF
405 work( 1 ) = lwmin
406*
407 IF( mm.LT.m ) THEN
408 info = -15
409 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
410 info = -18
411 END IF
412 END IF
413*
414 IF( info.NE.0 ) THEN
415 CALL xerbla( 'ZTGSNA', -info )
416 RETURN
417 ELSE IF( lquery ) THEN
418 RETURN
419 END IF
420*
421* Quick return if possible
422*
423 IF( n.EQ.0 )
424 $ RETURN
425*
426* Get machine constants
427*
428 eps = dlamch( 'P' )
429 smlnum = dlamch( 'S' ) / eps
430 bignum = one / smlnum
431 ks = 0
432 DO 20 k = 1, n
433*
434* Determine whether condition numbers are required for the k-th
435* eigenpair.
436*
437 IF( somcon ) THEN
438 IF( .NOT.SELECT( k ) )
439 $ GO TO 20
440 END IF
441*
442 ks = ks + 1
443*
444 IF( wants ) THEN
445*
446* Compute the reciprocal condition number of the k-th
447* eigenvalue.
448*
449 rnrm = dznrm2( n, vr( 1, ks ), 1 )
450 lnrm = dznrm2( n, vl( 1, ks ), 1 )
451 CALL zgemv( 'N', n, n, dcmplx( one, zero ), a, lda,
452 $ vr( 1, ks ), 1, dcmplx( zero, zero ), work, 1 )
453 yhax = zdotc( n, work, 1, vl( 1, ks ), 1 )
454 CALL zgemv( 'N', n, n, dcmplx( one, zero ), b, ldb,
455 $ vr( 1, ks ), 1, dcmplx( zero, zero ), work, 1 )
456 yhbx = zdotc( n, work, 1, vl( 1, ks ), 1 )
457 cond = dlapy2( abs( yhax ), abs( yhbx ) )
458 IF( cond.EQ.zero ) THEN
459 s( ks ) = -one
460 ELSE
461 s( ks ) = cond / ( rnrm*lnrm )
462 END IF
463 END IF
464*
465 IF( wantdf ) THEN
466 IF( n.EQ.1 ) THEN
467 dif( ks ) = dlapy2( abs( a( 1, 1 ) ), abs( b( 1,
468 $ 1 ) ) )
469 ELSE
470*
471* Estimate the reciprocal condition number of the k-th
472* eigenvectors.
473*
474* Copy the matrix (A, B) to the array WORK and move the
475* (k,k)th pair to the (1,1) position.
476*
477 CALL zlacpy( 'Full', n, n, a, lda, work, n )
478 CALL zlacpy( 'Full', n, n, b, ldb, work( n*n+1 ), n )
479 ifst = k
480 ilst = 1
481*
482 CALL ztgexc( .false., .false., n, work, n,
483 $ work( n*n+1 ),
484 $ n, dummy, 1, dummy1, 1, ifst, ilst, ierr )
485*
486 IF( ierr.GT.0 ) THEN
487*
488* Ill-conditioned problem - swap rejected.
489*
490 dif( ks ) = zero
491 ELSE
492*
493* Reordering successful, solve generalized Sylvester
494* equation for R and L,
495* A22 * R - L * A11 = A12
496* B22 * R - L * B11 = B12,
497* and compute estimate of Difl[(A11,B11), (A22, B22)].
498*
499 n1 = 1
500 n2 = n - n1
501 i = n*n + 1
502 CALL ztgsyl( 'N', idifjb, n2, n1,
503 $ work( n*n1+n1+1 ),
504 $ n, work, n, work( n1+1 ), n,
505 $ work( n*n1+n1+i ), n, work( i ), n,
506 $ work( n1+i ), n, scale, dif( ks ), dummy,
507 $ 1, iwork, ierr )
508 END IF
509 END IF
510 END IF
511*
512 20 CONTINUE
513 work( 1 ) = lwmin
514 RETURN
515*
516* End of ZTGSNA
517*
518 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:160
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
subroutine ztgexc(wantq, wantz, n, a, lda, b, ldb, q, ldq, z, ldz, ifst, ilst, info)
ZTGEXC
Definition ztgexc.f:198
subroutine ztgsna(job, howmny, select, n, a, lda, b, ldb, vl, ldvl, vr, ldvr, s, dif, mm, m, work, lwork, iwork, info)
ZTGSNA
Definition ztgsna.f:309
subroutine ztgsyl(trans, ijob, m, n, a, lda, b, ldb, c, ldc, d, ldd, e, lde, f, ldf, scale, dif, work, lwork, iwork, info)
ZTGSYL
Definition ztgsyl.f:294