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

◆ cgeevx()

subroutine cgeevx ( character  balanc,
character  jobvl,
character  jobvr,
character  sense,
integer  n,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( * )  w,
complex, dimension( ldvl, * )  vl,
integer  ldvl,
complex, dimension( ldvr, * )  vr,
integer  ldvr,
integer  ilo,
integer  ihi,
real, dimension( * )  scale,
real  abnrm,
real, dimension( * )  rconde,
real, dimension( * )  rcondv,
complex, dimension( * )  work,
integer  lwork,
real, dimension( * )  rwork,
integer  info 
)

CGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices

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

Purpose:
 CGEEVX computes for an N-by-N complex nonsymmetric matrix A, the
 eigenvalues and, optionally, the left and/or right eigenvectors.

 Optionally also, it computes a balancing transformation to improve
 the conditioning of the eigenvalues and eigenvectors (ILO, IHI,
 SCALE, and ABNRM), reciprocal condition numbers for the eigenvalues
 (RCONDE), and reciprocal condition numbers for the right
 eigenvectors (RCONDV).

 The right eigenvector v(j) of A satisfies
                  A * v(j) = lambda(j) * v(j)
 where lambda(j) is its eigenvalue.
 The left eigenvector u(j) of A satisfies
               u(j)**H * A = lambda(j) * u(j)**H
 where u(j)**H denotes the conjugate transpose of u(j).

 The computed eigenvectors are normalized to have Euclidean norm
 equal to 1 and largest component real.

 Balancing a matrix means permuting the rows and columns to make it
 more nearly upper triangular, and applying a diagonal similarity
 transformation D * A * D**(-1), where D is a diagonal matrix, to
 make its rows and columns closer in norm and the condition numbers
 of its eigenvalues and eigenvectors smaller.  The computed
 reciprocal condition numbers correspond to the balanced matrix.
 Permuting rows and columns will not change the condition numbers
 (in exact arithmetic) but diagonal scaling will.  For further
 explanation of balancing, see section 4.10.2 of the LAPACK
 Users' Guide.
Parameters
[in]BALANC
          BALANC is CHARACTER*1
          Indicates how the input matrix should be diagonally scaled
          and/or permuted to improve the conditioning of its
          eigenvalues.
          = 'N': Do not diagonally scale or permute;
          = 'P': Perform permutations to make the matrix more nearly
                 upper triangular. Do not diagonally scale;
          = 'S': Diagonally scale the matrix, ie. replace A by
                 D*A*D**(-1), where D is a diagonal matrix chosen
                 to make the rows and columns of A more equal in
                 norm. Do not permute;
          = 'B': Both diagonally scale and permute A.

          Computed reciprocal condition numbers will be for the matrix
          after balancing and/or permuting. Permuting does not change
          condition numbers (in exact arithmetic), but balancing does.
[in]JOBVL
          JOBVL is CHARACTER*1
          = 'N': left eigenvectors of A are not computed;
          = 'V': left eigenvectors of A are computed.
          If SENSE = 'E' or 'B', JOBVL must = 'V'.
[in]JOBVR
          JOBVR is CHARACTER*1
          = 'N': right eigenvectors of A are not computed;
          = 'V': right eigenvectors of A are computed.
          If SENSE = 'E' or 'B', JOBVR must = 'V'.
[in]SENSE
          SENSE is CHARACTER*1
          Determines which reciprocal condition numbers are computed.
          = 'N': None are computed;
          = 'E': Computed for eigenvalues only;
          = 'V': Computed for right eigenvectors only;
          = 'B': Computed for eigenvalues and right eigenvectors.

          If SENSE = 'E' or 'B', both left and right eigenvectors
          must also be computed (JOBVL = 'V' and JOBVR = 'V').
[in]N
          N is INTEGER
          The order of the matrix A. N >= 0.
[in,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the N-by-N matrix A.
          On exit, A has been overwritten.  If JOBVL = 'V' or
          JOBVR = 'V', A contains the Schur form of the balanced
          version of the matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]W
          W is COMPLEX array, dimension (N)
          W contains the computed eigenvalues.
[out]VL
          VL is COMPLEX array, dimension (LDVL,N)
          If JOBVL = 'V', the left eigenvectors u(j) are stored one
          after another in the columns of VL, in the same order
          as their eigenvalues.
          If JOBVL = 'N', VL is not referenced.
          u(j) = VL(:,j), the j-th column of VL.
[in]LDVL
          LDVL is INTEGER
          The leading dimension of the array VL.  LDVL >= 1; if
          JOBVL = 'V', LDVL >= N.
[out]VR
          VR is COMPLEX array, dimension (LDVR,N)
          If JOBVR = 'V', the right eigenvectors v(j) are stored one
          after another in the columns of VR, in the same order
          as their eigenvalues.
          If JOBVR = 'N', VR is not referenced.
          v(j) = VR(:,j), the j-th column of VR.
[in]LDVR
          LDVR is INTEGER
          The leading dimension of the array VR.  LDVR >= 1; if
          JOBVR = 'V', LDVR >= N.
[out]ILO
          ILO is INTEGER
[out]IHI
          IHI is INTEGER
          ILO and IHI are integer values determined when A was
          balanced.  The balanced A(i,j) = 0 if I > J and
          J = 1,...,ILO-1 or I = IHI+1,...,N.
[out]SCALE
          SCALE is REAL array, dimension (N)
          Details of the permutations and scaling factors applied
          when balancing A.  If P(j) is the index of the row and column
          interchanged with row and column j, and D(j) is the scaling
          factor applied to row and column j, then
          SCALE(J) = P(J),    for J = 1,...,ILO-1
                   = D(J),    for J = ILO,...,IHI
                   = P(J)     for J = IHI+1,...,N.
          The order in which the interchanges are made is N to IHI+1,
          then 1 to ILO-1.
[out]ABNRM
          ABNRM is REAL
          The one-norm of the balanced matrix (the maximum
          of the sum of absolute values of elements of any column).
[out]RCONDE
          RCONDE is REAL array, dimension (N)
          RCONDE(j) is the reciprocal condition number of the j-th
          eigenvalue.
[out]RCONDV
          RCONDV is REAL array, dimension (N)
          RCONDV(j) is the reciprocal condition number of the j-th
          right eigenvector.
[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.  If SENSE = 'N' or 'E',
          LWORK >= max(1,2*N), and if SENSE = 'V' or 'B',
          LWORK >= N*N+2*N.
          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]RWORK
          RWORK is REAL array, dimension (2*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, the QR algorithm failed to compute all the
                eigenvalues, and no eigenvectors or condition numbers
                have been computed; elements 1:ILO-1 and i+1:N of W
                contain eigenvalues which have converged.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 285 of file cgeevx.f.

288 implicit none
289*
290* -- LAPACK driver routine --
291* -- LAPACK is a software package provided by Univ. of Tennessee, --
292* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
293*
294* .. Scalar Arguments ..
295 CHARACTER BALANC, JOBVL, JOBVR, SENSE
296 INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
297 REAL ABNRM
298* ..
299* .. Array Arguments ..
300 REAL RCONDE( * ), RCONDV( * ), RWORK( * ),
301 $ SCALE( * )
302 COMPLEX A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ),
303 $ W( * ), WORK( * )
304* ..
305*
306* =====================================================================
307*
308* .. Parameters ..
309 REAL ZERO, ONE
310 parameter( zero = 0.0e0, one = 1.0e0 )
311* ..
312* .. Local Scalars ..
313 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR, WNTSNB, WNTSNE,
314 $ WNTSNN, WNTSNV
315 CHARACTER JOB, SIDE
316 INTEGER HSWORK, I, ICOND, IERR, ITAU, IWRK, K,
317 $ LWORK_TREVC, MAXWRK, MINWRK, NOUT
318 REAL ANRM, BIGNUM, CSCALE, EPS, SCL, SMLNUM
319 COMPLEX TMP
320* ..
321* .. Local Arrays ..
322 LOGICAL SELECT( 1 )
323 REAL DUM( 1 )
324* ..
325* .. External Subroutines ..
326 EXTERNAL slascl, xerbla, csscal, cgebak, cgebal,
328 $ ctrsna, cunghr
329* ..
330* .. External Functions ..
331 LOGICAL LSAME
332 INTEGER ISAMAX, ILAENV
333 REAL SLAMCH, SCNRM2, CLANGE, SROUNDUP_LWORK
334 EXTERNAL lsame, isamax, ilaenv, slamch, scnrm2, clange,
336* ..
337* .. Intrinsic Functions ..
338 INTRINSIC real, cmplx, conjg, aimag, max, sqrt
339* ..
340* .. Executable Statements ..
341*
342* Test the input arguments
343*
344 info = 0
345 lquery = ( lwork.EQ.-1 )
346 wantvl = lsame( jobvl, 'V' )
347 wantvr = lsame( jobvr, 'V' )
348 wntsnn = lsame( sense, 'N' )
349 wntsne = lsame( sense, 'E' )
350 wntsnv = lsame( sense, 'V' )
351 wntsnb = lsame( sense, 'B' )
352 IF( .NOT.( lsame( balanc, 'N' ) .OR. lsame( balanc, 'S' ) .OR.
353 $ lsame( balanc, 'P' ) .OR. lsame( balanc, 'B' ) ) ) THEN
354 info = -1
355 ELSE IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
356 info = -2
357 ELSE IF( ( .NOT.wantvr ) .AND. ( .NOT.lsame( jobvr, 'N' ) ) ) THEN
358 info = -3
359 ELSE IF( .NOT.( wntsnn .OR. wntsne .OR. wntsnb .OR. wntsnv ) .OR.
360 $ ( ( wntsne .OR. wntsnb ) .AND. .NOT.( wantvl .AND.
361 $ wantvr ) ) ) THEN
362 info = -4
363 ELSE IF( n.LT.0 ) THEN
364 info = -5
365 ELSE IF( lda.LT.max( 1, n ) ) THEN
366 info = -7
367 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) ) THEN
368 info = -10
369 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) ) THEN
370 info = -12
371 END IF
372*
373* Compute workspace
374* (Note: Comments in the code beginning "Workspace:" describe the
375* minimal amount of workspace needed at that point in the code,
376* as well as the preferred amount for good performance.
377* CWorkspace refers to complex workspace, and RWorkspace to real
378* workspace. NB refers to the optimal block size for the
379* immediately following subroutine, as returned by ILAENV.
380* HSWORK refers to the workspace preferred by CHSEQR, as
381* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
382* the worst case.)
383*
384 IF( info.EQ.0 ) THEN
385 IF( n.EQ.0 ) THEN
386 minwrk = 1
387 maxwrk = 1
388 ELSE
389 maxwrk = n + n*ilaenv( 1, 'CGEHRD', ' ', n, 1, n, 0 )
390*
391 IF( wantvl ) THEN
392 CALL ctrevc3( 'L', 'B', SELECT, n, a, lda,
393 $ vl, ldvl, vr, ldvr,
394 $ n, nout, work, -1, rwork, -1, ierr )
395 lwork_trevc = int( work(1) )
396 maxwrk = max( maxwrk, lwork_trevc )
397 CALL chseqr( 'S', 'V', n, 1, n, a, lda, w, vl, ldvl,
398 $ work, -1, info )
399 ELSE IF( wantvr ) THEN
400 CALL ctrevc3( 'R', 'B', SELECT, n, a, lda,
401 $ vl, ldvl, vr, ldvr,
402 $ n, nout, work, -1, rwork, -1, ierr )
403 lwork_trevc = int( work(1) )
404 maxwrk = max( maxwrk, lwork_trevc )
405 CALL chseqr( 'S', 'V', n, 1, n, a, lda, w, vr, ldvr,
406 $ work, -1, info )
407 ELSE
408 IF( wntsnn ) THEN
409 CALL chseqr( 'E', 'N', n, 1, n, a, lda, w, vr, ldvr,
410 $ work, -1, info )
411 ELSE
412 CALL chseqr( 'S', 'N', n, 1, n, a, lda, w, vr, ldvr,
413 $ work, -1, info )
414 END IF
415 END IF
416 hswork = int( work(1) )
417*
418 IF( ( .NOT.wantvl ) .AND. ( .NOT.wantvr ) ) THEN
419 minwrk = 2*n
420 IF( .NOT.( wntsnn .OR. wntsne ) )
421 $ minwrk = max( minwrk, n*n + 2*n )
422 maxwrk = max( maxwrk, hswork )
423 IF( .NOT.( wntsnn .OR. wntsne ) )
424 $ maxwrk = max( maxwrk, n*n + 2*n )
425 ELSE
426 minwrk = 2*n
427 IF( .NOT.( wntsnn .OR. wntsne ) )
428 $ minwrk = max( minwrk, n*n + 2*n )
429 maxwrk = max( maxwrk, hswork )
430 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1, 'CUNGHR',
431 $ ' ', n, 1, n, -1 ) )
432 IF( .NOT.( wntsnn .OR. wntsne ) )
433 $ maxwrk = max( maxwrk, n*n + 2*n )
434 maxwrk = max( maxwrk, 2*n )
435 END IF
436 maxwrk = max( maxwrk, minwrk )
437 END IF
438 work( 1 ) = sroundup_lwork(maxwrk)
439*
440 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
441 info = -20
442 END IF
443 END IF
444*
445 IF( info.NE.0 ) THEN
446 CALL xerbla( 'CGEEVX', -info )
447 RETURN
448 ELSE IF( lquery ) THEN
449 RETURN
450 END IF
451*
452* Quick return if possible
453*
454 IF( n.EQ.0 )
455 $ RETURN
456*
457* Get machine constants
458*
459 eps = slamch( 'P' )
460 smlnum = slamch( 'S' )
461 bignum = one / smlnum
462 smlnum = sqrt( smlnum ) / eps
463 bignum = one / smlnum
464*
465* Scale A if max element outside range [SMLNUM,BIGNUM]
466*
467 icond = 0
468 anrm = clange( 'M', n, n, a, lda, dum )
469 scalea = .false.
470 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
471 scalea = .true.
472 cscale = smlnum
473 ELSE IF( anrm.GT.bignum ) THEN
474 scalea = .true.
475 cscale = bignum
476 END IF
477 IF( scalea )
478 $ CALL clascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
479*
480* Balance the matrix and compute ABNRM
481*
482 CALL cgebal( balanc, n, a, lda, ilo, ihi, scale, ierr )
483 abnrm = clange( '1', n, n, a, lda, dum )
484 IF( scalea ) THEN
485 dum( 1 ) = abnrm
486 CALL slascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1, ierr )
487 abnrm = dum( 1 )
488 END IF
489*
490* Reduce to upper Hessenberg form
491* (CWorkspace: need 2*N, prefer N+N*NB)
492* (RWorkspace: none)
493*
494 itau = 1
495 iwrk = itau + n
496 CALL cgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
497 $ lwork-iwrk+1, ierr )
498*
499 IF( wantvl ) THEN
500*
501* Want left eigenvectors
502* Copy Householder vectors to VL
503*
504 side = 'L'
505 CALL clacpy( 'L', n, n, a, lda, vl, ldvl )
506*
507* Generate unitary matrix in VL
508* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
509* (RWorkspace: none)
510*
511 CALL cunghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
512 $ lwork-iwrk+1, ierr )
513*
514* Perform QR iteration, accumulating Schur vectors in VL
515* (CWorkspace: need 1, prefer HSWORK (see comments) )
516* (RWorkspace: none)
517*
518 iwrk = itau
519 CALL chseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vl, ldvl,
520 $ work( iwrk ), lwork-iwrk+1, info )
521*
522 IF( wantvr ) THEN
523*
524* Want left and right eigenvectors
525* Copy Schur vectors to VR
526*
527 side = 'B'
528 CALL clacpy( 'F', n, n, vl, ldvl, vr, ldvr )
529 END IF
530*
531 ELSE IF( wantvr ) THEN
532*
533* Want right eigenvectors
534* Copy Householder vectors to VR
535*
536 side = 'R'
537 CALL clacpy( 'L', n, n, a, lda, vr, ldvr )
538*
539* Generate unitary matrix in VR
540* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
541* (RWorkspace: none)
542*
543 CALL cunghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
544 $ lwork-iwrk+1, ierr )
545*
546* Perform QR iteration, accumulating Schur vectors in VR
547* (CWorkspace: need 1, prefer HSWORK (see comments) )
548* (RWorkspace: none)
549*
550 iwrk = itau
551 CALL chseqr( 'S', 'V', n, ilo, ihi, a, lda, w, vr, ldvr,
552 $ work( iwrk ), lwork-iwrk+1, info )
553*
554 ELSE
555*
556* Compute eigenvalues only
557* If condition numbers desired, compute Schur form
558*
559 IF( wntsnn ) THEN
560 job = 'E'
561 ELSE
562 job = 'S'
563 END IF
564*
565* (CWorkspace: need 1, prefer HSWORK (see comments) )
566* (RWorkspace: none)
567*
568 iwrk = itau
569 CALL chseqr( job, 'N', n, ilo, ihi, a, lda, w, vr, ldvr,
570 $ work( iwrk ), lwork-iwrk+1, info )
571 END IF
572*
573* If INFO .NE. 0 from CHSEQR, then quit
574*
575 IF( info.NE.0 )
576 $ GO TO 50
577*
578 IF( wantvl .OR. wantvr ) THEN
579*
580* Compute left and/or right eigenvectors
581* (CWorkspace: need 2*N, prefer N + 2*N*NB)
582* (RWorkspace: need N)
583*
584 CALL ctrevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
585 $ n, nout, work( iwrk ), lwork-iwrk+1,
586 $ rwork, n, ierr )
587 END IF
588*
589* Compute condition numbers if desired
590* (CWorkspace: need N*N+2*N unless SENSE = 'E')
591* (RWorkspace: need 2*N unless SENSE = 'E')
592*
593 IF( .NOT.wntsnn ) THEN
594 CALL ctrsna( sense, 'A', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
595 $ rconde, rcondv, n, nout, work( iwrk ), n, rwork,
596 $ icond )
597 END IF
598*
599 IF( wantvl ) THEN
600*
601* Undo balancing of left eigenvectors
602*
603 CALL cgebak( balanc, 'L', n, ilo, ihi, scale, n, vl, ldvl,
604 $ ierr )
605*
606* Normalize left eigenvectors and make largest component real
607*
608 DO 20 i = 1, n
609 scl = one / scnrm2( n, vl( 1, i ), 1 )
610 CALL csscal( n, scl, vl( 1, i ), 1 )
611 DO 10 k = 1, n
612 rwork( k ) = real( vl( k, i ) )**2 +
613 $ aimag( vl( k, i ) )**2
614 10 CONTINUE
615 k = isamax( n, rwork, 1 )
616 tmp = conjg( vl( k, i ) ) / sqrt( rwork( k ) )
617 CALL cscal( n, tmp, vl( 1, i ), 1 )
618 vl( k, i ) = cmplx( real( vl( k, i ) ), zero )
619 20 CONTINUE
620 END IF
621*
622 IF( wantvr ) THEN
623*
624* Undo balancing of right eigenvectors
625*
626 CALL cgebak( balanc, 'R', n, ilo, ihi, scale, n, vr, ldvr,
627 $ ierr )
628*
629* Normalize right eigenvectors and make largest component real
630*
631 DO 40 i = 1, n
632 scl = one / scnrm2( n, vr( 1, i ), 1 )
633 CALL csscal( n, scl, vr( 1, i ), 1 )
634 DO 30 k = 1, n
635 rwork( k ) = real( vr( k, i ) )**2 +
636 $ aimag( vr( k, i ) )**2
637 30 CONTINUE
638 k = isamax( n, rwork, 1 )
639 tmp = conjg( vr( k, i ) ) / sqrt( rwork( k ) )
640 CALL cscal( n, tmp, vr( 1, i ), 1 )
641 vr( k, i ) = cmplx( real( vr( k, i ) ), zero )
642 40 CONTINUE
643 END IF
644*
645* Undo scaling if necessary
646*
647 50 CONTINUE
648 IF( scalea ) THEN
649 CALL clascl( 'G', 0, 0, cscale, anrm, n-info, 1, w( info+1 ),
650 $ max( n-info, 1 ), ierr )
651 IF( info.EQ.0 ) THEN
652 IF( ( wntsnv .OR. wntsnb ) .AND. icond.EQ.0 )
653 $ CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, rcondv, n,
654 $ ierr )
655 ELSE
656 CALL clascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, w, n, ierr )
657 END IF
658 END IF
659*
660 work( 1 ) = sroundup_lwork(maxwrk)
661 RETURN
662*
663* End of CGEEVX
664*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
CGEBAK
Definition cgebak.f:131
subroutine cgebal(job, n, a, lda, ilo, ihi, scale, info)
CGEBAL
Definition cgebal.f:165
subroutine cgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
CGEHRD
Definition cgehrd.f:167
subroutine chseqr(job, compz, n, ilo, ihi, h, ldh, w, z, ldz, work, lwork, info)
CHSEQR
Definition chseqr.f:299
integer function isamax(n, sx, incx)
ISAMAX
Definition isamax.f:71
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
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:115
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:143
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function scnrm2(n, x, incx)
SCNRM2
Definition scnrm2.f90:90
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine csscal(n, sa, cx, incx)
CSSCAL
Definition csscal.f:78
subroutine cscal(n, ca, cx, incx)
CSCAL
Definition cscal.f:78
subroutine ctrevc3(side, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, mm, m, work, lwork, rwork, lrwork, info)
CTREVC3
Definition ctrevc3.f:244
subroutine ctrsna(job, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, s, sep, mm, m, work, ldwork, rwork, info)
CTRSNA
Definition ctrsna.f:249
subroutine cunghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
CUNGHR
Definition cunghr.f:126
Here is the call graph for this function:
Here is the caller graph for this function: