LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

◆ zgeevx()

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

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

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

Purpose:
!> !> ZGEEVX 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*16 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*16 array, dimension (N) !> W contains the computed eigenvalues. !>
[out]VL
!> VL is COMPLEX*16 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*16 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 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 COMPLEX*16 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 DOUBLE PRECISION 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 283 of file zgeevx.f.

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