LAPACK  3.10.0
LAPACK: Linear Algebra PACKage
zgges.f
Go to the documentation of this file.
1 *> \brief <b> ZGGES 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 *> \htmlonly
9 *> Download ZGGES + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgges.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgges.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgges.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZGGES( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, LDB,
22 * SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, WORK,
23 * LWORK, RWORK, BWORK, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER JOBVSL, JOBVSR, SORT
27 * INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
28 * ..
29 * .. Array Arguments ..
30 * LOGICAL BWORK( * )
31 * DOUBLE PRECISION RWORK( * )
32 * COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
33 * $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ),
34 * $ WORK( * )
35 * ..
36 * .. Function Arguments ..
37 * LOGICAL SELCTG
38 * EXTERNAL SELCTG
39 * ..
40 *
41 *
42 *> \par Purpose:
43 * =============
44 *>
45 *> \verbatim
46 *>
47 *> ZGGES computes for a pair of N-by-N complex nonsymmetric matrices
48 *> (A,B), the generalized eigenvalues, the generalized complex Schur
49 *> form (S, T), and optionally left and/or right Schur vectors (VSL
50 *> and VSR). This gives the generalized Schur factorization
51 *>
52 *> (A,B) = ( (VSL)*S*(VSR)**H, (VSL)*T*(VSR)**H )
53 *>
54 *> where (VSR)**H is the conjugate-transpose of VSR.
55 *>
56 *> Optionally, it also orders the eigenvalues so that a selected cluster
57 *> of eigenvalues appears in the leading diagonal blocks of the upper
58 *> triangular matrix S and the upper triangular matrix T. The leading
59 *> columns of VSL and VSR then form an unitary basis for the
60 *> corresponding left and right eigenspaces (deflating subspaces).
61 *>
62 *> (If only the generalized eigenvalues are needed, use the driver
63 *> ZGGEV instead, which is faster.)
64 *>
65 *> A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
66 *> or a ratio alpha/beta = w, such that A - w*B is singular. It is
67 *> usually represented as the pair (alpha,beta), as there is a
68 *> reasonable interpretation for beta=0, and even for both being zero.
69 *>
70 *> A pair of matrices (S,T) is in generalized complex Schur form if S
71 *> and T are upper triangular and, in addition, the diagonal elements
72 *> of T are non-negative real numbers.
73 *> \endverbatim
74 *
75 * Arguments:
76 * ==========
77 *
78 *> \param[in] JOBVSL
79 *> \verbatim
80 *> JOBVSL is CHARACTER*1
81 *> = 'N': do not compute the left Schur vectors;
82 *> = 'V': compute the left Schur vectors.
83 *> \endverbatim
84 *>
85 *> \param[in] JOBVSR
86 *> \verbatim
87 *> JOBVSR is CHARACTER*1
88 *> = 'N': do not compute the right Schur vectors;
89 *> = 'V': compute the right Schur vectors.
90 *> \endverbatim
91 *>
92 *> \param[in] SORT
93 *> \verbatim
94 *> SORT is CHARACTER*1
95 *> Specifies whether or not to order the eigenvalues on the
96 *> diagonal of the generalized Schur form.
97 *> = 'N': Eigenvalues are not ordered;
98 *> = 'S': Eigenvalues are ordered (see SELCTG).
99 *> \endverbatim
100 *>
101 *> \param[in] SELCTG
102 *> \verbatim
103 *> SELCTG is a LOGICAL FUNCTION of two COMPLEX*16 arguments
104 *> SELCTG must be declared EXTERNAL in the calling subroutine.
105 *> If SORT = 'N', SELCTG is not referenced.
106 *> If SORT = 'S', SELCTG is used to select eigenvalues to sort
107 *> to the top left of the Schur form.
108 *> An eigenvalue ALPHA(j)/BETA(j) is selected if
109 *> SELCTG(ALPHA(j),BETA(j)) is true.
110 *>
111 *> Note that a selected complex eigenvalue may no longer satisfy
112 *> SELCTG(ALPHA(j),BETA(j)) = .TRUE. after ordering, since
113 *> ordering may change the value of complex eigenvalues
114 *> (especially if the eigenvalue is ill-conditioned), in this
115 *> case INFO is set to N+2 (See INFO below).
116 *> \endverbatim
117 *>
118 *> \param[in] N
119 *> \verbatim
120 *> N is INTEGER
121 *> The order of the matrices A, B, VSL, and VSR. N >= 0.
122 *> \endverbatim
123 *>
124 *> \param[in,out] A
125 *> \verbatim
126 *> A is COMPLEX*16 array, dimension (LDA, N)
127 *> On entry, the first of the pair of matrices.
128 *> On exit, A has been overwritten by its generalized Schur
129 *> form S.
130 *> \endverbatim
131 *>
132 *> \param[in] LDA
133 *> \verbatim
134 *> LDA is INTEGER
135 *> The leading dimension of A. LDA >= max(1,N).
136 *> \endverbatim
137 *>
138 *> \param[in,out] B
139 *> \verbatim
140 *> B is COMPLEX*16 array, dimension (LDB, N)
141 *> On entry, the second of the pair of matrices.
142 *> On exit, B has been overwritten by its generalized Schur
143 *> form T.
144 *> \endverbatim
145 *>
146 *> \param[in] LDB
147 *> \verbatim
148 *> LDB is INTEGER
149 *> The leading dimension of B. LDB >= max(1,N).
150 *> \endverbatim
151 *>
152 *> \param[out] SDIM
153 *> \verbatim
154 *> SDIM is INTEGER
155 *> If SORT = 'N', SDIM = 0.
156 *> If SORT = 'S', SDIM = number of eigenvalues (after sorting)
157 *> for which SELCTG is true.
158 *> \endverbatim
159 *>
160 *> \param[out] ALPHA
161 *> \verbatim
162 *> ALPHA is COMPLEX*16 array, dimension (N)
163 *> \endverbatim
164 *>
165 *> \param[out] BETA
166 *> \verbatim
167 *> BETA is COMPLEX*16 array, dimension (N)
168 *> On exit, ALPHA(j)/BETA(j), j=1,...,N, will be the
169 *> generalized eigenvalues. ALPHA(j), j=1,...,N and BETA(j),
170 *> j=1,...,N are the diagonals of the complex Schur form (A,B)
171 *> output by ZGGES. The BETA(j) will be non-negative real.
172 *>
173 *> Note: the quotients ALPHA(j)/BETA(j) may easily over- or
174 *> underflow, and BETA(j) may even be zero. Thus, the user
175 *> should avoid naively computing the ratio alpha/beta.
176 *> However, ALPHA will be always less than and usually
177 *> comparable with norm(A) in magnitude, and BETA always less
178 *> than and usually comparable with norm(B).
179 *> \endverbatim
180 *>
181 *> \param[out] VSL
182 *> \verbatim
183 *> VSL is COMPLEX*16 array, dimension (LDVSL,N)
184 *> If JOBVSL = 'V', VSL will contain the left Schur vectors.
185 *> Not referenced if JOBVSL = 'N'.
186 *> \endverbatim
187 *>
188 *> \param[in] LDVSL
189 *> \verbatim
190 *> LDVSL is INTEGER
191 *> The leading dimension of the matrix VSL. LDVSL >= 1, and
192 *> if JOBVSL = 'V', LDVSL >= N.
193 *> \endverbatim
194 *>
195 *> \param[out] VSR
196 *> \verbatim
197 *> VSR is COMPLEX*16 array, dimension (LDVSR,N)
198 *> If JOBVSR = 'V', VSR will contain the right Schur vectors.
199 *> Not referenced if JOBVSR = 'N'.
200 *> \endverbatim
201 *>
202 *> \param[in] LDVSR
203 *> \verbatim
204 *> LDVSR is INTEGER
205 *> The leading dimension of the matrix VSR. LDVSR >= 1, and
206 *> if JOBVSR = 'V', LDVSR >= N.
207 *> \endverbatim
208 *>
209 *> \param[out] WORK
210 *> \verbatim
211 *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
212 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
213 *> \endverbatim
214 *>
215 *> \param[in] LWORK
216 *> \verbatim
217 *> LWORK is INTEGER
218 *> The dimension of the array WORK. LWORK >= max(1,2*N).
219 *> For good performance, LWORK must generally be larger.
220 *>
221 *> If LWORK = -1, then a workspace query is assumed; the routine
222 *> only calculates the optimal size of the WORK array, returns
223 *> this value as the first entry of the WORK array, and no error
224 *> message related to LWORK is issued by XERBLA.
225 *> \endverbatim
226 *>
227 *> \param[out] RWORK
228 *> \verbatim
229 *> RWORK is DOUBLE PRECISION array, dimension (8*N)
230 *> \endverbatim
231 *>
232 *> \param[out] BWORK
233 *> \verbatim
234 *> BWORK is LOGICAL array, dimension (N)
235 *> Not referenced if SORT = 'N'.
236 *> \endverbatim
237 *>
238 *> \param[out] INFO
239 *> \verbatim
240 *> INFO is INTEGER
241 *> = 0: successful exit
242 *> < 0: if INFO = -i, the i-th argument had an illegal value.
243 *> =1,...,N:
244 *> The QZ iteration failed. (A,B) are not in Schur
245 *> form, but ALPHA(j) and BETA(j) should be correct for
246 *> j=INFO+1,...,N.
247 *> > N: =N+1: other than QZ iteration failed in ZHGEQZ
248 *> =N+2: after reordering, roundoff changed values of
249 *> some complex eigenvalues so that leading
250 *> eigenvalues in the Generalized Schur form no
251 *> longer satisfy SELCTG=.TRUE. This could also
252 *> be caused due to scaling.
253 *> =N+3: reordering failed in ZTGSEN.
254 *> \endverbatim
255 *
256 * Authors:
257 * ========
258 *
259 *> \author Univ. of Tennessee
260 *> \author Univ. of California Berkeley
261 *> \author Univ. of Colorado Denver
262 *> \author NAG Ltd.
263 *
264 *> \ingroup complex16GEeigen
265 *
266 * =====================================================================
267  SUBROUTINE zgges( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, LDB,
268  $ SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, WORK,
269  $ LWORK, RWORK, BWORK, INFO )
270 *
271 * -- LAPACK driver routine --
272 * -- LAPACK is a software package provided by Univ. of Tennessee, --
273 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
274 *
275 * .. Scalar Arguments ..
276  CHARACTER JOBVSL, JOBVSR, SORT
277  INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
278 * ..
279 * .. Array Arguments ..
280  LOGICAL BWORK( * )
281  DOUBLE PRECISION RWORK( * )
282  COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
283  $ beta( * ), vsl( ldvsl, * ), vsr( ldvsr, * ),
284  $ work( * )
285 * ..
286 * .. Function Arguments ..
287  LOGICAL SELCTG
288  EXTERNAL SELCTG
289 * ..
290 *
291 * =====================================================================
292 *
293 * .. Parameters ..
294  DOUBLE PRECISION ZERO, ONE
295  PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
296  COMPLEX*16 CZERO, CONE
297  parameter( czero = ( 0.0d0, 0.0d0 ),
298  $ cone = ( 1.0d0, 0.0d0 ) )
299 * ..
300 * .. Local Scalars ..
301  LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
302  $ LQUERY, WANTST
303  INTEGER I, ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT,
304  $ ILO, IRIGHT, IROWS, IRWRK, ITAU, IWRK, LWKMIN,
305  $ lwkopt
306  DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PVSL,
307  $ PVSR, SMLNUM
308 * ..
309 * .. Local Arrays ..
310  INTEGER IDUM( 1 )
311  DOUBLE PRECISION DIF( 2 )
312 * ..
313 * .. External Subroutines ..
314  EXTERNAL dlabad, xerbla, zgeqrf, zggbak, zggbal, zgghrd,
316  $ zunmqr
317 * ..
318 * .. External Functions ..
319  LOGICAL LSAME
320  INTEGER ILAENV
321  DOUBLE PRECISION DLAMCH, ZLANGE
322  EXTERNAL lsame, ilaenv, dlamch, zlange
323 * ..
324 * .. Intrinsic Functions ..
325  INTRINSIC max, sqrt
326 * ..
327 * .. Executable Statements ..
328 *
329 * Decode the input arguments
330 *
331  IF( lsame( jobvsl, 'N' ) ) THEN
332  ijobvl = 1
333  ilvsl = .false.
334  ELSE IF( lsame( jobvsl, 'V' ) ) THEN
335  ijobvl = 2
336  ilvsl = .true.
337  ELSE
338  ijobvl = -1
339  ilvsl = .false.
340  END IF
341 *
342  IF( lsame( jobvsr, 'N' ) ) THEN
343  ijobvr = 1
344  ilvsr = .false.
345  ELSE IF( lsame( jobvsr, 'V' ) ) THEN
346  ijobvr = 2
347  ilvsr = .true.
348  ELSE
349  ijobvr = -1
350  ilvsr = .false.
351  END IF
352 *
353  wantst = lsame( sort, 'S' )
354 *
355 * Test the input arguments
356 *
357  info = 0
358  lquery = ( lwork.EQ.-1 )
359  IF( ijobvl.LE.0 ) THEN
360  info = -1
361  ELSE IF( ijobvr.LE.0 ) THEN
362  info = -2
363  ELSE IF( ( .NOT.wantst ) .AND. ( .NOT.lsame( sort, 'N' ) ) ) THEN
364  info = -3
365  ELSE IF( n.LT.0 ) THEN
366  info = -5
367  ELSE IF( lda.LT.max( 1, n ) ) THEN
368  info = -7
369  ELSE IF( ldb.LT.max( 1, n ) ) THEN
370  info = -9
371  ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
372  info = -14
373  ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
374  info = -16
375  END IF
376 *
377 * Compute workspace
378 * (Note: Comments in the code beginning "Workspace:" describe the
379 * minimal amount of workspace needed at that point in the code,
380 * as well as the preferred amount for good performance.
381 * NB refers to the optimal block size for the immediately
382 * following subroutine, as returned by ILAENV.)
383 *
384  IF( info.EQ.0 ) THEN
385  lwkmin = max( 1, 2*n )
386  lwkopt = max( 1, n + n*ilaenv( 1, 'ZGEQRF', ' ', n, 1, n, 0 ) )
387  lwkopt = max( lwkopt, n +
388  $ n*ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, -1 ) )
389  IF( ilvsl ) THEN
390  lwkopt = max( lwkopt, n +
391  $ n*ilaenv( 1, 'ZUNGQR', ' ', n, 1, n, -1 ) )
392  END IF
393  work( 1 ) = lwkopt
394 *
395  IF( lwork.LT.lwkmin .AND. .NOT.lquery )
396  $ info = -18
397  END IF
398 *
399  IF( info.NE.0 ) THEN
400  CALL xerbla( 'ZGGES ', -info )
401  RETURN
402  ELSE IF( lquery ) THEN
403  RETURN
404  END IF
405 *
406 * Quick return if possible
407 *
408  IF( n.EQ.0 ) THEN
409  sdim = 0
410  RETURN
411  END IF
412 *
413 * Get machine constants
414 *
415  eps = dlamch( 'P' )
416  smlnum = dlamch( 'S' )
417  bignum = one / smlnum
418  CALL dlabad( smlnum, bignum )
419  smlnum = sqrt( smlnum ) / eps
420  bignum = one / smlnum
421 *
422 * Scale A if max element outside range [SMLNUM,BIGNUM]
423 *
424  anrm = zlange( 'M', n, n, a, lda, rwork )
425  ilascl = .false.
426  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
427  anrmto = smlnum
428  ilascl = .true.
429  ELSE IF( anrm.GT.bignum ) THEN
430  anrmto = bignum
431  ilascl = .true.
432  END IF
433 *
434  IF( ilascl )
435  $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
436 *
437 * Scale B if max element outside range [SMLNUM,BIGNUM]
438 *
439  bnrm = zlange( 'M', n, n, b, ldb, rwork )
440  ilbscl = .false.
441  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
442  bnrmto = smlnum
443  ilbscl = .true.
444  ELSE IF( bnrm.GT.bignum ) THEN
445  bnrmto = bignum
446  ilbscl = .true.
447  END IF
448 *
449  IF( ilbscl )
450  $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
451 *
452 * Permute the matrix to make it more nearly triangular
453 * (Real Workspace: need 6*N)
454 *
455  ileft = 1
456  iright = n + 1
457  irwrk = iright + n
458  CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
459  $ rwork( iright ), rwork( irwrk ), ierr )
460 *
461 * Reduce B to triangular form (QR decomposition of B)
462 * (Complex Workspace: need N, prefer N*NB)
463 *
464  irows = ihi + 1 - ilo
465  icols = n + 1 - ilo
466  itau = 1
467  iwrk = itau + irows
468  CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
469  $ work( iwrk ), lwork+1-iwrk, ierr )
470 *
471 * Apply the orthogonal transformation to matrix A
472 * (Complex Workspace: need N, prefer N*NB)
473 *
474  CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
475  $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
476  $ lwork+1-iwrk, ierr )
477 *
478 * Initialize VSL
479 * (Complex Workspace: need N, prefer N*NB)
480 *
481  IF( ilvsl ) THEN
482  CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
483  IF( irows.GT.1 ) THEN
484  CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
485  $ vsl( ilo+1, ilo ), ldvsl )
486  END IF
487  CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
488  $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
489  END IF
490 *
491 * Initialize VSR
492 *
493  IF( ilvsr )
494  $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
495 *
496 * Reduce to generalized Hessenberg form
497 * (Workspace: none needed)
498 *
499  CALL zgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
500  $ ldvsl, vsr, ldvsr, ierr )
501 *
502  sdim = 0
503 *
504 * Perform QZ algorithm, computing Schur vectors if desired
505 * (Complex Workspace: need N)
506 * (Real Workspace: need N)
507 *
508  iwrk = itau
509  CALL zhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
510  $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwrk ),
511  $ lwork+1-iwrk, rwork( irwrk ), ierr )
512  IF( ierr.NE.0 ) THEN
513  IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
514  info = ierr
515  ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
516  info = ierr - n
517  ELSE
518  info = n + 1
519  END IF
520  GO TO 30
521  END IF
522 *
523 * Sort eigenvalues ALPHA/BETA if desired
524 * (Workspace: none needed)
525 *
526  IF( wantst ) THEN
527 *
528 * Undo scaling on eigenvalues before selecting
529 *
530  IF( ilascl )
531  $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, 1, alpha, n, ierr )
532  IF( ilbscl )
533  $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, 1, beta, n, ierr )
534 *
535 * Select eigenvalues
536 *
537  DO 10 i = 1, n
538  bwork( i ) = selctg( alpha( i ), beta( i ) )
539  10 CONTINUE
540 *
541  CALL ztgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb, alpha,
542  $ beta, vsl, ldvsl, vsr, ldvsr, sdim, pvsl, pvsr,
543  $ dif, work( iwrk ), lwork-iwrk+1, idum, 1, ierr )
544  IF( ierr.EQ.1 )
545  $ info = n + 3
546 *
547  END IF
548 *
549 * Apply back-permutation to VSL and VSR
550 * (Workspace: none needed)
551 *
552  IF( ilvsl )
553  $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
554  $ rwork( iright ), n, vsl, ldvsl, ierr )
555  IF( ilvsr )
556  $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
557  $ rwork( iright ), n, vsr, ldvsr, ierr )
558 *
559 * Undo scaling
560 *
561  IF( ilascl ) THEN
562  CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
563  CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
564  END IF
565 *
566  IF( ilbscl ) THEN
567  CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
568  CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
569  END IF
570 *
571  IF( wantst ) THEN
572 *
573 * Check if reordering is correct
574 *
575  lastsl = .true.
576  sdim = 0
577  DO 20 i = 1, n
578  cursl = selctg( alpha( i ), beta( i ) )
579  IF( cursl )
580  $ sdim = sdim + 1
581  IF( cursl .AND. .NOT.lastsl )
582  $ info = n + 2
583  lastsl = cursl
584  20 CONTINUE
585 *
586  END IF
587 *
588  30 CONTINUE
589 *
590  work( 1 ) = lwkopt
591 *
592  RETURN
593 *
594 * End of ZGGES
595 *
596  END
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine zggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
ZGGBAL
Definition: zggbal.f:177
subroutine zggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
ZGGBAK
Definition: zggbak.f:148
subroutine zhgeqz(JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, RWORK, INFO)
ZHGEQZ
Definition: zhgeqz.f:284
subroutine zgges(JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, WORK, LWORK, RWORK, BWORK, INFO)
ZGGES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE m...
Definition: zgges.f:270
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:143
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:103
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition: zlaset.f:106
subroutine ztgsen(IJOB, WANTQ, WANTZ, SELECT, N, A, LDA, B, LDB, ALPHA, BETA, Q, LDQ, Z, LDZ, M, PL, PR, DIF, WORK, LWORK, IWORK, LIWORK, INFO)
ZTGSEN
Definition: ztgsen.f:433
subroutine zungqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
ZUNGQR
Definition: zungqr.f:128
subroutine zgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
ZGGHRD
Definition: zgghrd.f:204
subroutine zunmqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMQR
Definition: zunmqr.f:167
subroutine zgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGEQRF VARIANT: left-looking Level 3 BLAS of the algorithm.
Definition: zgeqrf.f:151