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

◆ sgeesx()

subroutine sgeesx ( character jobvs,
character sort,
external select,
character sense,
integer n,
real, dimension( lda, * ) a,
integer lda,
integer sdim,
real, dimension( * ) wr,
real, dimension( * ) wi,
real, dimension( ldvs, * ) vs,
integer ldvs,
real rconde,
real rcondv,
real, dimension( * ) work,
integer lwork,
integer, dimension( * ) iwork,
integer liwork,
logical, dimension( * ) bwork,
integer info )

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

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

Purpose:
!>
!> SGEESX computes for an N-by-N real nonsymmetric matrix A, the
!> eigenvalues, the real Schur form T, and, optionally, the matrix of
!> Schur vectors Z.  This gives the Schur factorization A = Z*T*(Z**T).
!>
!> Optionally, it also orders the eigenvalues on the diagonal of the
!> real Schur form so that selected eigenvalues are at the top left;
!> computes a reciprocal condition number for the average of the
!> selected eigenvalues (RCONDE); and computes a reciprocal condition
!> number for the right invariant subspace corresponding to the
!> selected eigenvalues (RCONDV).  The leading columns of Z form an
!> orthonormal basis for this invariant subspace.
!>
!> For further explanation of the reciprocal condition numbers RCONDE
!> and RCONDV, see Section 4.10 of the LAPACK Users' Guide (where
!> these quantities are called s and sep respectively).
!>
!> A real matrix is in real Schur form if it is upper quasi-triangular
!> with 1-by-1 and 2-by-2 blocks. 2-by-2 blocks will be standardized in
!> the form
!>           [  a  b  ]
!>           [  c  a  ]
!>
!> where b*c < 0. The eigenvalues of such a block are a +- sqrt(bc).
!> 
Parameters
[in]JOBVS
!>          JOBVS is CHARACTER*1
!>          = 'N': Schur vectors are not computed;
!>          = 'V': Schur vectors are computed.
!> 
[in]SORT
!>          SORT is CHARACTER*1
!>          Specifies whether or not to order the eigenvalues on the
!>          diagonal of the Schur form.
!>          = 'N': Eigenvalues are not ordered;
!>          = 'S': Eigenvalues are ordered (see SELECT).
!> 
[in]SELECT
!>          SELECT is a LOGICAL FUNCTION of two REAL arguments
!>          SELECT must be declared EXTERNAL in the calling subroutine.
!>          If SORT = 'S', SELECT is used to select eigenvalues to sort
!>          to the top left of the Schur form.
!>          If SORT = 'N', SELECT is not referenced.
!>          An eigenvalue WR(j)+sqrt(-1)*WI(j) is selected if
!>          SELECT(WR(j),WI(j)) is true; i.e., if either one of a
!>          complex conjugate pair of eigenvalues is selected, then both
!>          are.  Note that a selected complex eigenvalue may no longer
!>          satisfy SELECT(WR(j),WI(j)) = .TRUE. after ordering, since
!>          ordering may change the value of complex eigenvalues
!>          (especially if the eigenvalue is ill-conditioned); in this
!>          case INFO may be set to N+3 (see INFO below).
!> 
[in]SENSE
!>          SENSE is CHARACTER*1
!>          Determines which reciprocal condition numbers are computed.
!>          = 'N': None are computed;
!>          = 'E': Computed for average of selected eigenvalues only;
!>          = 'V': Computed for selected right invariant subspace only;
!>          = 'B': Computed for both.
!>          If SENSE = 'E', 'V' or 'B', SORT must equal 'S'.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A. N >= 0.
!> 
[in,out]A
!>          A is REAL array, dimension (LDA, N)
!>          On entry, the N-by-N matrix A.
!>          On exit, A is overwritten by its real Schur form T.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]SDIM
!>          SDIM is INTEGER
!>          If SORT = 'N', SDIM = 0.
!>          If SORT = 'S', SDIM = number of eigenvalues (after sorting)
!>                         for which SELECT is true. (Complex conjugate
!>                         pairs for which SELECT is true for either
!>                         eigenvalue count as 2.)
!> 
[out]WR
!>          WR is REAL array, dimension (N)
!> 
[out]WI
!>          WI is REAL array, dimension (N)
!>          WR and WI contain the real and imaginary parts, respectively,
!>          of the computed eigenvalues, in the same order that they
!>          appear on the diagonal of the output Schur form T.  Complex
!>          conjugate pairs of eigenvalues appear consecutively with the
!>          eigenvalue having the positive imaginary part first.
!> 
[out]VS
!>          VS is REAL array, dimension (LDVS,N)
!>          If JOBVS = 'V', VS contains the orthogonal matrix Z of Schur
!>          vectors.
!>          If JOBVS = 'N', VS is not referenced.
!> 
[in]LDVS
!>          LDVS is INTEGER
!>          The leading dimension of the array VS.  LDVS >= 1, and if
!>          JOBVS = 'V', LDVS >= N.
!> 
[out]RCONDE
!>          RCONDE is REAL
!>          If SENSE = 'E' or 'B', RCONDE contains the reciprocal
!>          condition number for the average of the selected eigenvalues.
!>          Not referenced if SENSE = 'N' or 'V'.
!> 
[out]RCONDV
!>          RCONDV is REAL
!>          If SENSE = 'V' or 'B', RCONDV contains the reciprocal
!>          condition number for the selected right invariant subspace.
!>          Not referenced if SENSE = 'N' or 'E'.
!> 
[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.  LWORK >= max(1,3*N).
!>          Also, if SENSE = 'E' or 'V' or 'B',
!>          LWORK >= N+2*SDIM*(N-SDIM), where SDIM is the number of
!>          selected eigenvalues computed by this routine.  Note that
!>          N+2*SDIM*(N-SDIM) <= N+N*N/2. Note also that an error is only
!>          returned if LWORK < max(1,3*N), but if SENSE = 'E' or 'V' or
!>          'B' this may not be large enough.
!>          For good performance, LWORK must generally be larger.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates upper bounds on the optimal sizes of the
!>          arrays WORK and IWORK, returns these values as the first
!>          entries of the WORK and IWORK arrays, and no error messages
!>          related to LWORK or LIWORK are issued by XERBLA.
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
!>          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
!> 
[in]LIWORK
!>          LIWORK is INTEGER
!>          The dimension of the array IWORK.
!>          LIWORK >= 1; if SENSE = 'V' or 'B', LIWORK >= SDIM*(N-SDIM).
!>          Note that SDIM*(N-SDIM) <= N*N/4. Note also that an error is
!>          only returned if LIWORK < 1, but if SENSE = 'V' or 'B' this
!>          may not be large enough.
!>
!>          If LIWORK = -1, then a workspace query is assumed; the
!>          routine only calculates upper bounds on the optimal sizes of
!>          the arrays WORK and IWORK, returns these values as the first
!>          entries of the WORK and IWORK arrays, and no error messages
!>          related to LWORK or LIWORK are 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.
!>          > 0: if INFO = i, and i is
!>             <= N: the QR algorithm failed to compute all the
!>                   eigenvalues; elements 1:ILO-1 and i+1:N of WR and WI
!>                   contain those eigenvalues which have converged; if
!>                   JOBVS = 'V', VS contains the transformation which
!>                   reduces A to its partially converged Schur form.
!>             = N+1: the eigenvalues could not be reordered because some
!>                   eigenvalues were too close to separate (the problem
!>                   is very ill-conditioned);
!>             = N+2: after reordering, roundoff changed values of some
!>                   complex eigenvalues so that leading eigenvalues in
!>                   the Schur form no longer satisfy SELECT=.TRUE.  This
!>                   could also be caused by underflow due to scaling.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 276 of file sgeesx.f.

279*
280* -- LAPACK driver routine --
281* -- LAPACK is a software package provided by Univ. of Tennessee, --
282* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
283*
284* .. Scalar Arguments ..
285 CHARACTER JOBVS, SENSE, SORT
286 INTEGER INFO, LDA, LDVS, LIWORK, LWORK, N, SDIM
287 REAL RCONDE, RCONDV
288* ..
289* .. Array Arguments ..
290 LOGICAL BWORK( * )
291 INTEGER IWORK( * )
292 REAL A( LDA, * ), VS( LDVS, * ), WI( * ), WORK( * ),
293 $ WR( * )
294* ..
295* .. Function Arguments ..
296 LOGICAL SELECT
297 EXTERNAL SELECT
298* ..
299*
300* =====================================================================
301*
302* .. Parameters ..
303 REAL ZERO, ONE
304 parameter( zero = 0.0e0, one = 1.0e0 )
305* ..
306* .. Local Scalars ..
307 LOGICAL CURSL, LASTSL, LQUERY, LST2SL, SCALEA, WANTSB,
308 $ WANTSE, WANTSN, WANTST, WANTSV, WANTVS
309 INTEGER HSWORK, I, I1, I2, IBAL, ICOND, IERR, IEVAL,
310 $ IHI, ILO, INXT, IP, ITAU, IWRK, LWRK, LIWRK,
311 $ MAXWRK, MINWRK
312 REAL ANRM, BIGNUM, CSCALE, EPS, SMLNUM
313* ..
314* .. Local Arrays ..
315 REAL DUM( 1 )
316* ..
317* .. External Subroutines ..
318 EXTERNAL scopy, sgebak, sgebal, sgehrd,
319 $ shseqr,
321* ..
322* .. External Functions ..
323 LOGICAL LSAME
324 INTEGER ILAENV
325 REAL SLAMCH, SLANGE, SROUNDUP_LWORK
326 EXTERNAL lsame, ilaenv, slamch,
328* ..
329* .. Intrinsic Functions ..
330 INTRINSIC max, sqrt
331* ..
332* .. Executable Statements ..
333*
334* Test the input arguments
335*
336 info = 0
337 wantvs = lsame( jobvs, 'V' )
338 wantst = lsame( sort, 'S' )
339 wantsn = lsame( sense, 'N' )
340 wantse = lsame( sense, 'E' )
341 wantsv = lsame( sense, 'V' )
342 wantsb = lsame( sense, 'B' )
343 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
344*
345 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
346 info = -1
347 ELSE IF( ( .NOT.wantst ) .AND.
348 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
349 info = -2
350 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
351 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
352 info = -4
353 ELSE IF( n.LT.0 ) THEN
354 info = -5
355 ELSE IF( lda.LT.max( 1, n ) ) THEN
356 info = -7
357 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
358 info = -12
359 END IF
360*
361* Compute workspace
362* (Note: Comments in the code beginning "RWorkspace:" describe the
363* minimal amount of real workspace needed at that point in the
364* code, as well as the preferred amount for good performance.
365* IWorkspace refers to integer workspace.
366* NB refers to the optimal block size for the immediately
367* following subroutine, as returned by ILAENV.
368* HSWORK refers to the workspace preferred by SHSEQR, as
369* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
370* the worst case.
371* If SENSE = 'E', 'V' or 'B', then the amount of workspace needed
372* depends on SDIM, which is computed by the routine STRSEN later
373* in the code.)
374*
375 IF( info.EQ.0 ) THEN
376 liwrk = 1
377 IF( n.EQ.0 ) THEN
378 minwrk = 1
379 lwrk = 1
380 ELSE
381 maxwrk = 2*n + n*ilaenv( 1, 'SGEHRD', ' ', n, 1, n, 0 )
382 minwrk = 3*n
383*
384 CALL shseqr( 'S', jobvs, n, 1, n, a, lda, wr, wi, vs,
385 $ ldvs,
386 $ work, -1, ieval )
387 hswork = int( work( 1 ) )
388*
389 IF( .NOT.wantvs ) THEN
390 maxwrk = max( maxwrk, n + hswork )
391 ELSE
392 maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
393 $ 'SORGHR', ' ', n, 1, n, -1 ) )
394 maxwrk = max( maxwrk, n + hswork )
395 END IF
396 lwrk = maxwrk
397 IF( .NOT.wantsn )
398 $ lwrk = max( lwrk, n + ( n*n )/2 )
399 IF( wantsv .OR. wantsb )
400 $ liwrk = ( n*n )/4
401 END IF
402 iwork( 1 ) = liwrk
403 work( 1 ) = sroundup_lwork(lwrk)
404*
405 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
406 info = -16
407 ELSE IF( liwork.LT.1 .AND. .NOT.lquery ) THEN
408 info = -18
409 END IF
410 END IF
411*
412 IF( info.NE.0 ) THEN
413 CALL xerbla( 'SGEESX', -info )
414 RETURN
415 ELSE IF( lquery ) THEN
416 RETURN
417 END IF
418*
419* Quick return if possible
420*
421 IF( n.EQ.0 ) THEN
422 sdim = 0
423 RETURN
424 END IF
425*
426* Get machine constants
427*
428 eps = slamch( 'P' )
429 smlnum = slamch( 'S' )
430 bignum = one / smlnum
431 smlnum = sqrt( smlnum ) / eps
432 bignum = one / smlnum
433*
434* Scale A if max element outside range [SMLNUM,BIGNUM]
435*
436 anrm = slange( 'M', n, n, a, lda, dum )
437 scalea = .false.
438 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
439 scalea = .true.
440 cscale = smlnum
441 ELSE IF( anrm.GT.bignum ) THEN
442 scalea = .true.
443 cscale = bignum
444 END IF
445 IF( scalea )
446 $ CALL slascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
447*
448* Permute the matrix to make it more nearly triangular
449* (RWorkspace: need N)
450*
451 ibal = 1
452 CALL sgebal( 'P', n, a, lda, ilo, ihi, work( ibal ), ierr )
453*
454* Reduce to upper Hessenberg form
455* (RWorkspace: need 3*N, prefer 2*N+N*NB)
456*
457 itau = n + ibal
458 iwrk = n + itau
459 CALL sgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
460 $ lwork-iwrk+1, ierr )
461*
462 IF( wantvs ) THEN
463*
464* Copy Householder vectors to VS
465*
466 CALL slacpy( 'L', n, n, a, lda, vs, ldvs )
467*
468* Generate orthogonal matrix in VS
469* (RWorkspace: need 3*N-1, prefer 2*N+(N-1)*NB)
470*
471 CALL sorghr( n, ilo, ihi, vs, ldvs, work( itau ),
472 $ work( iwrk ),
473 $ lwork-iwrk+1, ierr )
474 END IF
475*
476 sdim = 0
477*
478* Perform QR iteration, accumulating Schur vectors in VS if desired
479* (RWorkspace: need N+1, prefer N+HSWORK (see comments) )
480*
481 iwrk = itau
482 CALL shseqr( 'S', jobvs, n, ilo, ihi, a, lda, wr, wi, vs, ldvs,
483 $ work( iwrk ), lwork-iwrk+1, ieval )
484 IF( ieval.GT.0 )
485 $ info = ieval
486*
487* Sort eigenvalues if desired
488*
489 IF( wantst .AND. info.EQ.0 ) THEN
490 IF( scalea ) THEN
491 CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, wr, n, ierr )
492 CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, wi, n, ierr )
493 END IF
494 DO 10 i = 1, n
495 bwork( i ) = SELECT( wr( i ), wi( i ) )
496 10 CONTINUE
497*
498* Reorder eigenvalues, transform Schur vectors, and compute
499* reciprocal condition numbers
500* (RWorkspace: if SENSE is not 'N', need N+2*SDIM*(N-SDIM)
501* otherwise, need N )
502* (IWorkspace: if SENSE is 'V' or 'B', need SDIM*(N-SDIM)
503* otherwise, need 0 )
504*
505 CALL strsen( sense, jobvs, bwork, n, a, lda, vs, ldvs, wr,
506 $ wi,
507 $ sdim, rconde, rcondv, work( iwrk ), lwork-iwrk+1,
508 $ iwork, liwork, icond )
509 IF( .NOT.wantsn )
510 $ maxwrk = max( maxwrk, n+2*sdim*( n-sdim ) )
511 IF( icond.EQ.-15 ) THEN
512*
513* Not enough real workspace
514*
515 info = -16
516 ELSE IF( icond.EQ.-17 ) THEN
517*
518* Not enough integer workspace
519*
520 info = -18
521 ELSE IF( icond.GT.0 ) THEN
522*
523* STRSEN failed to reorder or to restore standard Schur form
524*
525 info = icond + n
526 END IF
527 END IF
528*
529 IF( wantvs ) THEN
530*
531* Undo balancing
532* (RWorkspace: need N)
533*
534 CALL sgebak( 'P', 'R', n, ilo, ihi, work( ibal ), n, vs,
535 $ ldvs,
536 $ ierr )
537 END IF
538*
539 IF( scalea ) THEN
540*
541* Undo scaling for the Schur form of A
542*
543 CALL slascl( 'H', 0, 0, cscale, anrm, n, n, a, lda, ierr )
544 CALL scopy( n, a, lda+1, wr, 1 )
545 IF( ( wantsv .OR. wantsb ) .AND. info.EQ.0 ) THEN
546 dum( 1 ) = rcondv
547 CALL slascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1,
548 $ ierr )
549 rcondv = dum( 1 )
550 END IF
551 IF( cscale.EQ.smlnum ) THEN
552*
553* If scaling back towards underflow, adjust WI if an
554* offdiagonal element of a 2-by-2 block in the Schur form
555* underflows.
556*
557 IF( ieval.GT.0 ) THEN
558 i1 = ieval + 1
559 i2 = ihi - 1
560 CALL slascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
561 $ ierr )
562 ELSE IF( wantst ) THEN
563 i1 = 1
564 i2 = n - 1
565 ELSE
566 i1 = ilo
567 i2 = ihi - 1
568 END IF
569 inxt = i1 - 1
570 DO 20 i = i1, i2
571 IF( i.LT.inxt )
572 $ GO TO 20
573 IF( wi( i ).EQ.zero ) THEN
574 inxt = i + 1
575 ELSE
576 IF( a( i+1, i ).EQ.zero ) THEN
577 wi( i ) = zero
578 wi( i+1 ) = zero
579 ELSE IF( a( i+1, i ).NE.zero .AND. a( i, i+1 ).EQ.
580 $ zero ) THEN
581 wi( i ) = zero
582 wi( i+1 ) = zero
583 IF( i.GT.1 )
584 $ CALL sswap( i-1, a( 1, i ), 1, a( 1, i+1 ),
585 $ 1 )
586 IF( n.GT.i+1 )
587 $ CALL sswap( n-i-1, a( i, i+2 ), lda,
588 $ a( i+1, i+2 ), lda )
589 IF( wantvs ) THEN
590 CALL sswap( n, vs( 1, i ), 1, vs( 1, i+1 ),
591 $ 1 )
592 END IF
593 a( i, i+1 ) = a( i+1, i )
594 a( i+1, i ) = zero
595 END IF
596 inxt = i + 2
597 END IF
598 20 CONTINUE
599 END IF
600 CALL slascl( 'G', 0, 0, cscale, anrm, n-ieval, 1,
601 $ wi( ieval+1 ), max( n-ieval, 1 ), ierr )
602 END IF
603*
604 IF( wantst .AND. info.EQ.0 ) THEN
605*
606* Check if reordering successful
607*
608 lastsl = .true.
609 lst2sl = .true.
610 sdim = 0
611 ip = 0
612 DO 30 i = 1, n
613 cursl = SELECT( wr( i ), wi( i ) )
614 IF( wi( i ).EQ.zero ) THEN
615 IF( cursl )
616 $ sdim = sdim + 1
617 ip = 0
618 IF( cursl .AND. .NOT.lastsl )
619 $ info = n + 2
620 ELSE
621 IF( ip.EQ.1 ) THEN
622*
623* Last eigenvalue of conjugate pair
624*
625 cursl = cursl .OR. lastsl
626 lastsl = cursl
627 IF( cursl )
628 $ sdim = sdim + 2
629 ip = -1
630 IF( cursl .AND. .NOT.lst2sl )
631 $ info = n + 2
632 ELSE
633*
634* First eigenvalue of conjugate pair
635*
636 ip = 1
637 END IF
638 END IF
639 lst2sl = lastsl
640 lastsl = cursl
641 30 CONTINUE
642 END IF
643*
644 work( 1 ) = sroundup_lwork(maxwrk)
645 IF( wantsv .OR. wantsb ) THEN
646 iwork( 1 ) = sdim*(n-sdim)
647 ELSE
648 iwork( 1 ) = 1
649 END IF
650*
651 RETURN
652*
653* End of SGEESX
654*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
SGEBAK
Definition sgebak.f:128
subroutine sgebal(job, n, a, lda, ilo, ihi, scale, info)
SGEBAL
Definition sgebal.f:161
subroutine sgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
SGEHRD
Definition sgehrd.f:166
subroutine shseqr(job, compz, n, ilo, ihi, h, ldh, wr, wi, z, ldz, work, lwork, info)
SHSEQR
Definition shseqr.f:314
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
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
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine strsen(job, compq, select, n, t, ldt, q, ldq, wr, wi, m, s, sep, work, lwork, iwork, liwork, info)
STRSEN
Definition strsen.f:313
subroutine sorghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
SORGHR
Definition sorghr.f:125
Here is the call graph for this function:
Here is the caller graph for this function: