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

◆ sgges3()

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

SGGES3 computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices (blocked algorithm)

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

Purpose:
!>
!> SGGES3 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
!> SGGEV 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 REAL 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 REAL 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 REAL 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 REAL array, dimension (N)
!> 
[out]ALPHAI
!>          ALPHAI is REAL array, dimension (N)
!> 
[out]BETA
!>          BETA is REAL 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 REAL 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 REAL 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 REAL 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 >= 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 SLAQZ0.
!>                =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 STGSEN.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 279 of file sgges3.f.

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