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

◆ dgeevx()

subroutine dgeevx ( character  balanc,
character  jobvl,
character  jobvr,
character  sense,
integer  n,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( * )  wr,
double precision, dimension( * )  wi,
double precision, dimension( ldvl, * )  vl,
integer  ldvl,
double precision, dimension( ldvr, * )  vr,
integer  ldvr,
integer  ilo,
integer  ihi,
double precision, dimension( * )  scale,
double precision  abnrm,
double precision, dimension( * )  rconde,
double precision, dimension( * )  rcondv,
double precision, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  info 
)

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

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

Purpose:
 DGEEVX computes for an N-by-N real 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, i.e. 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 DOUBLE PRECISION 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 real Schur form of the balanced
          version of the input matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]WR
          WR is DOUBLE PRECISION array, dimension (N)
[out]WI
          WI is DOUBLE PRECISION array, dimension (N)
          WR and WI contain the real and imaginary parts,
          respectively, of the computed eigenvalues.  Complex
          conjugate pairs of eigenvalues will appear consecutively
          with the eigenvalue having the positive imaginary part
          first.
[out]VL
          VL is DOUBLE PRECISION 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.
          If the j-th eigenvalue is real, then u(j) = VL(:,j),
          the j-th column of VL.
          If the j-th and (j+1)-st eigenvalues form a complex
          conjugate pair, then u(j) = VL(:,j) + i*VL(:,j+1) and
          u(j+1) = VL(:,j) - i*VL(:,j+1).
[in]LDVL
          LDVL is INTEGER
          The leading dimension of the array VL.  LDVL >= 1; if
          JOBVL = 'V', LDVL >= N.
[out]VR
          VR is DOUBLE PRECISION 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.
          If the j-th eigenvalue is real, then v(j) = VR(:,j),
          the j-th column of VR.
          If the j-th and (j+1)-st eigenvalues form a complex
          conjugate pair, then v(j) = VR(:,j) + i*VR(:,j+1) and
          v(j+1) = VR(:,j) - i*VR(:,j+1).
[in]LDVR
          LDVR is INTEGER
          The leading dimension of the array VR.  LDVR >= 1, and 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 DOUBLE PRECISION 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 DOUBLE PRECISION
          The one-norm of the balanced matrix (the maximum
          of the sum of absolute values of elements of any column).
[out]RCONDE
          RCONDE is DOUBLE PRECISION array, dimension (N)
          RCONDE(j) is the reciprocal condition number of the j-th
          eigenvalue.
[out]RCONDV
          RCONDV is DOUBLE PRECISION array, dimension (N)
          RCONDV(j) is the reciprocal condition number of the j-th
          right eigenvector.
[out]WORK
          WORK is DOUBLE PRECISION 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 JOBVL = 'V' or JOBVR = 'V',
          LWORK >= 3*N.  If SENSE = 'V' or 'B', LWORK >= N*(N+6).
          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]IWORK
          IWORK is INTEGER array, dimension (2*N-2)
          If SENSE = 'N' or 'E', not referenced.
