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