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

◆ dgges()

subroutine dgges ( character jobvsl,
character jobvsr,
character sort,
external selctg,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldb, * ) b,
integer ldb,
integer sdim,
double precision, dimension( * ) alphar,
double precision, dimension( * ) alphai,
double precision, dimension( * ) beta,
double precision, dimension( ldvsl, * ) vsl,
integer ldvsl,
double precision, dimension( ldvsr, * ) vsr,
integer ldvsr,
double precision, dimension( * ) work,
integer lwork,
logical, dimension( * ) bwork,
integer info )

DGGES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices

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

Purpose:
!>
!> DGGES computes for a pair of N-by-N real nonsymmetric matrices (A,B),
!> the generalized eigenvalues, the generalized real Schur form (S,T),
!> optionally, the left and/or right matrices of Schur vectors (VSL and
!> VSR). This gives the generalized Schur factorization
!>
!>          (A,B) = ( (VSL)*S*(VSR)**T, (VSL)*T*(VSR)**T )
!>
!> Optionally, it also orders the eigenvalues so that a selected cluster
!> of eigenvalues appears in the leading diagonal blocks of the upper
!> quasi-triangular matrix S and the upper triangular matrix T.The
!> leading columns of VSL and VSR then form an orthonormal basis for the
!> corresponding left and right eigenspaces (deflating subspaces).
!>
!> (If only the generalized eigenvalues are needed, use the driver
!> DGGEV instead, which is faster.)
!>
!> A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
!> or a ratio alpha/beta = w, such that  A - w*B is singular.  It is
!> usually represented as the pair (alpha,beta), as there is a
!> reasonable interpretation for beta=0 or both being zero.
!>
!> A pair of matrices (S,T) is in generalized real Schur form if T is
!> upper triangular with non-negative diagonal and S is block upper
!> triangular with 1-by-1 and 2-by-2 blocks.  1-by-1 blocks correspond
!> to real generalized eigenvalues, while 2-by-2 blocks of S will be
!>  by making the corresponding elements of T have the
!> form:
!>         [  a  0  ]
!>         [  0  b  ]
!>
!> and the pair of corresponding 2-by-2 blocks in S and T will have a
!> complex conjugate pair of generalized eigenvalues.
!>
!> 
Parameters
[in]JOBVSL
!>          JOBVSL is CHARACTER*1
!>          = 'N':  do not compute the left Schur vectors;
!>          = 'V':  compute the left Schur vectors.
!> 
[in]JOBVSR
!>          JOBVSR is CHARACTER*1
!>          = 'N':  do not compute the right Schur vectors;
!>          = 'V':  compute the right Schur vectors.
!> 
[in]SORT
!>          SORT is CHARACTER*1
!>          Specifies whether or not to order the eigenvalues on the
!>          diagonal of the generalized Schur form.
!>          = 'N':  Eigenvalues are not ordered;
!>          = 'S':  Eigenvalues are ordered (see SELCTG);
!> 
[in]SELCTG
!>          SELCTG is a LOGICAL FUNCTION of three DOUBLE PRECISION arguments
!>          SELCTG must be declared EXTERNAL in the calling subroutine.
!>          If SORT = 'N', SELCTG is not referenced.
!>          If SORT = 'S', SELCTG is used to select eigenvalues to sort
!>          to the top left of the Schur form.
!>          An eigenvalue (ALPHAR(j)+ALPHAI(j))/BETA(j) is selected if
!>          SELCTG(ALPHAR(j),ALPHAI(j),BETA(j)) is true; i.e. if either
!>          one of a complex conjugate pair of eigenvalues is selected,
!>          then both complex eigenvalues are selected.
!>
!>          Note that in the ill-conditioned case, a selected complex
!>          eigenvalue may no longer satisfy SELCTG(ALPHAR(j),ALPHAI(j),
!>          BETA(j)) = .TRUE. after ordering. INFO is to be set to N+2
!>          in this case.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices A, B, VSL, and VSR.  N >= 0.
!> 
[in,out]A
!>          A is DOUBLE PRECISION array, dimension (LDA, N)
!>          On entry, the first of the pair of matrices.
!>          On exit, A has been overwritten by its generalized Schur
!>          form S.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of A.  LDA >= max(1,N).
!> 
[in,out]B
!>          B is DOUBLE PRECISION array, dimension (LDB, N)
!>          On entry, the second of the pair of matrices.
!>          On exit, B has been overwritten by its generalized Schur
!>          form T.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of B.  LDB >= max(1,N).
!> 
[out]SDIM
!>          SDIM is INTEGER
!>          If SORT = 'N', SDIM = 0.
!>          If SORT = 'S', SDIM = number of eigenvalues (after sorting)
!>          for which SELCTG is true.  (Complex conjugate pairs for which
!>          SELCTG is true for either eigenvalue count as 2.)
!> 
[out]ALPHAR
!>          ALPHAR is DOUBLE PRECISION array, dimension (N)
!> 
[out]ALPHAI
!>          ALPHAI is DOUBLE PRECISION array, dimension (N)
!> 
[out]BETA
!>          BETA is DOUBLE PRECISION array, dimension (N)
!>          On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
!>          be the generalized eigenvalues.  ALPHAR(j) + ALPHAI(j)*i,
!>          and  BETA(j),j=1,...,N are the diagonals of the complex Schur
!>          form (S,T) that would result if the 2-by-2 diagonal blocks of
!>          the real Schur form of (A,B) were further reduced to
!>          triangular form using 2-by-2 complex unitary transformations.
!>          If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
!>          positive, then the j-th and (j+1)-st eigenvalues are a
!>          complex conjugate pair, with ALPHAI(j+1) negative.
!>
!>          Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j)
!>          may easily over- or underflow, and BETA(j) may even be zero.
!>          Thus, the user should avoid naively computing the ratio.
!>          However, ALPHAR and ALPHAI will be always less than and
!>          usually comparable with norm(A) in magnitude, and BETA always
!>          less than and usually comparable with norm(B).
!> 
[out]VSL
!>          VSL is DOUBLE PRECISION array, dimension (LDVSL,N)
!>          If JOBVSL = 'V', VSL will contain the left Schur vectors.
!>          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 DOUBLE PRECISION array, dimension (LDVSR,N)
!>          If JOBVSR = 'V', VSR will contain the right Schur vectors.
!>          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 DOUBLE PRECISION 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.
!>          If N = 0, LWORK >= 1, else LWORK >= MAX(8*N,6*N+16).
!>          For good performance, LWORK must generally be larger.
!>
!>          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]BWORK
!>          BWORK is LOGICAL array, dimension (N)
!>          Not referenced if SORT = '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 ALPHAR(j), ALPHAI(j), and BETA(j) should
!>                be correct for j=INFO+1,...,N.
!>          > N:  =N+1: other than QZ iteration failed in DHGEQZ.
!>                =N+2: after reordering, roundoff changed values of
!>                      some complex eigenvalues so that leading
!>                      eigenvalues in the Generalized Schur form no
!>                      longer satisfy SELCTG=.TRUE.  This could also
!>                      be caused due to scaling.
!>                =N+3: reordering failed in DTGSEN.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 279 of file dgges.f.

