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

◆ cgegv()

subroutine cgegv ( character jobvl,
character jobvr,
integer n,
complex, dimension( lda, * ) a,
integer lda,
complex, dimension( ldb, * ) b,
integer ldb,
complex, dimension( * ) alpha,
complex, dimension( * ) beta,
complex, dimension( ldvl, * ) vl,
integer ldvl,
complex, dimension( ldvr, * ) vr,
integer ldvr,
complex, dimension( * ) work,
integer lwork,
real, dimension( * ) rwork,
integer info )

CGEGV computes the eigenvalues and, optionally, the left and/or right eigenvectors of a complex matrix pair (A,B).

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

Purpose:
!>
!> This routine is deprecated and has been replaced by routine CGGEV.
!>
!> CGEGV computes the eigenvalues and, optionally, the left and/or right
!> eigenvectors of a complex matrix pair (A,B).
!> Given two square matrices A and B,
!> the generalized nonsymmetric eigenvalue problem (GNEP) is to find the
!> eigenvalues lambda and corresponding (non-zero) eigenvectors x such
!> that
!>    A*x = lambda*B*x.
!>
!> An alternate form is to find the eigenvalues mu and corresponding
!> eigenvectors y such that
!>    mu*A*y = B*y.
!>
!> These two forms are equivalent with mu = 1/lambda and x = y if
!> neither lambda nor mu is zero.  In order to deal with the case that
!> lambda or mu is zero or small, two values alpha and beta are returned
!> for each eigenvalue, such that lambda = alpha/beta and
!> mu = beta/alpha.
!>
!> The vectors x and y in the above equations are right eigenvectors of
!> the matrix pair (A,B).  Vectors u and v satisfying
!>    u**H*A = lambda*u**H*B  or  mu*v**H*A = v**H*B
!> are left eigenvectors of (A,B).
!>
!> Note: this routine performs  on A and B
!> 
Parameters
[in]JOBVL
!>          JOBVL is CHARACTER*1
!>          = 'N':  do not compute the left generalized eigenvectors;
!>          = 'V':  compute the left generalized eigenvectors (returned
!>                  in VL).
!> 
[in]JOBVR
!>          JOBVR is CHARACTER*1
!>          = 'N':  do not compute the right generalized eigenvectors;
!>          = 'V':  compute the right generalized eigenvectors (returned
!>                  in VR).
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrices A, B, VL, and VR.  N >= 0.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA, N)
!>          On entry, the matrix A.
!>          If JOBVL = 'V' or JOBVR = 'V', then on exit A
!>          contains the Schur form of A from the generalized Schur
!>          factorization of the pair (A,B) after balancing.  If no
!>          eigenvectors were computed, then only the diagonal elements
!>          of the Schur form will be correct.  See CGGHRD and CHGEQZ
!>          for details.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of A.  LDA >= max(1,N).
!> 
[in,out]B
!>          B is COMPLEX array, dimension (LDB, N)
!>          On entry, the matrix B.
!>          If JOBVL = 'V' or JOBVR = 'V', then on exit B contains the
!>          upper triangular matrix obtained from B in the generalized
!>          Schur factorization of the pair (A,B) after balancing.
!>          If no eigenvectors were computed, then only the diagonal
!>          elements of B will be correct.  See CGGHRD and CHGEQZ for
!>          details.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of B.  LDB >= max(1,N).
!> 
[out]ALPHA
!>          ALPHA is COMPLEX array, dimension (N)
!>          The complex scalars alpha that define the eigenvalues of
!>          GNEP.
!> 
[out]BETA
!>          BETA is COMPLEX array, dimension (N)
!>          The complex scalars beta that define the eigenvalues of GNEP.
!>
!>          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]VL
!>          VL is COMPLEX array, dimension (LDVL,N)
!>          If JOBVL = 'V', the left eigenvectors u(j) are stored
!>          in the columns of VL, in the same order as their eigenvalues.
!>          Each eigenvector is scaled so that its largest component has
!>          abs(real part) + abs(imag. part) = 1, except for eigenvectors
!>          corresponding to an eigenvalue with alpha = beta = 0, which
!>          are set to zero.
!>          Not referenced if JOBVL = 'N'.
!> 
[in]LDVL
!>          LDVL is INTEGER
!>          The leading dimension of the matrix VL. LDVL >= 1, and
!>          if JOBVL = 'V', LDVL >= N.
!> 
[out]VR
!>          VR is COMPLEX array, dimension (LDVR,N)
!>          If JOBVR = 'V', the right eigenvectors x(j) are stored
!>          in the columns of VR, in the same order as their eigenvalues.
!>          Each eigenvector is scaled so that its largest component has
!>          abs(real part) + abs(imag. part) = 1, except for eigenvectors
!>          corresponding to an eigenvalue with alpha = beta = 0, which
!>          are set to zero.
!>          Not referenced if JOBVR = 'N'.
!> 
[in]LDVR
!>          LDVR is INTEGER
!>          The leading dimension of the matrix VR. LDVR >= 1, and
!>          if JOBVR = 'V', LDVR >= N.
!> 
[out]WORK
!>          WORK is COMPLEX 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 CGEQRF, CUNMQR, and CUNGQR.)  Then compute:
!>          NB  -- MAX of the blocksizes for CGEQRF, CUNMQR, and CUNGQR;
!>          The optimal LWORK is  MAX( 2*N, 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 REAL array, dimension (8*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.  No eigenvectors have been
!>                calculated, 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 CGGBAL
!>                =N+2: error return from CGEQRF
!>                =N+3: error return from CUNMQR
!>                =N+4: error return from CUNGQR
!>                =N+5: error return from CGGHRD
!>                =N+6: error return from CHGEQZ (other than failed
!>                                               iteration)
!>                =N+7: error return from CTGEVC
!>                =N+8: error return from CGGBAK (computing VL)
!>                =N+9: error return from CGGBAK (computing VR)
!>                =N+10: error return from CLASCL (various calls)
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Balancing
!>  ---------
!>
!>  This driver calls CGGBAL to both permute and scale rows and columns
!>  of A and B.  The permutations PL and PR are chosen so that PL*A*PR
!>  and PL*B*R will be upper triangular except for the diagonal blocks
!>  A(i:j,i:j) and B(i:j,i:j), with i and j as close together as
!>  possible.  The diagonal scaling matrices DL and DR are chosen so
!>  that the pair  DL*PL*A*PR*DR, DL*PL*B*PR*DR have elements close to
!>  one (except for the elements that start out zero.)
!>
!>  After the eigenvalues and eigenvectors of the balanced matrices
!>  have been computed, CGGBAK transforms the eigenvectors back to what
!>  they would have been (in perfect arithmetic) if they had not been
!>  balanced.
!>
!>  Contents of A and B on Exit
!>  -------- -- - --- - -- ----
!>
!>  If any eigenvectors are computed (either JOBVL='V' or JOBVR='V' or
!>  both), then on exit the arrays A and B will contain the complex Schur
!>  form[*] of the  versions of A and B.  If no eigenvectors
!>  are computed, then only the diagonal blocks will be correct.
!>
!>  [*] In other words, upper triangular form.
!> 

Definition at line 278 of file cgegv.f.

280*
281* -- LAPACK driver routine --
282* -- LAPACK is a software package provided by Univ. of Tennessee, --
283* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
284*
285* .. Scalar Arguments ..
286 CHARACTER JOBVL, JOBVR
287 INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, N
288* ..
289* .. Array Arguments ..
290 REAL RWORK( * )
291 COMPLEX A( LDA, * ), ALPHA( * ), B( LDB, * ),
292 $ BETA( * ), VL( LDVL, * ), VR( LDVR, * ),
293 $ WORK( * )
294* ..
295*
296* =====================================================================
297*
298* .. Parameters ..
299 REAL ZERO, ONE
300 parameter( zero = 0.0e0, one = 1.0e0 )
301 COMPLEX CZERO, CONE
302 parameter( czero = ( 0.0e0, 0.0e0 ),
303 $ cone = ( 1.0e0, 0.0e0 ) )
304* ..
305* .. Local Scalars ..
306 LOGICAL ILIMIT, ILV, ILVL, ILVR, LQUERY
307 CHARACTER CHTEMP
308 INTEGER ICOLS, IHI, IINFO, IJOBVL, IJOBVR, ILEFT, ILO,
309 $ IN, IRIGHT, IROWS, IRWORK, ITAU, IWORK, JC, JR,
310 $ LOPT, LWKMIN, LWKOPT, NB, NB1, NB2, NB3
311 REAL ABSAI, ABSAR, ABSB, ANRM, ANRM1, ANRM2, BNRM,
312 $ BNRM1, BNRM2, EPS, SAFMAX, SAFMIN, SALFAI,
313 $ SALFAR, SBETA, SCALE, TEMP
314 COMPLEX X
315* ..
316* .. Local Arrays ..
317 LOGICAL LDUMMA( 1 )
318* ..
319* .. External Subroutines ..
320 EXTERNAL cgeqrf, cggbak, cggbal, cgghrd, chgeqz, clacpy,
322* ..
323* .. External Functions ..
324 LOGICAL LSAME
325 INTEGER ILAENV
326 REAL CLANGE, SLAMCH
327 EXTERNAL ilaenv, lsame, clange, slamch
328* ..
329* .. Intrinsic Functions ..
330 INTRINSIC abs, aimag, cmplx, int, max, real
331* ..
332* .. Statement Functions ..
333 REAL ABS1
334* ..
335* .. Statement Function definitions ..
336 abs1( x ) = abs( real( x ) ) + abs( aimag( x ) )
337* ..
338* .. Executable Statements ..
339*
340* Decode the input arguments
341*
342 IF( lsame( jobvl, 'N' ) ) THEN
343 ijobvl = 1
344 ilvl = .false.
345 ELSE IF( lsame( jobvl, 'V' ) ) THEN
346 ijobvl = 2
347 ilvl = .true.
348 ELSE
349 ijobvl = -1
350 ilvl = .false.
351 END IF
352*
353 IF( lsame( jobvr, 'N' ) ) THEN
354 ijobvr = 1
355 ilvr = .false.
356 ELSE IF( lsame( jobvr, 'V' ) ) THEN
357 ijobvr = 2
358 ilvr = .true.
359 ELSE
360 ijobvr = -1
361 ilvr = .false.
362 END IF
363 ilv = ilvl .OR. ilvr
364*
365* Test the input arguments
366*
367 lwkmin = max( 2*n, 1 )
368 lwkopt = lwkmin
369 work( 1 ) = lwkopt
370 lquery = ( lwork.EQ.-1 )
371 info = 0
372 IF( ijobvl.LE.0 ) THEN
373 info = -1
374 ELSE IF( ijobvr.LE.0 ) THEN
375 info = -2
376 ELSE IF( n.LT.0 ) THEN
377 info = -3
378 ELSE IF( lda.LT.max( 1, n ) ) THEN
379 info = -5
380 ELSE IF( ldb.LT.max( 1, n ) ) THEN
381 info = -7
382 ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) ) THEN
383 info = -11
384 ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) ) THEN
385 info = -13
386 ELSE IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
387 info = -15
388 END IF
389*
390 IF( info.EQ.0 ) THEN
391 nb1 = ilaenv( 1, 'CGEQRF', ' ', n, n, -1, -1 )
392 nb2 = ilaenv( 1, 'CUNMQR', ' ', n, n, n, -1 )
393 nb3 = ilaenv( 1, 'CUNGQR', ' ', n, n, n, -1 )
394 nb = max( nb1, nb2, nb3 )
395 lopt = max( 2*n, n*(nb+1) )
396 work( 1 ) = lopt
397 END IF
398*
399 IF( info.NE.0 ) THEN
400 CALL xerbla( 'CGEGV ', -info )
401 RETURN
402 ELSE IF( lquery ) THEN
403 RETURN
404 END IF
405*
406* Quick return if possible
407*
408 IF( n.EQ.0 )
409 $ RETURN
410*
411* Get machine constants
412*
413 eps = slamch( 'E' )*slamch( 'B' )
414 safmin = slamch( 'S' )
415 safmin = safmin + safmin
416 safmax = one / safmin
417*
418* Scale A
419*
420 anrm = clange( 'M', n, n, a, lda, rwork )
421 anrm1 = anrm
422 anrm2 = one
423 IF( anrm.LT.one ) THEN
424 IF( safmax*anrm.LT.one ) THEN
425 anrm1 = safmin
426 anrm2 = safmax*anrm
427 END IF
428 END IF
429*
430 IF( anrm.GT.zero ) THEN
431 CALL clascl( 'G', -1, -1, anrm, one, n, n, a, lda, iinfo )
432 IF( iinfo.NE.0 ) THEN
433 info = n + 10
434 RETURN
435 END IF
436 END IF
437*
438* Scale B
439*
440 bnrm = clange( 'M', n, n, b, ldb, rwork )
441 bnrm1 = bnrm
442 bnrm2 = one
443 IF( bnrm.LT.one ) THEN
444 IF( safmax*bnrm.LT.one ) THEN
445 bnrm1 = safmin
446 bnrm2 = safmax*bnrm
447 END IF
448 END IF
449*
450 IF( bnrm.GT.zero ) THEN
451 CALL clascl( 'G', -1, -1, bnrm, one, n, n, b, ldb, iinfo )
452 IF( iinfo.NE.0 ) THEN
453 info = n + 10
454 RETURN
455 END IF
456 END IF
457*
458* Permute the matrix to make it more nearly triangular
459* Also "balance" the matrix.
460*
461 ileft = 1
462 iright = n + 1
463 irwork = iright + n
464 CALL cggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
465 $ rwork( iright ), rwork( irwork ), iinfo )
466 IF( iinfo.NE.0 ) THEN
467 info = n + 1
468 GO TO 80
469 END IF
470*
471* Reduce B to triangular form, and initialize VL and/or VR
472*
473 irows = ihi + 1 - ilo
474 IF( ilv ) THEN
475 icols = n + 1 - ilo
476 ELSE
477 icols = irows
478 END IF
479 itau = 1
480 iwork = itau + irows
481 CALL cgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
482 $ work( iwork ), lwork+1-iwork, iinfo )
483 IF( iinfo.GE.0 )
484 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
485 IF( iinfo.NE.0 ) THEN
486 info = n + 2
487 GO TO 80
488 END IF
489*
490 CALL cunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
491 $ work( itau ), a( ilo, ilo ), lda, work( iwork ),
492 $ lwork+1-iwork, iinfo )
493 IF( iinfo.GE.0 )
494 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
495 IF( iinfo.NE.0 ) THEN
496 info = n + 3
497 GO TO 80
498 END IF
499*
500 IF( ilvl ) THEN
501 CALL claset( 'Full', n, n, czero, cone, vl, ldvl )
502 CALL clacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
503 $ vl( ilo+1, ilo ), ldvl )
504 CALL cungqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
505 $ work( itau ), work( iwork ), lwork+1-iwork,
506 $ iinfo )
507 IF( iinfo.GE.0 )
508 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
509 IF( iinfo.NE.0 ) THEN
510 info = n + 4
511 GO TO 80
512 END IF
513 END IF
514*
515 IF( ilvr )
516 $ CALL claset( 'Full', n, n, czero, cone, vr, ldvr )
517*
518* Reduce to generalized Hessenberg form
519*
520 IF( ilv ) THEN
521*
522* Eigenvectors requested -- work on whole matrix.
523*
524 CALL cgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
525 $ ldvl, vr, ldvr, iinfo )
526 ELSE
527 CALL cgghrd( 'N', 'N', irows, 1, irows, a( ilo, ilo ), lda,
528 $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, iinfo )
529 END IF
530 IF( iinfo.NE.0 ) THEN
531 info = n + 5
532 GO TO 80
533 END IF
534*
535* Perform QZ algorithm
536*
537 iwork = itau
538 IF( ilv ) THEN
539 chtemp = 'S'
540 ELSE
541 chtemp = 'E'
542 END IF
543 CALL chgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
544 $ alpha, beta, vl, ldvl, vr, ldvr, work( iwork ),
545 $ lwork+1-iwork, rwork( irwork ), iinfo )
546 IF( iinfo.GE.0 )
547 $ lwkopt = max( lwkopt, int( work( iwork ) )+iwork-1 )
548 IF( iinfo.NE.0 ) THEN
549 IF( iinfo.GT.0 .AND. iinfo.LE.n ) THEN
550 info = iinfo
551 ELSE IF( iinfo.GT.n .AND. iinfo.LE.2*n ) THEN
552 info = iinfo - n
553 ELSE
554 info = n + 6
555 END IF
556 GO TO 80
557 END IF
558*
559 IF( ilv ) THEN
560*
561* Compute Eigenvectors
562*
563 IF( ilvl ) THEN
564 IF( ilvr ) THEN
565 chtemp = 'B'
566 ELSE
567 chtemp = 'L'
568 END IF
569 ELSE
570 chtemp = 'R'
571 END IF
572*
573 CALL ctgevc( chtemp, 'B', ldumma, n, a, lda, b, ldb, vl, ldvl,
574 $ vr, ldvr, n, in, work( iwork ), rwork( irwork ),
575 $ iinfo )
576 IF( iinfo.NE.0 ) THEN
577 info = n + 7
578 GO TO 80
579 END IF
580*
581* Undo balancing on VL and VR, rescale
582*
583 IF( ilvl ) THEN
584 CALL cggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
585 $ rwork( iright ), n, vl, ldvl, iinfo )
586 IF( iinfo.NE.0 ) THEN
587 info = n + 8
588 GO TO 80
589 END IF
590 DO 30 jc = 1, n
591 temp = zero
592 DO 10 jr = 1, n
593 temp = max( temp, abs1( vl( jr, jc ) ) )
594 10 CONTINUE
595 IF( temp.LT.safmin )
596 $ GO TO 30
597 temp = one / temp
598 DO 20 jr = 1, n
599 vl( jr, jc ) = vl( jr, jc )*temp
600 20 CONTINUE
601 30 CONTINUE
602 END IF
603 IF( ilvr ) THEN
604 CALL cggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
605 $ rwork( iright ), n, vr, ldvr, iinfo )
606 IF( iinfo.NE.0 ) THEN
607 info = n + 9
608 GO TO 80
609 END IF
610 DO 60 jc = 1, n
611 temp = zero
612 DO 40 jr = 1, n
613 temp = max( temp, abs1( vr( jr, jc ) ) )
614 40 CONTINUE
615 IF( temp.LT.safmin )
616 $ GO TO 60
617 temp = one / temp
618 DO 50 jr = 1, n
619 vr( jr, jc ) = vr( jr, jc )*temp
620 50 CONTINUE
621 60 CONTINUE
622 END IF
623*
624* End of eigenvector calculation
625*
626 END IF
627*
628* Undo scaling in alpha, beta
629*
630* Note: this does not give the alpha and beta for the unscaled
631* problem.
632*
633* Un-scaling is limited to avoid underflow in alpha and beta
634* if they are significant.
635*
636 DO 70 jc = 1, n
637 absar = abs( real( alpha( jc ) ) )
638 absai = abs( aimag( alpha( jc ) ) )
639 absb = abs( real( beta( jc ) ) )
640 salfar = anrm*real( alpha( jc ) )
641 salfai = anrm*aimag( alpha( jc ) )
642 sbeta = bnrm*real( beta( jc ) )
643 ilimit = .false.
644 scale = one
645*
646* Check for significant underflow in imaginary part of ALPHA
647*
648 IF( abs( salfai ).LT.safmin .AND. absai.GE.
649 $ max( safmin, eps*absar, eps*absb ) ) THEN
650 ilimit = .true.
651 scale = ( safmin / anrm1 ) / max( safmin, anrm2*absai )
652 END IF
653*
654* Check for significant underflow in real part of ALPHA
655*
656 IF( abs( salfar ).LT.safmin .AND. absar.GE.
657 $ max( safmin, eps*absai, eps*absb ) ) THEN
658 ilimit = .true.
659 scale = max( scale, ( safmin / anrm1 ) /
660 $ max( safmin, anrm2*absar ) )
661 END IF
662*
663* Check for significant underflow in BETA
664*
665 IF( abs( sbeta ).LT.safmin .AND. absb.GE.
666 $ max( safmin, eps*absar, eps*absai ) ) THEN
667 ilimit = .true.
668 scale = max( scale, ( safmin / bnrm1 ) /
669 $ max( safmin, bnrm2*absb ) )
670 END IF
671*
672* Check for possible overflow when limiting scaling
673*
674 IF( ilimit ) THEN
675 temp = ( scale*safmin )*max( abs( salfar ), abs( salfai ),
676 $ abs( sbeta ) )
677 IF( temp.GT.one )
678 $ scale = scale / temp
679 IF( scale.LT.one )
680 $ ilimit = .false.
681 END IF
682*
683* Recompute un-scaled ALPHA, BETA if necessary.
684*
685 IF( ilimit ) THEN
686 salfar = ( scale*real( alpha( jc ) ) )*anrm
687 salfai = ( scale*aimag( alpha( jc ) ) )*anrm
688 sbeta = ( scale*beta( jc ) )*bnrm
689 END IF
690 alpha( jc ) = cmplx( salfar, salfai )
691 beta( jc ) = sbeta
692 70 CONTINUE
693*
694 80 CONTINUE
695 work( 1 ) = lwkopt
696*
697 RETURN
698*
699* End of CGEGV
700*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgeqrf(m, n, a, lda, tau, work, lwork, info)
CGEQRF
Definition cgeqrf.f:144
subroutine cggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
CGGBAK
Definition cggbak.f:147
subroutine cggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
CGGBAL
Definition cggbal.f:175
subroutine cgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
CGGHRD
Definition cgghrd.f:203
subroutine chgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
CHGEQZ
Definition chgeqz.f:283
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clange(norm, m, n, a, lda, work)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition clange.f:113
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition clascl.f:142
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:104
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine ctgevc(side, howmny, select, n, s, lds, p, ldp, vl, ldvl, vr, ldvr, mm, m, work, rwork, info)
CTGEVC
Definition ctgevc.f:217
subroutine cungqr(m, n, k, a, lda, tau, work, lwork, info)
CUNGQR
Definition cungqr.f:126
subroutine cunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
CUNMQR
Definition cunmqr.f:166
Here is the call graph for this function: