LAPACK 3.12.1
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 301 of file dgeevx.f.

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