283*
284* -- LAPACK driver routine --
285* -- LAPACK is a software package provided by Univ. of Tennessee, --
286* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
287*
288* .. Scalar Arguments ..
289 CHARACTER JOBVSL, JOBVSR, SORT
290 INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
291* ..
292* .. Array Arguments ..
293 LOGICAL BWORK( * )
294 DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
295 $ B( LDB, * ), BETA( * ), VSL( LDVSL, * ),
296 $ VSR( LDVSR, * ), WORK( * )
297* ..
298* .. Function Arguments ..
299 LOGICAL SELCTG
300 EXTERNAL selctg
301* ..
302*
303* =====================================================================
304*
305* .. Parameters ..
306 DOUBLE PRECISION ZERO, ONE
307 parameter( zero = 0.0d+0, one = 1.0d+0 )
308* ..
309* .. Local Scalars ..
310 LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
311 $ LQUERY, LST2SL, WANTST
312 INTEGER I, ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT,
313 $ ILO, IP, IRIGHT, IROWS, ITAU, IWRK, MAXWRK,
314 $ MINWRK
315 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PVSL,
316 $ PVSR, SAFMAX, SAFMIN, SMLNUM
317* ..
318* .. Local Arrays ..
319 INTEGER IDUM( 1 )
320 DOUBLE PRECISION DIF( 2 )
321* ..
322* .. External Subroutines ..
323 EXTERNAL dgeqrf, dggbak, dggbal, dgghrd, dhgeqz,
324 $ dlacpy,
326* ..
327* .. External Functions ..
328 LOGICAL LSAME
329 INTEGER ILAENV
330 DOUBLE PRECISION DLAMCH, DLANGE
331 EXTERNAL lsame, ilaenv, dlamch, dlange
332* ..
333* .. Intrinsic Functions ..
334 INTRINSIC abs, max, sqrt
335* ..
336* .. Executable Statements ..
337*
338* Decode the input arguments
339*
340 IF( lsame( jobvsl, 'N' ) ) THEN
341 ijobvl = 1
342 ilvsl = .false.
343 ELSE IF( lsame( jobvsl, 'V' ) ) THEN
344 ijobvl = 2
345 ilvsl = .true.
346 ELSE
347 ijobvl = -1
348 ilvsl = .false.
349 END IF
350*
351 IF( lsame( jobvsr, 'N' ) ) THEN
352 ijobvr = 1
353 ilvsr = .false.
354 ELSE IF( lsame( jobvsr, 'V' ) ) THEN
355 ijobvr = 2
356 ilvsr = .true.
357 ELSE
358 ijobvr = -1
359 ilvsr = .false.
360 END IF
361*
362 wantst = lsame( sort, 'S' )
363*
364* Test the input arguments
365*
366 info = 0
367 lquery = ( lwork.EQ.-1 )
368 IF( ijobvl.LE.0 ) THEN
369 info = -1
370 ELSE IF( ijobvr.LE.0 ) THEN
371 info = -2
372 ELSE IF( ( .NOT.wantst ) .AND.
373 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
374 info = -3
375 ELSE IF( n.LT.0 ) THEN
376 info = -5
377 ELSE IF( lda.LT.max( 1, n ) ) THEN
378 info = -7
379 ELSE IF( ldb.LT.max( 1, n ) ) THEN
380 info = -9
381 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
382 info = -15
383 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
384 info = -17
385 END IF
386*
387* Compute workspace
388* (Note: Comments in the code beginning "Workspace:" describe the
389* minimal amount of workspace needed at that point in the code,
390* as well as the preferred amount for good performance.
391* NB refers to the optimal block size for the immediately
392* following subroutine, as returned by ILAENV.)
393*
394 IF( info.EQ.0 ) THEN
395 IF( n.GT.0 )THEN
396 minwrk = max( 8*n, 6*n + 16 )
397 maxwrk = minwrk - n +
398 $ n*ilaenv( 1, 'DGEQRF', ' ', n, 1, n, 0 )
399 maxwrk = max( maxwrk, minwrk - n +
400 $ n*ilaenv( 1, 'DORMQR', ' ', n, 1, n, -1 ) )
401 IF( ilvsl ) THEN
402 maxwrk = max( maxwrk, minwrk - n +
403 $ n*ilaenv( 1, 'DORGQR', ' ', n, 1, n,
404 $ -1 ) )
405 END IF
406 ELSE
407 minwrk = 1
408 maxwrk = 1
409 END IF
410 work( 1 ) = maxwrk
411*
412 IF( lwork.LT.minwrk .AND. .NOT.lquery )
413 $ info = -19
414 END IF
415*
416 IF( info.NE.0 ) THEN
417 CALL xerbla( 'DGGES ', -info )
418 RETURN
419 ELSE IF( lquery ) THEN
420 RETURN
421 END IF
422*
423* Quick return if possible
424*
425 IF( n.EQ.0 ) THEN
426 sdim = 0
427 RETURN
428 END IF
429*
430* Get machine constants
431*
432 eps = dlamch( 'P' )
433 safmin = dlamch( 'S' )
434 safmax = one / safmin
435 smlnum = sqrt( safmin ) / eps
436 bignum = one / smlnum
437*
438* Scale A if max element outside range [SMLNUM,BIGNUM]
439*
440 anrm = dlange( 'M', n, n, a, lda, work )
441 ilascl = .false.
442 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
443 anrmto = smlnum
444 ilascl = .true.
445 ELSE IF( anrm.GT.bignum ) THEN
446 anrmto = bignum
447 ilascl = .true.
448 END IF
449 IF( ilascl )
450 $ CALL dlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
451*
452* Scale B if max element outside range [SMLNUM,BIGNUM]
453*
454 bnrm = dlange( 'M', n, n, b, ldb, work )
455 ilbscl = .false.
456 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
457 bnrmto = smlnum
458 ilbscl = .true.
459 ELSE IF( bnrm.GT.bignum ) THEN
460 bnrmto = bignum
461 ilbscl = .true.
462 END IF
463 IF( ilbscl )
464 $ CALL dlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
465*
466* Permute the matrix to make it more nearly triangular
467* (Workspace: need 6*N + 2*N space for storing balancing factors)
468*
469 ileft = 1
470 iright = n + 1
471 iwrk = iright + n
472 CALL dggbal( 'P', n, a, lda, b, ldb, ilo, ihi, work( ileft ),
473 $ work( iright ), work( iwrk ), ierr )
474*
475* Reduce B to triangular form (QR decomposition of B)
476* (Workspace: need N, prefer N*NB)
477*
478 irows = ihi + 1 - ilo
479 icols = n + 1 - ilo
480 itau = iwrk
481 iwrk = itau + irows
482 CALL dgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
483 $ work( iwrk ), lwork+1-iwrk, ierr )
484*
485* Apply the orthogonal transformation to matrix A
486* (Workspace: need N, prefer N*NB)
487*
488 CALL dormqr( 'L', 'T', irows, icols, irows, b( ilo, ilo ), ldb,
489 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
490 $ lwork+1-iwrk, ierr )
491*
492* Initialize VSL
493* (Workspace: need N, prefer N*NB)
494*
495 IF( ilvsl ) THEN
496 CALL dlaset( 'Full', n, n, zero, one, vsl, ldvsl )
497 IF( irows.GT.1 ) THEN
498 CALL dlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
499 $ vsl( ilo+1, ilo ), ldvsl )
500 END IF
501 CALL dorgqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
502 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
503 END IF
504*
505* Initialize VSR
506*
507 IF( ilvsr )
508 $ CALL dlaset( 'Full', n, n, zero, one, vsr, ldvsr )
509*
510* Reduce to generalized Hessenberg form
511* (Workspace: none needed)
512*
513 CALL dgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
514 $ ldvsl, vsr, ldvsr, ierr )
515*
516* Perform QZ algorithm, computing Schur vectors if desired
517* (Workspace: need N)
518*
519 iwrk = itau
520 CALL dhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
521 $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
522 $ work( iwrk ), lwork+1-iwrk, ierr )
523 IF( ierr.NE.0 ) THEN
524 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
525 info = ierr
526 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
527 info = ierr - n
528 ELSE
529 info = n + 1
530 END IF
531 GO TO 50
532 END IF
533*
534* Sort eigenvalues ALPHA/BETA if desired
535* (Workspace: need 4*N+16 )
536*
537 sdim = 0
538 IF( wantst ) THEN
539*
540* Undo scaling on eigenvalues before SELCTGing
541*
542 IF( ilascl ) THEN
543 CALL dlascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n,
544 $ ierr )
545 CALL dlascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n,
546 $ ierr )
547 END IF
548 IF( ilbscl )
549 $ CALL dlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n,
550 $ ierr )
551*
552* Select eigenvalues
553*
554 DO 10 i = 1, n
555 bwork( i ) = selctg( alphar( i ), alphai( i ),
556 $ beta( i ) )
557 10 CONTINUE
558*
559 CALL dtgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
560 $ alphar,
561 $ alphai, beta, vsl, ldvsl, vsr, ldvsr, sdim, pvsl,
562 $ pvsr, dif, work( iwrk ), lwork-iwrk+1, idum, 1,
563 $ ierr )
564 IF( ierr.EQ.1 )
565 $ info = n + 3
566*
567 END IF
568*
569* Apply back-permutation to VSL and VSR
570* (Workspace: none needed)
571*
572 IF( ilvsl )
573 $ CALL dggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
574 $ work( iright ), n, vsl, ldvsl, ierr )
575*
576 IF( ilvsr )
577 $ CALL dggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
578 $ work( iright ), n, vsr, ldvsr, ierr )
579*
580* Check if unscaling would cause over/underflow, if so, rescale
581* (ALPHAR(I),ALPHAI(I),BETA(I)) so BETA(I) is on the order of
582* B(I,I) and ALPHAR(I) and ALPHAI(I) are on the order of A(I,I)
583*
584 IF( ilascl ) THEN
585 DO 20 i = 1, n
586 IF( alphai( i ).NE.zero ) THEN
587 IF( ( alphar( i ) / safmax ).GT.( anrmto / anrm ) .OR.
588 $ ( safmin / alphar( i ) ).GT.( anrm / anrmto ) ) THEN
589 work( 1 ) = abs( a( i, i ) / alphar( i ) )
590 beta( i ) = beta( i )*work( 1 )
591 alphar( i ) = alphar( i )*work( 1 )
592 alphai( i ) = alphai( i )*work( 1 )
593 ELSE IF( ( alphai( i ) / safmax ).GT.
594 $ ( anrmto / anrm ) .OR.
595 $ ( safmin / alphai( i ) ).GT.( anrm / anrmto ) )
596 $ THEN
597 work( 1 ) = abs( a( i, i+1 ) / alphai( i ) )
598 beta( i ) = beta( i )*work( 1 )
599 alphar( i ) = alphar( i )*work( 1 )
600 alphai( i ) = alphai( i )*work( 1 )
601 END IF
602 END IF
603 20 CONTINUE
604 END IF
605*
606 IF( ilbscl ) THEN
607 DO 30 i = 1, n
608 IF( alphai( i ).NE.zero ) THEN
609 IF( ( beta( i ) / safmax ).GT.( bnrmto / bnrm ) .OR.
610 $ ( safmin / beta( i ) ).GT.( bnrm / bnrmto ) ) THEN
611 work( 1 ) = abs( b( i, i ) / beta( i ) )
612 beta( i ) = beta( i )*work( 1 )
613 alphar( i ) = alphar( i )*work( 1 )
614 alphai( i ) = alphai( i )*work( 1 )
615 END IF
616 END IF
617 30 CONTINUE
618 END IF
619*
620* Undo scaling
621*
622 IF( ilascl ) THEN
623 CALL dlascl( 'H', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
624 CALL dlascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n,
625 $ ierr )
626 CALL dlascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n,
627 $ ierr )
628 END IF
629*
630 IF( ilbscl ) THEN
631 CALL dlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
632 CALL dlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
633 END IF
634*
635 IF( wantst ) THEN
636*
637* Check if reordering is correct
638*
639 lastsl = .true.
640 lst2sl = .true.
641 sdim = 0
642 ip = 0
643 DO 40 i = 1, n
644 cursl = selctg( alphar( i ), alphai( i ), beta( i ) )
645 IF( alphai( i ).EQ.zero ) THEN
646 IF( cursl )
647 $ sdim = sdim + 1
648 ip = 0
649 IF( cursl .AND. .NOT.lastsl )
650 $ info = n + 2
651 ELSE
652 IF( ip.EQ.1 ) THEN
653*
654* Last eigenvalue of conjugate pair
655*
656 cursl = cursl .OR. lastsl
657 lastsl = cursl
658 IF( cursl )
659 $ sdim = sdim + 2
660 ip = -1
661 IF( cursl .AND. .NOT.lst2sl )
662 $ info = n + 2
663 ELSE
664*
665* First eigenvalue of conjugate pair
666*
667 ip = 1
668 END IF
669 END IF
670 lst2sl = lastsl
671 lastsl = cursl
672 40 CONTINUE
673*
674 END IF
675*
676 50 CONTINUE
677*
678 work( 1 ) = maxwrk
679*
680 RETURN
681*
682* End of DGGES
683*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgeqrf(m, n, a, lda, tau, work, lwork, info)
DGEQRF
Definition dgeqrf.f:144
subroutine dggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
DGGBAK
Definition dggbak.f:146
subroutine dggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
DGGBAL
Definition dggbal.f:175
subroutine dgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
DGGHRD
Definition dgghrd.f:206
subroutine dhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alphar, alphai, beta, q, ldq, z, ldz, work, lwork, info)
DHGEQZ
Definition dhgeqz.f:303
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:112
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 dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:108
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dtgsen(ijob, wantq, wantz, select, n, a, lda, b, ldb, alphar, alphai, beta, q, ldq, z, ldz, m, pl, pr, dif, work, lwork, iwork, liwork, info)
DTGSEN
Definition dtgsen.f:450
subroutine dorgqr(m, n, k, a, lda, tau, work, lwork, info)
DORGQR
Definition dorgqr.f:126
subroutine dormqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
DORMQR
Definition dormqr.f:165
Here is the call graph for this function:
Here is the caller graph for this function: