LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
zgeesx.f
Go to the documentation of this file.
1*> \brief <b> ZGEESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGEESX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgeesx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgeesx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgeesx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGEESX( JOBVS, SORT, SELECT, SENSE, N, A, LDA, SDIM, W,
20* VS, LDVS, RCONDE, RCONDV, WORK, LWORK, RWORK,
21* BWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER JOBVS, SENSE, SORT
25* INTEGER INFO, LDA, LDVS, LWORK, N, SDIM
26* DOUBLE PRECISION RCONDE, RCONDV
27* ..
28* .. Array Arguments ..
29* LOGICAL BWORK( * )
30* DOUBLE PRECISION RWORK( * )
31* COMPLEX*16 A( LDA, * ), VS( LDVS, * ), W( * ), WORK( * )
32* ..
33* .. Function Arguments ..
34* LOGICAL SELECT
35* EXTERNAL SELECT
36* ..
37*
38*
39*> \par Purpose:
40* =============
41*>
42*> \verbatim
43*>
44*> ZGEESX computes for an N-by-N complex nonsymmetric matrix A, the
45*> eigenvalues, the Schur form T, and, optionally, the matrix of Schur
46*> vectors Z. This gives the Schur factorization A = Z*T*(Z**H).
47*>
48*> Optionally, it also orders the eigenvalues on the diagonal of the
49*> Schur form so that selected eigenvalues are at the top left;
50*> computes a reciprocal condition number for the average of the
51*> selected eigenvalues (RCONDE); and computes a reciprocal condition
52*> number for the right invariant subspace corresponding to the
53*> selected eigenvalues (RCONDV). The leading columns of Z form an
54*> orthonormal basis for this invariant subspace.
55*>
56*> For further explanation of the reciprocal condition numbers RCONDE
57*> and RCONDV, see Section 4.10 of the LAPACK Users' Guide (where
58*> these quantities are called s and sep respectively).
59*>
60*> A complex matrix is in Schur form if it is upper triangular.
61*> \endverbatim
62*
63* Arguments:
64* ==========
65*
66*> \param[in] JOBVS
67*> \verbatim
68*> JOBVS is CHARACTER*1
69*> = 'N': Schur vectors are not computed;
70*> = 'V': Schur vectors are computed.
71*> \endverbatim
72*>
73*> \param[in] SORT
74*> \verbatim
75*> SORT is CHARACTER*1
76*> Specifies whether or not to order the eigenvalues on the
77*> diagonal of the Schur form.
78*> = 'N': Eigenvalues are not ordered;
79*> = 'S': Eigenvalues are ordered (see SELECT).
80*> \endverbatim
81*>
82*> \param[in] SELECT
83*> \verbatim
84*> SELECT is a LOGICAL FUNCTION of one COMPLEX*16 argument
85*> SELECT must be declared EXTERNAL in the calling subroutine.
86*> If SORT = 'S', SELECT is used to select eigenvalues to order
87*> to the top left of the Schur form.
88*> If SORT = 'N', SELECT is not referenced.
89*> An eigenvalue W(j) is selected if SELECT(W(j)) is true.
90*> \endverbatim
91*>
92*> \param[in] SENSE
93*> \verbatim
94*> SENSE is CHARACTER*1
95*> Determines which reciprocal condition numbers are computed.
96*> = 'N': None are computed;
97*> = 'E': Computed for average of selected eigenvalues only;
98*> = 'V': Computed for selected right invariant subspace only;
99*> = 'B': Computed for both.
100*> If SENSE = 'E', 'V' or 'B', SORT must equal 'S'.
101*> \endverbatim
102*>
103*> \param[in] N
104*> \verbatim
105*> N is INTEGER
106*> The order of the matrix A. N >= 0.
107*> \endverbatim
108*>
109*> \param[in,out] A
110*> \verbatim
111*> A is COMPLEX*16 array, dimension (LDA, N)
112*> On entry, the N-by-N matrix A.
113*> On exit, A is overwritten by its Schur form T.
114*> \endverbatim
115*>
116*> \param[in] LDA
117*> \verbatim
118*> LDA is INTEGER
119*> The leading dimension of the array A. LDA >= max(1,N).
120*> \endverbatim
121*>
122*> \param[out] SDIM
123*> \verbatim
124*> SDIM is INTEGER
125*> If SORT = 'N', SDIM = 0.
126*> If SORT = 'S', SDIM = number of eigenvalues for which
127*> SELECT is true.
128*> \endverbatim
129*>
130*> \param[out] W
131*> \verbatim
132*> W is COMPLEX*16 array, dimension (N)
133*> W contains the computed eigenvalues, in the same order
134*> that they appear on the diagonal of the output Schur form T.
135*> \endverbatim
136*>
137*> \param[out] VS
138*> \verbatim
139*> VS is COMPLEX*16 array, dimension (LDVS,N)
140*> If JOBVS = 'V', VS contains the unitary matrix Z of Schur
141*> vectors.
142*> If JOBVS = 'N', VS is not referenced.
143*> \endverbatim
144*>
145*> \param[in] LDVS
146*> \verbatim
147*> LDVS is INTEGER
148*> The leading dimension of the array VS. LDVS >= 1, and if
149*> JOBVS = 'V', LDVS >= N.
150*> \endverbatim
151*>
152*> \param[out] RCONDE
153*> \verbatim
154*> RCONDE is DOUBLE PRECISION
155*> If SENSE = 'E' or 'B', RCONDE contains the reciprocal
156*> condition number for the average of the selected eigenvalues.
157*> Not referenced if SENSE = 'N' or 'V'.
158*> \endverbatim
159*>
160*> \param[out] RCONDV
161*> \verbatim
162*> RCONDV is DOUBLE PRECISION
163*> If SENSE = 'V' or 'B', RCONDV contains the reciprocal
164*> condition number for the selected right invariant subspace.
165*> Not referenced if SENSE = 'N' or 'E'.
166*> \endverbatim
167*>
168*> \param[out] WORK
169*> \verbatim
170*> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
171*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
172*> \endverbatim
173*>
174*> \param[in] LWORK
175*> \verbatim
176*> LWORK is INTEGER
177*> The dimension of the array WORK. LWORK >= max(1,2*N).
178*> Also, if SENSE = 'E' or 'V' or 'B', LWORK >= 2*SDIM*(N-SDIM),
179*> where SDIM is the number of selected eigenvalues computed by
180*> this routine. Note that 2*SDIM*(N-SDIM) <= N*N/2. Note also
181*> that an error is only returned if LWORK < max(1,2*N), but if
182*> SENSE = 'E' or 'V' or 'B' this may not be large enough.
183*> For good performance, LWORK must generally be larger.
184*>
185*> If LWORK = -1, then a workspace query is assumed; the routine
186*> only calculates upper bound on the optimal size of the
187*> array WORK, returns this value as the first entry of the WORK
188*> array, and no error message related to LWORK is issued by
189*> XERBLA.
190*> \endverbatim
191*>
192*> \param[out] RWORK
193*> \verbatim
194*> RWORK is DOUBLE PRECISION array, dimension (N)
195*> \endverbatim
196*>
197*> \param[out] BWORK
198*> \verbatim
199*> BWORK is LOGICAL array, dimension (N)
200*> Not referenced if SORT = 'N'.
201*> \endverbatim
202*>
203*> \param[out] INFO
204*> \verbatim
205*> INFO is INTEGER
206*> = 0: successful exit
207*> < 0: if INFO = -i, the i-th argument had an illegal value.
208*> > 0: if INFO = i, and i is
209*> <= N: the QR algorithm failed to compute all the
210*> eigenvalues; elements 1:ILO-1 and i+1:N of W
211*> contain those eigenvalues which have converged; if
212*> JOBVS = 'V', VS contains the transformation which
213*> reduces A to its partially converged Schur form.
214*> = N+1: the eigenvalues could not be reordered because some
215*> eigenvalues were too close to separate (the problem
216*> is very ill-conditioned);
217*> = N+2: after reordering, roundoff changed values of some
218*> complex eigenvalues so that leading eigenvalues in
219*> the Schur form no longer satisfy SELECT=.TRUE. This
220*> could also be caused by underflow due to scaling.
221*> \endverbatim
222*
223* Authors:
224* ========
225*
226*> \author Univ. of Tennessee
227*> \author Univ. of California Berkeley
228*> \author Univ. of Colorado Denver
229*> \author NAG Ltd.
230*
231*> \ingroup geesx
232*
233* =====================================================================
234 SUBROUTINE zgeesx( JOBVS, SORT, SELECT, SENSE, N, A, LDA, SDIM,
235 $ W,
236 $ VS, LDVS, RCONDE, RCONDV, WORK, LWORK, RWORK,
237 $ BWORK, INFO )
238*
239* -- LAPACK driver routine --
240* -- LAPACK is a software package provided by Univ. of Tennessee, --
241* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
242*
243* .. Scalar Arguments ..
244 CHARACTER JOBVS, SENSE, SORT
245 INTEGER INFO, LDA, LDVS, LWORK, N, SDIM
246 DOUBLE PRECISION RCONDE, RCONDV
247* ..
248* .. Array Arguments ..
249 LOGICAL BWORK( * )
250 DOUBLE PRECISION RWORK( * )
251 COMPLEX*16 A( LDA, * ), VS( LDVS, * ), W( * ), WORK( * )
252* ..
253* .. Function Arguments ..
254 LOGICAL SELECT
255 EXTERNAL SELECT
256* ..
257*
258* =====================================================================
259*
260* .. Parameters ..
261 DOUBLE PRECISION ZERO, ONE
262 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
263* ..
264* .. Local Scalars ..
265 LOGICAL LQUERY, SCALEA, WANTSB, WANTSE, WANTSN, WANTST,
266 $ WANTSV, WANTVS
267 INTEGER HSWORK, I, IBAL, ICOND, IERR, IEVAL, IHI, ILO,
268 $ itau, iwrk, lwrk, maxwrk, minwrk
269 DOUBLE PRECISION ANRM, BIGNUM, CSCALE, EPS, SMLNUM
270* ..
271* .. Local Arrays ..
272 DOUBLE PRECISION DUM( 1 )
273* ..
274* .. External Subroutines ..
275 EXTERNAL dlascl, xerbla, zcopy, zgebak, zgebal,
276 $ zgehrd,
278* ..
279* .. External Functions ..
280 LOGICAL LSAME
281 INTEGER ILAENV
282 DOUBLE PRECISION DLAMCH, ZLANGE
283 EXTERNAL lsame, ilaenv, dlamch, zlange
284* ..
285* .. Intrinsic Functions ..
286 INTRINSIC max, sqrt
287* ..
288* .. Executable Statements ..
289*
290* Test the input arguments
291*
292 info = 0
293 wantvs = lsame( jobvs, 'V' )
294 wantst = lsame( sort, 'S' )
295 wantsn = lsame( sense, 'N' )
296 wantse = lsame( sense, 'E' )
297 wantsv = lsame( sense, 'V' )
298 wantsb = lsame( sense, 'B' )
299 lquery = ( lwork.EQ.-1 )
300*
301 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
302 info = -1
303 ELSE IF( ( .NOT.wantst ) .AND.
304 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
305 info = -2
306 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
307 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
308 info = -4
309 ELSE IF( n.LT.0 ) THEN
310 info = -5
311 ELSE IF( lda.LT.max( 1, n ) ) THEN
312 info = -7
313 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
314 info = -11
315 END IF
316*
317* Compute workspace
318* (Note: Comments in the code beginning "Workspace:" describe the
319* minimal amount of real workspace needed at that point in the
320* code, as well as the preferred amount for good performance.
321* CWorkspace refers to complex workspace, and RWorkspace to real
322* workspace. NB refers to the optimal block size for the
323* immediately following subroutine, as returned by ILAENV.
324* HSWORK refers to the workspace preferred by ZHSEQR, as
325* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
326* the worst case.
327* If SENSE = 'E', 'V' or 'B', then the amount of workspace needed
328* depends on SDIM, which is computed by the routine ZTRSEN later
329* in the code.)
330*
331 IF( info.EQ.0 ) THEN
332 IF( n.EQ.0 ) THEN
333 minwrk = 1
334 lwrk = 1
335 ELSE
336 maxwrk = n + n*ilaenv( 1, 'ZGEHRD', ' ', n, 1, n, 0 )
337 minwrk = 2*n
338*
339 CALL zhseqr( 'S', jobvs, n, 1, n, a, lda, w, vs, ldvs,
340 $ work, -1, ieval )
341 hswork = int( work( 1 ) )
342*
343 IF( .NOT.wantvs ) THEN
344 maxwrk = max( maxwrk, hswork )
345 ELSE
346 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
347 $ 'ZUNGHR',
348 $ ' ', n, 1, n, -1 ) )
349 maxwrk = max( maxwrk, hswork )
350 END IF
351 lwrk = maxwrk
352 IF( .NOT.wantsn )
353 $ lwrk = max( lwrk, ( n*n )/2 )
354 END IF
355 work( 1 ) = lwrk
356*
357 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
358 info = -15
359 END IF
360 END IF
361*
362 IF( info.NE.0 ) THEN
363 CALL xerbla( 'ZGEESX', -info )
364 RETURN
365 ELSE IF( lquery ) THEN
366 RETURN
367 END IF
368*
369* Quick return if possible
370*
371 IF( n.EQ.0 ) THEN
372 sdim = 0
373 RETURN
374 END IF
375*
376* Get machine constants
377*
378 eps = dlamch( 'P' )
379 smlnum = dlamch( 'S' )
380 bignum = one / smlnum
381 smlnum = sqrt( smlnum ) / eps
382 bignum = one / smlnum
383*
384* Scale A if max element outside range [SMLNUM,BIGNUM]
385*
386 anrm = zlange( 'M', n, n, a, lda, dum )
387 scalea = .false.
388 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
389 scalea = .true.
390 cscale = smlnum
391 ELSE IF( anrm.GT.bignum ) THEN
392 scalea = .true.
393 cscale = bignum
394 END IF
395 IF( scalea )
396 $ CALL zlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
397*
398*
399* Permute the matrix to make it more nearly triangular
400* (CWorkspace: none)
401* (RWorkspace: need N)
402*
403 ibal = 1
404 CALL zgebal( 'P', n, a, lda, ilo, ihi, rwork( ibal ), ierr )
405*
406* Reduce to upper Hessenberg form
407* (CWorkspace: need 2*N, prefer N+N*NB)
408* (RWorkspace: none)
409*
410 itau = 1
411 iwrk = n + itau
412 CALL zgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
413 $ lwork-iwrk+1, ierr )
414*
415 IF( wantvs ) THEN
416*
417* Copy Householder vectors to VS
418*
419 CALL zlacpy( 'L', n, n, a, lda, vs, ldvs )
420*
421* Generate unitary matrix in VS
422* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
423* (RWorkspace: none)
424*
425 CALL zunghr( n, ilo, ihi, vs, ldvs, work( itau ),
426 $ work( iwrk ),
427 $ lwork-iwrk+1, ierr )
428 END IF
429*
430 sdim = 0
431*
432* Perform QR iteration, accumulating Schur vectors in VS if desired
433* (CWorkspace: need 1, prefer HSWORK (see comments) )
434* (RWorkspace: none)
435*
436 iwrk = itau
437 CALL zhseqr( 'S', jobvs, n, ilo, ihi, a, lda, w, vs, ldvs,
438 $ work( iwrk ), lwork-iwrk+1, ieval )
439 IF( ieval.GT.0 )
440 $ info = ieval
441*
442* Sort eigenvalues if desired
443*
444 IF( wantst .AND. info.EQ.0 ) THEN
445 IF( scalea )
446 $ CALL zlascl( 'G', 0, 0, cscale, anrm, n, 1, w, n, ierr )
447 DO 10 i = 1, n
448 bwork( i ) = SELECT( w( i ) )
449 10 CONTINUE
450*
451* Reorder eigenvalues, transform Schur vectors, and compute
452* reciprocal condition numbers
453* (CWorkspace: if SENSE is not 'N', need 2*SDIM*(N-SDIM)
454* otherwise, need none )
455* (RWorkspace: none)
456*
457 CALL ztrsen( sense, jobvs, bwork, n, a, lda, vs, ldvs, w,
458 $ sdim,
459 $ rconde, rcondv, work( iwrk ), lwork-iwrk+1,
460 $ icond )
461 IF( .NOT.wantsn )
462 $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
463 IF( icond.EQ.-14 ) THEN
464*
465* Not enough complex workspace
466*
467 info = -15
468 END IF
469 END IF
470*
471 IF( wantvs ) THEN
472*
473* Undo balancing
474* (CWorkspace: none)
475* (RWorkspace: need N)
476*
477 CALL zgebak( 'P', 'R', n, ilo, ihi, rwork( ibal ), n, vs,
478 $ ldvs,
479 $ ierr )
480 END IF
481*
482 IF( scalea ) THEN
483*
484* Undo scaling for the Schur form of A
485*
486 CALL zlascl( 'U', 0, 0, cscale, anrm, n, n, a, lda, ierr )
487 CALL zcopy( n, a, lda+1, w, 1 )
488 IF( ( wantsv .OR. wantsb ) .AND. info.EQ.0 ) THEN
489 dum( 1 ) = rcondv
490 CALL dlascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1,
491 $ ierr )
492 rcondv = dum( 1 )
493 END IF
494 END IF
495*
496 work( 1 ) = maxwrk
497 RETURN
498*
499* End of ZGEESX
500*
501 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
ZGEBAK
Definition zgebak.f:129
subroutine zgebal(job, n, a, lda, ilo, ihi, scale, info)
ZGEBAL
Definition zgebal.f:163
subroutine zgeesx(jobvs, sort, select, sense, n, a, lda, sdim, w, vs, ldvs, rconde, rcondv, work, lwork, rwork, bwork, info)
ZGEESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE ...
Definition zgeesx.f:238
subroutine zgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZGEHRD
Definition zgehrd.f:166
subroutine zhseqr(job, compz, n, ilo, ihi, h, ldh, w, z, ldz, work, lwork, info)
ZHSEQR
Definition zhseqr.f:297
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 zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition zlascl.f:142
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:142
subroutine ztrsen(job, compq, select, n, t, ldt, q, ldq, w, m, s, sep, work, lwork, info)
ZTRSEN
Definition ztrsen.f:263
subroutine zunghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
ZUNGHR
Definition zunghr.f:125