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

◆ sgges()

subroutine sgges ( 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 
)

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

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

Purpose:
 SGGES 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
 "standardized" 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 >= 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 SHGEQZ.
                =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 281 of file sgges.f.

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