LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgges3.f
Go to the documentation of this file.
1*> \brief <b> ZGGES3 computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices (blocked algorithm)</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGGES3 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgges3.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgges3.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgges3.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGGES3( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B,
20* $ LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR,
21* $ WORK, 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* DOUBLE PRECISION RWORK( * )
30* COMPLEX*16 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*> ZGGES3 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*> ZGGEV 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*16 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*16 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*16 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*16 array, dimension (N)
161*> \endverbatim
162*>
163*> \param[out] BETA
164*> \verbatim
165*> BETA is COMPLEX*16 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 ZGGES3. 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*16 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*16 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*16 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 DOUBLE PRECISION 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 ZLAQZ0
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 ZTGSEN.
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 gges3
263*
264* =====================================================================
265 SUBROUTINE zgges3( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B,
266 $ LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR,
267 $ WORK, LWORK, RWORK, BWORK, INFO )
268*
269* -- LAPACK driver routine --
270* -- LAPACK is a software package provided by Univ. of Tennessee, --
271* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
272*
273* .. Scalar Arguments ..
274 CHARACTER JOBVSL, JOBVSR, SORT
275 INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
276* ..
277* .. Array Arguments ..
278 LOGICAL BWORK( * )
279 DOUBLE PRECISION RWORK( * )
280 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
281 $ beta( * ), vsl( ldvsl, * ), vsr( ldvsr, * ),
282 $ work( * )
283* ..
284* .. Function Arguments ..
285 LOGICAL SELCTG
286 EXTERNAL SELCTG
287* ..
288*
289* =====================================================================
290*
291* .. Parameters ..
292 DOUBLE PRECISION ZERO, ONE
293 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
294 COMPLEX*16 CZERO, CONE
295 parameter( czero = ( 0.0d0, 0.0d0 ),
296 $ cone = ( 1.0d0, 0.0d0 ) )
297* ..
298* .. Local Scalars ..
299 LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
300 $ LQUERY, WANTST
301 INTEGER I, ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT,
302 $ ILO, IRIGHT, IROWS, IRWRK, ITAU, IWRK, LWKOPT,
303 $ lwkmin
304 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PVSL,
305 $ PVSR, SMLNUM
306* ..
307* .. Local Arrays ..
308 INTEGER IDUM( 1 )
309 DOUBLE PRECISION DIF( 2 )
310* ..
311* .. External Subroutines ..
312 EXTERNAL xerbla, zgeqrf, zggbak, zggbal, zgghd3,
313 $ zlaqz0,
315* ..
316* .. External Functions ..
317 LOGICAL LSAME
318 DOUBLE PRECISION DLAMCH, ZLANGE
319 EXTERNAL lsame, dlamch, zlange
320* ..
321* .. Intrinsic Functions ..
322 INTRINSIC max, sqrt
323* ..
324* .. Executable Statements ..
325*
326* Decode the input arguments
327*
328 IF( lsame( jobvsl, 'N' ) ) THEN
329 ijobvl = 1
330 ilvsl = .false.
331 ELSE IF( lsame( jobvsl, 'V' ) ) THEN
332 ijobvl = 2
333 ilvsl = .true.
334 ELSE
335 ijobvl = -1
336 ilvsl = .false.
337 END IF
338*
339 IF( lsame( jobvsr, 'N' ) ) THEN
340 ijobvr = 1
341 ilvsr = .false.
342 ELSE IF( lsame( jobvsr, 'V' ) ) THEN
343 ijobvr = 2
344 ilvsr = .true.
345 ELSE
346 ijobvr = -1
347 ilvsr = .false.
348 END IF
349*
350 wantst = lsame( sort, 'S' )
351*
352* Test the input arguments
353*
354 info = 0
355 lquery = ( lwork.EQ.-1 )
356 lwkmin = max( 1, 2*n )
357*
358 IF( ijobvl.LE.0 ) THEN
359 info = -1
360 ELSE IF( ijobvr.LE.0 ) THEN
361 info = -2
362 ELSE IF( ( .NOT.wantst ) .AND.
363 $ ( .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 ELSE IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
376 info = -18
377 END IF
378*
379* Compute workspace
380*
381 IF( info.EQ.0 ) THEN
382 CALL zgeqrf( n, n, b, ldb, work, work, -1, ierr )
383 lwkopt = max( lwkmin, n + int( work( 1 ) ) )
384 CALL zunmqr( 'L', 'C', n, n, n, b, ldb, work, a, lda, work,
385 $ -1, ierr )
386 lwkopt = max( lwkopt, n + int( work( 1 ) ) )
387 IF( ilvsl ) THEN
388 CALL zungqr( n, n, n, vsl, ldvsl, work, work, -1, ierr )
389 lwkopt = max( lwkopt, n + int( work( 1 ) ) )
390 END IF
391 CALL zgghd3( jobvsl, jobvsr, n, 1, n, a, lda, b, ldb, vsl,
392 $ ldvsl, vsr, ldvsr, work, -1, ierr )
393 lwkopt = max( lwkopt, n + int( work( 1 ) ) )
394 CALL zlaqz0( 'S', jobvsl, jobvsr, n, 1, n, a, lda, b, ldb,
395 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work, -1,
396 $ rwork, 0, ierr )
397 lwkopt = max( lwkopt, int( work( 1 ) ) )
398 IF( wantst ) THEN
399 CALL ztgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
400 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, sdim,
401 $ pvsl, pvsr, dif, work, -1, idum, 1, ierr )
402 lwkopt = max( lwkopt, int( work( 1 ) ) )
403 END IF
404 IF( n.EQ.0 ) THEN
405 work( 1 ) = 1
406 ELSE
407 work( 1 ) = dcmplx( lwkopt )
408 END IF
409 END IF
410*
411 IF( info.NE.0 ) THEN
412 CALL xerbla( 'ZGGES3 ', -info )
413 RETURN
414 ELSE IF( lquery ) THEN
415 RETURN
416 END IF
417*
418* Quick return if possible
419*
420 IF( n.EQ.0 ) THEN
421 sdim = 0
422 RETURN
423 END IF
424*
425* Get machine constants
426*
427 eps = dlamch( 'P' )
428 smlnum = dlamch( 'S' )
429 bignum = one / smlnum
430 smlnum = sqrt( smlnum ) / eps
431 bignum = one / smlnum
432*
433* Scale A if max element outside range [SMLNUM,BIGNUM]
434*
435 anrm = zlange( 'M', n, n, a, lda, rwork )
436 ilascl = .false.
437 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
438 anrmto = smlnum
439 ilascl = .true.
440 ELSE IF( anrm.GT.bignum ) THEN
441 anrmto = bignum
442 ilascl = .true.
443 END IF
444*
445 IF( ilascl )
446 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
447*
448* Scale B if max element outside range [SMLNUM,BIGNUM]
449*
450 bnrm = zlange( 'M', n, n, b, ldb, rwork )
451 ilbscl = .false.
452 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
453 bnrmto = smlnum
454 ilbscl = .true.
455 ELSE IF( bnrm.GT.bignum ) THEN
456 bnrmto = bignum
457 ilbscl = .true.
458 END IF
459*
460 IF( ilbscl )
461 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
462*
463* Permute the matrix to make it more nearly triangular
464*
465 ileft = 1
466 iright = n + 1
467 irwrk = iright + n
468 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
469 $ rwork( iright ), rwork( irwrk ), ierr )
470*
471* Reduce B to triangular form (QR decomposition of B)
472*
473 irows = ihi + 1 - ilo
474 icols = n + 1 - ilo
475 itau = 1
476 iwrk = itau + irows
477 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
478 $ work( iwrk ), lwork+1-iwrk, ierr )
479*
480* Apply the orthogonal transformation to matrix A
481*
482 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
483 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
484 $ lwork+1-iwrk, ierr )
485*
486* Initialize VSL
487*
488 IF( ilvsl ) THEN
489 CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
490 IF( irows.GT.1 ) THEN
491 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
492 $ vsl( ilo+1, ilo ), ldvsl )
493 END IF
494 CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
495 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
496 END IF
497*
498* Initialize VSR
499*
500 IF( ilvsr )
501 $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
502*
503* Reduce to generalized Hessenberg form
504*
505 CALL zgghd3( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
506 $ ldvsl, vsr, ldvsr, work( iwrk ), lwork+1-iwrk, ierr )
507*
508 sdim = 0
509*
510* Perform QZ algorithm, computing Schur vectors if desired
511*
512 iwrk = itau
513 CALL zlaqz0( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
514 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwrk ),
515 $ lwork+1-iwrk, rwork( irwrk ), 0, ierr )
516 IF( ierr.NE.0 ) THEN
517 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
518 info = ierr
519 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
520 info = ierr - n
521 ELSE
522 info = n + 1
523 END IF
524 GO TO 30
525 END IF
526*
527* Sort eigenvalues ALPHA/BETA if desired
528*
529 IF( wantst ) THEN
530*
531* Undo scaling on eigenvalues before selecting
532*
533 IF( ilascl )
534 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, 1, alpha, n,
535 $ ierr )
536 IF( ilbscl )
537 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, 1, beta, n,
538 $ ierr )
539*
540* Select eigenvalues
541*
542 DO 10 i = 1, n
543 bwork( i ) = selctg( alpha( i ), beta( i ) )
544 10 CONTINUE
545*
546 CALL ztgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
547 $ alpha,
548 $ beta, vsl, ldvsl, vsr, ldvsr, sdim, pvsl, pvsr,
549 $ dif, work( iwrk ), lwork-iwrk+1, idum, 1, ierr )
550 IF( ierr.EQ.1 )
551 $ info = n + 3
552*
553 END IF
554*
555* Apply back-permutation to VSL and VSR
556*
557 IF( ilvsl )
558 $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
559 $ rwork( iright ), n, vsl, ldvsl, ierr )
560 IF( ilvsr )
561 $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
562 $ rwork( iright ), n, vsr, ldvsr, ierr )
563*
564* Undo scaling
565*
566 IF( ilascl ) THEN
567 CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
568 CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
569 END IF
570*
571 IF( ilbscl ) THEN
572 CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
573 CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
574 END IF
575*
576 IF( wantst ) THEN
577*
578* Check if reordering is correct
579*
580 lastsl = .true.
581 sdim = 0
582 DO 20 i = 1, n
583 cursl = selctg( alpha( i ), beta( i ) )
584 IF( cursl )
585 $ sdim = sdim + 1
586 IF( cursl .AND. .NOT.lastsl )
587 $ info = n + 2
588 lastsl = cursl
589 20 CONTINUE
590*
591 END IF
592*
593 30 CONTINUE
594*
595 work( 1 ) = dcmplx( lwkopt )
596*
597 RETURN
598*
599* End of ZGGES3
600*
601 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqrf(m, n, a, lda, tau, work, lwork, info)
ZGEQRF
Definition zgeqrf.f:144
subroutine zggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
ZGGBAK
Definition zggbak.f:147
subroutine zggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
ZGGBAL
Definition zggbal.f:175
subroutine zgges3(jobvsl, jobvsr, sort, selctg, n, a, lda, b, ldb, sdim, alpha, beta, vsl, ldvsl, vsr, ldvsr, work, lwork, rwork, bwork, info)
ZGGES3 computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE ...
Definition zgges3.f:268
subroutine zgghd3(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, work, lwork, info)
ZGGHD3
Definition zgghd3.f:226
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
recursive subroutine zlaqz0(wants, wantq, wantz, n, ilo, ihi, a, lda, b, ldb, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, rec, info)
ZLAQZ0
Definition zlaqz0.f:283
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 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:104
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:432
subroutine zungqr(m, n, k, a, lda, tau, work, lwork, info)
ZUNGQR
Definition zungqr.f:126
subroutine zunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
ZUNMQR
Definition zunmqr.f:165