[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 WR
                and WI contain eigenvalues which have converged.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 303 of file dgeevx.f.

306 implicit none
307*
308* -- LAPACK driver routine --
309* -- LAPACK is a software package provided by Univ. of Tennessee, --
310* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
311*
312* .. Scalar Arguments ..
313 CHARACTER BALANC, JOBVL, JOBVR, SENSE
314 INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N
315 DOUBLE PRECISION ABNRM
316* ..
317* .. Array Arguments ..
318 INTEGER IWORK( * )
319 DOUBLE PRECISION A( LDA, * ), RCONDE( * ), RCONDV( * ),
320 $ SCALE( * ), VL( LDVL, * ), VR( LDVR, * ),
321 $ WI( * ), WORK( * ), WR( * )
322* ..
323*
324* =====================================================================
325*
326* .. Parameters ..
327 DOUBLE PRECISION ZERO, ONE
328 parameter( zero = 0.0d0, one = 1.0d0 )
329* ..
330* .. Local Scalars ..
331 LOGICAL LQUERY, SCALEA, WANTVL, WANTVR, WNTSNB, WNTSNE,
332 $ WNTSNN, WNTSNV
333 CHARACTER JOB, SIDE
334 INTEGER HSWORK, I, ICOND, IERR, ITAU, IWRK, K,
335 $ LWORK_TREVC, MAXWRK, MINWRK, NOUT
336 DOUBLE PRECISION ANRM, BIGNUM, CS, CSCALE, EPS, R, SCL, SMLNUM,
337 $ SN
338* ..
339* .. Local Arrays ..
340 LOGICAL SELECT( 1 )
341 DOUBLE PRECISION DUM( 1 )
342* ..
343* .. External Subroutines ..
344 EXTERNAL dgebak, dgebal, dgehrd, dhseqr, dlacpy, dlartg,
346 $ xerbla
347* ..
348* .. External Functions ..
349 LOGICAL LSAME
350 INTEGER IDAMAX, ILAENV
351 DOUBLE PRECISION DLAMCH, DLANGE, DLAPY2, DNRM2
352 EXTERNAL lsame, idamax, ilaenv, dlamch, dlange, dlapy2,
353 $ dnrm2
354* ..
355* .. Intrinsic Functions ..
356 INTRINSIC max, sqrt
357* ..
358* .. Executable Statements ..
359*
360* Test the input arguments
361*
362 info = 0
363 lquery = ( lwork.EQ.-1 )
364 wantvl = lsame( jobvl, 'V' )
365 wantvr = lsame( jobvr, 'V' )
366 wntsnn = lsame( sense, 'N' )
367 wntsne = lsame( sense, 'E' )
368 wntsnv = lsame( sense, 'V' )
369 wntsnb = lsame( sense, 'B' )
370 IF( .NOT.( lsame( balanc, 'N' ) .OR. lsame( balanc, 'S' )
371 $ .OR. lsame( balanc, 'P' ) .OR. lsame( balanc, 'B' ) ) )
372 $ THEN
373 info = -1
374 ELSE IF( ( .NOT.wantvl ) .AND. ( .NOT.lsame( jobvl, 'N' ) ) ) THEN
375 info = -2
376 ELSE IF( ( .NOT.wantvr ) .AND. ( .NOT.lsame( jobvr, 'N' ) ) ) THEN
377 info = -3
378 ELSE IF( .NOT.( wntsnn .OR. wntsne .OR. wntsnb .OR. wntsnv ) .OR.
379 $ ( ( wntsne .OR. wntsnb ) .AND. .NOT.( wantvl .AND.
380 $ wantvr ) ) ) THEN
381 info = -4
382 ELSE IF( n.LT.0 ) THEN
383 info = -5
384 ELSE IF( lda.LT.max( 1, n ) ) THEN
385 info = -7
386 ELSE IF( ldvl.LT.1 .OR. ( wantvl .AND. ldvl.LT.n ) ) THEN
387 info = -11
388 ELSE IF( ldvr.LT.1 .OR. ( wantvr .AND. ldvr.LT.n ) ) THEN
389 info = -13
390 END IF
391*
392* Compute workspace
393* (Note: Comments in the code beginning "Workspace:" describe the
394* minimal amount of workspace needed at that point in the code,
395* as well as the preferred amount for good performance.
396* NB refers to the optimal block size for the immediately
397* following subroutine, as returned by ILAENV.
398* HSWORK refers to the workspace preferred by DHSEQR, as
399* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
400* the worst case.)
401*
402 IF( info.EQ.0 ) THEN
403 IF( n.EQ.0 ) THEN
404 minwrk = 1
405 maxwrk = 1
406 ELSE
407 maxwrk = n + n*ilaenv( 1, 'DGEHRD', ' ', n, 1, n, 0 )
408*
409 IF( wantvl ) THEN
410 CALL dtrevc3( 'L', 'B', SELECT, n, a, lda,
411 $ vl, ldvl, vr, ldvr,
412 $ n, nout, work, -1, ierr )
413 lwork_trevc = int( work(1) )
414 maxwrk = max( maxwrk, n + lwork_trevc )
415 CALL dhseqr( 'S', 'V', n, 1, n, a, lda, wr, wi, vl, ldvl,
416 $ work, -1, info )
417 ELSE IF( wantvr ) THEN
418 CALL dtrevc3( 'R', 'B', SELECT, n, a, lda,
419 $ vl, ldvl, vr, ldvr,
420 $ n, nout, work, -1, ierr )
421 lwork_trevc = int( work(1) )
422 maxwrk = max( maxwrk, n + lwork_trevc )
423 CALL dhseqr( 'S', 'V', n, 1, n, a, lda, wr, wi, vr, ldvr,
424 $ work, -1, info )
425 ELSE
426 IF( wntsnn ) THEN
427 CALL dhseqr( 'E', 'N', n, 1, n, a, lda, wr, wi, vr,
428 $ ldvr, work, -1, info )
429 ELSE
430 CALL dhseqr( 'S', 'N', n, 1, n, a, lda, wr, wi, vr,
431 $ ldvr, work, -1, info )
432 END IF
433 END IF
434 hswork = int( work(1) )
435*
436 IF( ( .NOT.wantvl ) .AND. ( .NOT.wantvr ) ) THEN
437 minwrk = 2*n
438 IF( .NOT.wntsnn )
439 $ minwrk = max( minwrk, n*n+6*n )
440 maxwrk = max( maxwrk, hswork )
441 IF( .NOT.wntsnn )
442 $ maxwrk = max( maxwrk, n*n + 6*n )
443 ELSE
444 minwrk = 3*n
445 IF( ( .NOT.wntsnn ) .AND. ( .NOT.wntsne ) )
446 $ minwrk = max( minwrk, n*n + 6*n )
447 maxwrk = max( maxwrk, hswork )
448 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1, 'DORGHR',
449 $ ' ', n, 1, n, -1 ) )
450 IF( ( .NOT.wntsnn ) .AND. ( .NOT.wntsne ) )
451 $ maxwrk = max( maxwrk, n*n + 6*n )
452 maxwrk = max( maxwrk, 3*n )
453 END IF
454 maxwrk = max( maxwrk, minwrk )
455 END IF
456 work( 1 ) = maxwrk
457*
458 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
459 info = -21
460 END IF
461 END IF
462*
463 IF( info.NE.0 ) THEN
464 CALL xerbla( 'DGEEVX', -info )
465 RETURN
466 ELSE IF( lquery ) THEN
467 RETURN
468 END IF
469*
470* Quick return if possible
471*
472 IF( n.EQ.0 )
473 $ RETURN
474*
475* Get machine constants
476*
477 eps = dlamch( 'P' )
478 smlnum = dlamch( 'S' )
479 bignum = one / smlnum
480 smlnum = sqrt( smlnum ) / eps
481 bignum = one / smlnum
482*
483* Scale A if max element outside range [SMLNUM,BIGNUM]
484*
485 icond = 0
486 anrm = dlange( 'M', n, n, a, lda, dum )
487 scalea = .false.
488 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
489 scalea = .true.
490 cscale = smlnum
491 ELSE IF( anrm.GT.bignum ) THEN
492 scalea = .true.
493 cscale = bignum
494 END IF
495 IF( scalea )
496 $ CALL dlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
497*
498* Balance the matrix and compute ABNRM
499*
500 CALL dgebal( balanc, n, a, lda, ilo, ihi, scale, ierr )
501 abnrm = dlange( '1', n, n, a, lda, dum )
502 IF( scalea ) THEN
503 dum( 1 ) = abnrm
504 CALL dlascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1, ierr )
505 abnrm = dum( 1 )
506 END IF
507*
508* Reduce to upper Hessenberg form
509* (Workspace: need 2*N, prefer N+N*NB)
510*
511 itau = 1
512 iwrk = itau + n
513 CALL dgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
514 $ lwork-iwrk+1, ierr )
515*
516 IF( wantvl ) THEN
517*
518* Want left eigenvectors
519* Copy Householder vectors to VL
520*
521 side = 'L'
522 CALL dlacpy( 'L', n, n, a, lda, vl, ldvl )
523*
524* Generate orthogonal matrix in VL
525* (Workspace: need 2*N-1, prefer N+(N-1)*NB)
526*
527 CALL dorghr( n, ilo, ihi, vl, ldvl, work( itau ), work( iwrk ),
528 $ lwork-iwrk+1, ierr )
529*
530* Perform QR iteration, accumulating Schur vectors in VL
531* (Workspace: need 1, prefer HSWORK (see comments) )
532*
533 iwrk = itau
534 CALL dhseqr( 'S', 'V', n, ilo, ihi, a, lda, wr, wi, vl, ldvl,
535 $ work( iwrk ), lwork-iwrk+1, info )
536*
537 IF( wantvr ) THEN
538*
539* Want left and right eigenvectors
540* Copy Schur vectors to VR
541*
542 side = 'B'
543 CALL dlacpy( 'F', n, n, vl, ldvl, vr, ldvr )
544 END IF
545*
546 ELSE IF( wantvr ) THEN
547*
548* Want right eigenvectors
549* Copy Householder vectors to VR
550*
551 side = 'R'
552 CALL dlacpy( 'L', n, n, a, lda, vr, ldvr )
553*
554* Generate orthogonal matrix in VR
555* (Workspace: need 2*N-1, prefer N+(N-1)*NB)
556*
557 CALL dorghr( n, ilo, ihi, vr, ldvr, work( itau ), work( iwrk ),
558 $ lwork-iwrk+1, ierr )
559*
560* Perform QR iteration, accumulating Schur vectors in VR
561* (Workspace: need 1, prefer HSWORK (see comments) )
562*
563 iwrk = itau
564 CALL dhseqr( 'S', 'V', n, ilo, ihi, a, lda, wr, wi, vr, ldvr,
565 $ work( iwrk ), lwork-iwrk+1, info )
566*
567 ELSE
568*
569* Compute eigenvalues only
570* If condition numbers desired, compute Schur form
571*
572 IF( wntsnn ) THEN
573 job = 'E'
574 ELSE
575 job = 'S'
576 END IF
577*
578* (Workspace: need 1, prefer HSWORK (see comments) )
579*
580 iwrk = itau
581 CALL dhseqr( job, 'N', n, ilo, ihi, a, lda, wr, wi, vr, ldvr,
582 $ work( iwrk ), lwork-iwrk+1, info )
583 END IF
584*
585* If INFO .NE. 0 from DHSEQR, then quit
586*
587 IF( info.NE.0 )
588 $ GO TO 50
589*
590 IF( wantvl .OR. wantvr ) THEN
591*
592* Compute left and/or right eigenvectors
593* (Workspace: need 3*N, prefer N + 2*N*NB)
594*
595 CALL dtrevc3( side, 'B', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
596 $ n, nout, work( iwrk ), lwork-iwrk+1, ierr )
597 END IF
598*
599* Compute condition numbers if desired
600* (Workspace: need N*N+6*N unless SENSE = 'E')
601*
602 IF( .NOT.wntsnn ) THEN
603 CALL dtrsna( sense, 'A', SELECT, n, a, lda, vl, ldvl, vr, ldvr,
604 $ rconde, rcondv, n, nout, work( iwrk ), n, iwork,
605 $ icond )
606 END IF
607*
608 IF( wantvl ) THEN
609*
610* Undo balancing of left eigenvectors
611*
612 CALL dgebak( balanc, 'L', n, ilo, ihi, scale, n, vl, ldvl,
613 $ ierr )
614*
615* Normalize left eigenvectors and make largest component real
616*
617 DO 20 i = 1, n
618 IF( wi( i ).EQ.zero ) THEN
619 scl = one / dnrm2( n, vl( 1, i ), 1 )
620 CALL dscal( n, scl, vl( 1, i ), 1 )
621 ELSE IF( wi( i ).GT.zero ) THEN
622 scl = one / dlapy2( dnrm2( n, vl( 1, i ), 1 ),
623 $ dnrm2( n, vl( 1, i+1 ), 1 ) )
624 CALL dscal( n, scl, vl( 1, i ), 1 )
625 CALL dscal( n, scl, vl( 1, i+1 ), 1 )
626 DO 10 k = 1, n
627 work( k ) = vl( k, i )**2 + vl( k, i+1 )**2
628 10 CONTINUE
629 k = idamax( n, work, 1 )
630 CALL dlartg( vl( k, i ), vl( k, i+1 ), cs, sn, r )
631 CALL drot( n, vl( 1, i ), 1, vl( 1, i+1 ), 1, cs, sn )
632 vl( k, i+1 ) = zero
633 END IF
634 20 CONTINUE
635 END IF
636*
637 IF( wantvr ) THEN
638*
639* Undo balancing of right eigenvectors
640*
641 CALL dgebak( balanc, 'R', n, ilo, ihi, scale, n, vr, ldvr,
642 $ ierr )
643*
644* Normalize right eigenvectors and make largest component real
645*
646 DO 40 i = 1, n
647 IF( wi( i ).EQ.zero ) THEN
648 scl = one / dnrm2( n, vr( 1, i ), 1 )
649 CALL dscal( n, scl, vr( 1, i ), 1 )
650 ELSE IF( wi( i ).GT.zero ) THEN
651 scl = one / dlapy2( dnrm2( n, vr( 1, i ), 1 ),
652 $ dnrm2( n, vr( 1, i+1 ), 1 ) )
653 CALL dscal( n, scl, vr( 1, i ), 1 )
654 CALL dscal( n, scl, vr( 1, i+1 ), 1 )
655 DO 30 k = 1, n
656 work( k ) = vr( k, i )**2 + vr( k, i+1 )**2
657 30 CONTINUE
658 k = idamax( n, work, 1 )
659 CALL dlartg( vr( k, i ), vr( k, i+1 ), cs, sn, r )
660 CALL drot( n, vr( 1, i ), 1, vr( 1, i+1 ), 1, cs, sn )
661 vr( k, i+1 ) = zero
662 END IF
663 40 CONTINUE
664 END IF
665*
666* Undo scaling if necessary
667*
668 50 CONTINUE
669 IF( scalea ) THEN
670 CALL dlascl( 'G', 0, 0, cscale, anrm, n-info, 1, wr( info+1 ),
671 $ max( n-info, 1 ), ierr )
672 CALL dlascl( 'G', 0, 0, cscale, anrm, n-info, 1, wi( info+1 ),
673 $ max( n-info, 1 ), ierr )
674 IF( info.EQ.0 ) THEN
675 IF( ( wntsnv .OR. wntsnb ) .AND. icond.EQ.0 )
676 $ CALL dlascl( 'G', 0, 0, cscale, anrm, n, 1, rcondv, n,
677 $ ierr )
678 ELSE
679 CALL dlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wr, n,
680 $ ierr )
681 CALL dlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
682 $ ierr )
683 END IF
684 END IF
685*
686 work( 1 ) = maxwrk
687 RETURN
688*
689* End of DGEEVX
690*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
DGEBAK
Definition dgebak.f:130
subroutine dgebal(job, n, a, lda, ilo, ihi, scale, info)
DGEBAL
Definition dgebal.f:163
subroutine dgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
DGEHRD
Definition dgehrd.f:167
subroutine dhseqr(job, compz, n, ilo, ihi, h, ldh, wr, wi, z, ldz, work, lwork, info)
DHSEQR
Definition dhseqr.f:316
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
double precision function dlapy2(x, y)
DLAPY2 returns sqrt(x2+y2).
Definition dlapy2.f:63
subroutine dlartg(f, g, c, s, r)
DLARTG generates a plane rotation with real cosine and real sine.
Definition dlartg.f90:111
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function dnrm2(n, x, incx)
DNRM2
Definition dnrm2.f90:89
subroutine drot(n, dx, incx, dy, incy, c, s)
DROT
Definition drot.f:92
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dtrevc3(side, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, mm, m, work, lwork, info)
DTREVC3
Definition dtrevc3.f:237
subroutine dtrsna(job, howmny, select, n, t, ldt, vl, ldvl, vr, ldvr, s, sep, mm, m, work, ldwork, iwork, info)
DTRSNA
Definition dtrsna.f:265
subroutine dorghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
DORGHR
Definition dorghr.f:126
Here is the call graph for this function:
Here is the caller graph for this function: