LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zgegs()

subroutine zgegs ( character jobvsl,
character jobvsr,
integer n,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( * ) alpha,
complex*16, dimension( * ) beta,
complex*16, dimension( ldvsl, * ) vsl,
integer ldvsl,
complex*16, dimension( ldvsr, * ) vsr,
integer ldvsr,
complex*16, dimension( * ) work,
integer lwork,
double precision, dimension( * ) rwork,
integer info )

ZGEGS computes the eigenvalues, Schur form, and, optionally, the left and or/right Schur vectors of a complex matrix pair (A,B)

Download ZGEGS + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> This routine is deprecated and has been replaced by routine ZGGES.
!>
!> ZGEGS computes the eigenvalues, Schur form, and, optionally, the
!> left and or/right Schur vectors of a complex matrix pair (A,B).
!> Given two square matrices A and B, the generalized Schur
!> factorization has the form
!>
!>    A = Q*S*Z**H,  B = Q*T*Z**H
!>
!> where Q and Z are unitary matrices and S and T are upper triangular.
!> The columns of Q are the left Schur vectors
!> and the columns of Z are the right Schur vectors.
!>
!> If only the eigenvalues of (A,B) are needed, the driver routine
!> ZGEGV should be used instead.  See ZGEGV for a description of the
!> eigenvalues of the generalized nonsymmetric eigenvalue problem
!> (GNEP).
!> 
Parameters
[in]JOBVSL
!>          JOBVSL is CHARACTER*1
!>          = 'N':  do not compute the left Schur vectors;
!>          = 'V':  compute the left Schur vectors (returned in VSL).
!> 
[in]JOBVSR
!>          JOBVSR is CHARACTER*1
!>          = 'N':  do not compute the right Schur vectors;
!>          = 'V':  compute the right Schur vectors (returned in VSR).
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices A, B, VSL, and VSR.  N >= 0.
!> 
[in,out]A
!>          A is COMPLEX*16 array, dimension (LDA, N)
!>          On entry, the matrix A.
!>          On exit, the upper triangular matrix S from the generalized
!>          Schur factorization.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of A.  LDA >= max(1,N).
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension (LDB, N)
!>          On entry, the matrix B.
!>          On exit, the upper triangular matrix T from the generalized
!>          Schur factorization.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of B.  LDB >= max(1,N).
!> 
[out]ALPHA
!>          ALPHA is COMPLEX*16 array, dimension (N)
!>          The complex scalars alpha that define the eigenvalues of
!>          GNEP.  ALPHA(j) = S(j,j), the diagonal element of the Schur
!>          form of A.
!> 
[out]BETA
!>          BETA is COMPLEX*16 array, dimension (N)
!>          The non-negative real scalars beta that define the
!>          eigenvalues of GNEP.  BETA(j) = T(j,j), the diagonal element
!>          of the triangular factor T.
!>
!>          Together, the quantities alpha = ALPHA(j) and beta = BETA(j)
!>          represent the j-th eigenvalue of the matrix pair (A,B), in
!>          one of the forms lambda = alpha/beta or mu = beta/alpha.
!>          Since either lambda or mu may overflow, they should not,
!>          in general, be computed.
!> 
[out]VSL
!>          VSL is COMPLEX*16 array, dimension (LDVSL,N)
!>          If JOBVSL = 'V', the matrix of left Schur vectors Q.
!>          Not referenced if JOBVSL = 'N'.
!> 
[in]LDVSL
!>          LDVSL is INTEGER
!>          The leading dimension of the matrix VSL. LDVSL >= 1, and
!>          if JOBVSL = 'V', LDVSL >= N.
!> 
[out]VSR
!>          VSR is COMPLEX*16 array, dimension (LDVSR,N)
!>          If JOBVSR = 'V', the matrix of right Schur vectors Z.
!>          Not referenced if JOBVSR = 'N'.
!> 
[in]LDVSR
!>          LDVSR is INTEGER
!>          The leading dimension of the matrix VSR. LDVSR >= 1, and
!>          if JOBVSR = 'V', LDVSR >= N.
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.  LWORK >= max(1,2*N).
!>          For good performance, LWORK must generally be larger.
!>          To compute the optimal value of LWORK, call ILAENV to get
!>          blocksizes (for ZGEQRF, ZUNMQR, and CUNGQR.)  Then compute:
!>          NB  -- MAX of the blocksizes for ZGEQRF, ZUNMQR, and CUNGQR;
!>          the optimal LWORK is N*(NB+1).
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (3*N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!>          =1,...,N:
!>                The QZ iteration failed.  (A,B) are not in Schur
!>                form, but ALPHA(j) and BETA(j) should be correct for
!>                j=INFO+1,...,N.
!>          > N:  errors that usually indicate LAPACK problems:
!>                =N+1: error return from ZGGBAL
!>                =N+2: error return from ZGEQRF
!>                =N+3: error return from ZUNMQR
!>                =N+4: error return from ZUNGQR
!>                =N+5: error return from ZGGHRD
!>                =N+6: error return from ZHGEQZ (other than failed
!>                                               iteration)
!>                =N+7: error return from ZGGBAK (computing VSL)
!>                =N+8: error return from ZGGBAK (computing VSR)
!>                =N+9: error return from ZLASCL (various places)
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 220 of file zgegs.f.

223*
224* -- LAPACK driver routine --
225* -- LAPACK is a software package provided by Univ. of Tennessee, --
226* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
227*
228* .. Scalar Arguments ..
229 CHARACTER JOBVSL, JOBVSR
230 INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N
231* ..
232* .. Array Arguments ..
233 DOUBLE PRECISION RWORK( * )
234 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
235 $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ),
236 $ WORK( * )
237* ..
238*
239* =====================================================================
240*
241* .. Parameters ..
242 DOUBLE PRECISION ZERO, ONE
243 parameter( zero = 0.0d0, one = 1.0d0 )
244 COMPLEX*16 CZERO, CONE
245 parameter( czero = ( 0.0d0, 0.0d0 ),
246 $ cone = ( 1.0d0, 0.0d0 ) )
247* ..
248* .. Local Scalars ..
249 LOGICAL ILASCL, ILBSCL, ILVSL, ILVSR, LQUERY
250 INTEGER ICOLS, IHI, IINFO, IJOBVL, IJOBVR, ILEFT, ILO,
251 $ IRIGHT, IROWS, IRWORK, ITAU, IWORK, LOPT,
252 $ LWKMIN, LWKOPT, NB, NB1, NB2, NB3
253 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
254 $ SAFMIN, SMLNUM
255* ..
256* .. External Subroutines ..
257 EXTERNAL xerbla, zgeqrf, zggbak, zggbal, zgghrd, zhgeqz,
259* ..
260* .. External Functions ..
261 LOGICAL LSAME
262 INTEGER ILAENV
263 DOUBLE PRECISION DLAMCH, ZLANGE
264 EXTERNAL lsame, ilaenv, dlamch, zlange
265* ..
266* .. Intrinsic Functions ..
267 INTRINSIC int, max
268* ..
269* .. Executable Statements ..
270*
271* Decode the input arguments
272*
273 IF( lsame( jobvsl, 'N' ) ) THEN
274 ijobvl = 1
275 ilvsl = .false.
276 ELSE IF( lsame( jobvsl, 'V' ) ) THEN
277 ijobvl = 2
278 ilvsl = .true.
279 ELSE
280 ijobvl = -1
281 ilvsl = .false.
282 END IF
283*
284 IF( lsame( jobvsr, 'N' ) ) THEN
285 ijobvr = 1
286 ilvsr = .false.
287 ELSE IF( lsame( jobvsr, 'V' ) ) THEN
288 ijobvr = 2
289 ilvsr = .true.
290 ELSE
291 ijobvr = -1
292 ilvsr = .false.
293 END IF
294*
295* Test the input arguments
296*
297 lwkmin = max( 2*n, 1 )
298 lwkopt = lwkmin
299 work( 1 ) = lwkopt
300 lquery = ( lwork.EQ.-1 )
301 info = 0
302 IF( ijobvl.LE.0 ) THEN
303 info = -1
304 ELSE IF( ijobvr.LE.0 ) THEN
305 info = -2
306 ELSE IF( n.LT.0 ) THEN
307 info = -3
308 ELSE IF( lda.LT.max( 1, n ) ) THEN
309 info = -5
310 ELSE IF( ldb.LT.max( 1, n ) ) THEN
311 info = -7
312 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
313 info = -11
314 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
315 info = -13
316 ELSE IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
317 info = -15
318 END IF
319*
320 IF( info.EQ.0 ) THEN
321 nb1 = ilaenv( 1, 'ZGEQRF', ' ', n, n, -1, -1 )
322 nb2 = ilaenv( 1, 'ZUNMQR', ' ', n, n, n, -1 )
323 nb3 = ilaenv( 1, 'ZUNGQR', ' ', n, n, n, -1 )
324 nb = max( nb1, nb2, nb3 )
325 lopt = n*( nb+1 )
326 work( 1 ) = lopt
327 END IF
328*
329 IF( info.NE.0 ) THEN
330 CALL xerbla( 'ZGEGS ', -info )
331 RETURN
332 ELSE IF( lquery ) THEN
333 RETURN
334 END IF
335*
336* Quick return if possible
337*
338 IF( n.EQ.0 )
339 $ RETURN
340*
341* Get machine constants
342*
343 eps = dlamch( 'E' )*dlamch( 'B' )
344 safmin = dlamch( 'S' )
345 smlnum = n*safmin / eps
346 bignum = one / smlnum
347*
348* Scale A if max element outside range [SMLNUM,BIGNUM]
349*
350 anrm = zlange( 'M', n, n, a, lda, rwork )
351 ilascl = .false.
352 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
353 anrmto = smlnum
354 ilascl = .true.
355 ELSE IF( anrm.GT.bignum ) THEN
356 anrmto = bignum
357 ilascl = .true.
358 END IF
359*
360 IF( ilascl ) THEN
361 CALL zlascl( 'G', -1, -1, anrm, anrmto, n, n, a, lda,
362 $ iinfo )
363 IF( iinfo.NE.0 ) THEN
364 info = n + 9
365 RETURN
366 END IF
367 END IF
368*
369* Scale B if max element outside range [SMLNUM,BIGNUM]
370*
371 bnrm = zlange( 'M', n, n, b, ldb, rwork )
372 ilbscl = .false.
373 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
374 bnrmto = smlnum
375 ilbscl = .true.
376 ELSE IF( bnrm.GT.bignum ) THEN
377 bnrmto = bignum
378 ilbscl = .true.
379 END IF
380*
381 IF( ilbscl ) THEN
382 CALL zlascl( 'G', -1, -1, bnrm, bnrmto, n, n, b, ldb,
383 $ iinfo )
384 IF( iinfo.NE.0 ) THEN
385 info = n + 9
386 RETURN
387 END IF
388 END IF
389*
390* Permute the matrix to make it more nearly triangular
391*
392 ileft = 1
393 iright = n + 1
394 irwork = iright + n
395 iwork = 1
396 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
397 $ rwork( iright ), rwork( irwork ), iinfo )
398 IF( iinfo.NE.0 ) THEN
399 info = n + 1
400 GO TO 10
401 END IF
402*
403* Reduce B to triangular form, and initialize VSL and/or VSR
404*
405 irows = ihi + 1 - ilo
406 icols = n + 1 - ilo
407 itau = iwork
408 iwork = itau + irows
409 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
410 $ work( iwork ), lwork+1-iwork, iinfo )
411 IF( iinfo.GE.0 )
412 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
413 IF( iinfo.NE.0 ) THEN
414 info = n + 2
415 GO TO 10
416 END IF
417*
418 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
419 $ work( itau ), a( ilo, ilo ), lda, work( iwork ),
420 $ lwork+1-iwork, iinfo )
421 IF( iinfo.GE.0 )
422 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
423 IF( iinfo.NE.0 ) THEN
424 info = n + 3
425 GO TO 10
426 END IF
427*
428 IF( ilvsl ) THEN
429 CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
430 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
431 $ vsl( ilo+1, ilo ), ldvsl )
432 CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
433 $ work( itau ), work( iwork ), lwork+1-iwork,
434 $ iinfo )
435 IF( iinfo.GE.0 )
436 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
437 IF( iinfo.NE.0 ) THEN
438 info = n + 4
439 GO TO 10
440 END IF
441 END IF
442*
443 IF( ilvsr )
444 $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
445*
446* Reduce to generalized Hessenberg form
447*
448 CALL zgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
449 $ ldvsl, vsr, ldvsr, iinfo )
450 IF( iinfo.NE.0 ) THEN
451 info = n + 5
452 GO TO 10
453 END IF
454*
455* Perform QZ algorithm, computing Schur vectors if desired
456*
457 iwork = itau
458 CALL zhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
459 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwork ),
460 $ lwork+1-iwork, rwork( irwork ), iinfo )
461 IF( iinfo.GE.0 )
462 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
463 IF( iinfo.NE.0 ) THEN
464 IF( iinfo.GT.0 .AND. iinfo.LE.n ) THEN
465 info = iinfo
466 ELSE IF( iinfo.GT.n .AND. iinfo.LE.2*n ) THEN
467 info = iinfo - n
468 ELSE
469 info = n + 6
470 END IF
471 GO TO 10
472 END IF
473*
474* Apply permutation to VSL and VSR
475*
476 IF( ilvsl ) THEN
477 CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
478 $ rwork( iright ), n, vsl, ldvsl, iinfo )
479 IF( iinfo.NE.0 ) THEN
480 info = n + 7
481 GO TO 10
482 END IF
483 END IF
484 IF( ilvsr ) THEN
485 CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
486 $ rwork( iright ), n, vsr, ldvsr, iinfo )
487 IF( iinfo.NE.0 ) THEN
488 info = n + 8
489 GO TO 10
490 END IF
491 END IF
492*
493* Undo scaling
494*
495 IF( ilascl ) THEN
496 CALL zlascl( 'U', -1, -1, anrmto, anrm, n, n, a, lda,
497 $ iinfo )
498 IF( iinfo.NE.0 ) THEN
499 info = n + 9
500 RETURN
501 END IF
502 CALL zlascl( 'G', -1, -1, anrmto, anrm, n, 1, alpha, n,
503 $ iinfo )
504 IF( iinfo.NE.0 ) THEN
505 info = n + 9
506 RETURN
507 END IF
508 END IF
509*
510 IF( ilbscl ) THEN
511 CALL zlascl( 'U', -1, -1, bnrmto, bnrm, n, n, b, ldb,
512 $ iinfo )
513 IF( iinfo.NE.0 ) THEN
514 info = n + 9
515 RETURN
516 END IF
517 CALL zlascl( 'G', -1, -1, bnrmto, bnrm, n, 1, beta, n,
518 $ iinfo )
519 IF( iinfo.NE.0 ) THEN
520 info = n + 9
521 RETURN
522 END IF
523 END IF
524*
525 10 CONTINUE
526 work( 1 ) = lwkopt
527*
528 RETURN
529*
530* End of ZGEGS
531*
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 zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
ZGGHRD
Definition zgghrd.f:203
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:283
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:113
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
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
Here is the call graph for this function: