378 SUBROUTINE zdrges( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH,
379 $ NOUNIT, A, LDA, B, S, T, Q, LDQ, Z, ALPHA,
380 $ BETA, WORK, LWORK, RWORK, RESULT, BWORK, INFO )
387 INTEGER INFO, LDA, LDQ, LWORK, NOUNIT, NSIZES, NTYPES
388 DOUBLE PRECISION THRESH
391 LOGICAL BWORK( * ), DOTYPE( * )
392 INTEGER ISEED( 4 ), NN( * )
393 DOUBLE PRECISION RESULT( 13 ), RWORK( * )
394 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDA, * ),
395 $ beta( * ), q( ldq, * ), s( lda, * ),
396 $ t( lda, * ), work( * ), z( ldq, * )
402 DOUBLE PRECISION ZERO, ONE
403 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
404 COMPLEX*16 CZERO, CONE
405 parameter( czero = ( 0.0d+0, 0.0d+0 ),
406 $ cone = ( 1.0d+0, 0.0d+0 ) )
408 PARAMETER ( MAXTYP = 26 )
411 LOGICAL BADNN, ILABAD
413 INTEGER I, IADD, IINFO, IN, ISORT, J, JC, JR, JSIZE,
414 $ jtype, knteig, maxwrk, minwrk, mtypes, n, n1,
415 $ nb, nerrs, nmats, nmax, ntest, ntestt, rsub,
417 DOUBLE PRECISION SAFMAX, SAFMIN, TEMP1, TEMP2, ULP, ULPINV
421 LOGICAL LASIGN( MAXTYP ), LBSIGN( MAXTYP )
422 INTEGER IOLDSD( 4 ), KADD( 6 ), KAMAGN( MAXTYP ),
423 $ katype( maxtyp ), kazero( maxtyp ),
424 $ kbmagn( maxtyp ), kbtype( maxtyp ),
425 $ kbzero( maxtyp ), kclass( maxtyp ),
426 $ ktrian( maxtyp ), kz1( 6 ), kz2( 6 )
427 DOUBLE PRECISION RMAGN( 0: 3 )
432 DOUBLE PRECISION DLAMCH
434 EXTERNAL zlctes, ilaenv, dlamch, zlarnd
441 INTRINSIC abs, dble, dconjg, dimag, max, min, sign
444 DOUBLE PRECISION ABS1
447 abs1( x ) = abs( dble( x ) ) + abs( dimag( x ) )
450 DATA kclass / 15*1, 10*2, 1*3 /
451 DATA kz1 / 0, 1, 2, 1, 3, 3 /
452 DATA kz2 / 0, 0, 1, 2, 1, 1 /
453 DATA kadd / 0, 0, 0, 0, 3, 2 /
454 DATA katype / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4,
455 $ 4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0 /
456 DATA kbtype / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4,
457 $ 1, 1, -4, 2, -4, 8*8, 0 /
458 DATA kazero / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3,
460 DATA kbzero / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4,
462 DATA kamagn / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3,
464 DATA kbmagn / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3,
466 DATA ktrian / 16*0, 10*1 /
467 DATA lasign / 6*.false., .true., .false., 2*.true.,
468 $ 2*.false., 3*.true., .false., .true.,
469 $ 3*.false., 5*.true., .false. /
470 DATA lbsign / 7*.false., .true., 2*.false.,
471 $ 2*.true., 2*.false., .true., .false., .true.,
483 nmax = max( nmax, nn( j ) )
488 IF( nsizes.LT.0 )
THEN
490 ELSE IF( badnn )
THEN
492 ELSE IF( ntypes.LT.0 )
THEN
494 ELSE IF( thresh.LT.zero )
THEN
496 ELSE IF( lda.LE.1 .OR. lda.LT.nmax )
THEN
498 ELSE IF( ldq.LE.1 .OR. ldq.LT.nmax )
THEN
510 IF( info.EQ.0 .AND. lwork.GE.1 )
THEN
512 nb = max( 1, ilaenv( 1,
'ZGEQRF',
' ', nmax, nmax, -1, -1 ),
513 $ ilaenv( 1,
'ZUNMQR',
'LC', nmax, nmax, nmax, -1 ),
514 $ ilaenv( 1,
'ZUNGQR',
' ', nmax, nmax, nmax, -1 ) )
515 maxwrk = max( nmax+nmax*nb, 3*nmax*nmax )
519 IF( lwork.LT.minwrk )
523 CALL xerbla(
'ZDRGES', -info )
529 IF( nsizes.EQ.0 .OR. ntypes.EQ.0 )
532 ulp = dlamch(
'Precision' )
533 safmin = dlamch(
'Safe minimum' )
534 safmin = safmin / ulp
535 safmax = one / safmin
536 CALL dlabad( safmin, safmax )
550 DO 190 jsize = 1, nsizes
553 rmagn( 2 ) = safmax*ulp / dble( n1 )
554 rmagn( 3 ) = safmin*ulpinv*dble( n1 )
556 IF( nsizes.NE.1 )
THEN
557 mtypes = min( maxtyp, ntypes )
559 mtypes = min( maxtyp+1, ntypes )
564 DO 180 jtype = 1, mtypes
565 IF( .NOT.dotype( jtype ) )
573 ioldsd( j ) = iseed( j )
603 IF( mtypes.GT.maxtyp )
606 IF( kclass( jtype ).LT.3 )
THEN
610 IF( abs( katype( jtype ) ).EQ.3 )
THEN
611 in = 2*( ( n-1 ) / 2 ) + 1
613 $
CALL zlaset(
'Full', n, n, czero, czero, a, lda )
617 CALL zlatm4( katype( jtype ), in, kz1( kazero( jtype ) ),
618 $ kz2( kazero( jtype ) ), lasign( jtype ),
619 $ rmagn( kamagn( jtype ) ), ulp,
620 $ rmagn( ktrian( jtype )*kamagn( jtype ) ), 2,
622 iadd = kadd( kazero( jtype ) )
623 IF( iadd.GT.0 .AND. iadd.LE.n )
624 $ a( iadd, iadd ) = rmagn( kamagn( jtype ) )
628 IF( abs( kbtype( jtype ) ).EQ.3 )
THEN
629 in = 2*( ( n-1 ) / 2 ) + 1
631 $
CALL zlaset(
'Full', n, n, czero, czero, b, lda )
635 CALL zlatm4( kbtype( jtype ), in, kz1( kbzero( jtype ) ),
636 $ kz2( kbzero( jtype ) ), lbsign( jtype ),
637 $ rmagn( kbmagn( jtype ) ), one,
638 $ rmagn( ktrian( jtype )*kbmagn( jtype ) ), 2,
640 iadd = kadd( kbzero( jtype ) )
641 IF( iadd.NE.0 .AND. iadd.LE.n )
642 $ b( iadd, iadd ) = rmagn( kbmagn( jtype ) )
644 IF( kclass( jtype ).EQ.2 .AND. n.GT.0 )
THEN
653 q( jr, jc ) = zlarnd( 3, iseed )
654 z( jr, jc ) = zlarnd( 3, iseed )
656 CALL zlarfg( n+1-jc, q( jc, jc ), q( jc+1, jc ), 1,
658 work( 2*n+jc ) = sign( one, dble( q( jc, jc ) ) )
660 CALL zlarfg( n+1-jc, z( jc, jc ), z( jc+1, jc ), 1,
662 work( 3*n+jc ) = sign( one, dble( z( jc, jc ) ) )
665 ctemp = zlarnd( 3, iseed )
668 work( 3*n ) = ctemp / abs( ctemp )
669 ctemp = zlarnd( 3, iseed )
672 work( 4*n ) = ctemp / abs( ctemp )
678 a( jr, jc ) = work( 2*n+jr )*
679 $ dconjg( work( 3*n+jc ) )*
681 b( jr, jc ) = work( 2*n+jr )*
682 $ dconjg( work( 3*n+jc ) )*
686 CALL zunm2r(
'L',
'N', n, n, n-1, q, ldq, work, a,
687 $ lda, work( 2*n+1 ), iinfo )
690 CALL zunm2r(
'R',
'C', n, n, n-1, z, ldq, work( n+1 ),
691 $ a, lda, work( 2*n+1 ), iinfo )
694 CALL zunm2r(
'L',
'N', n, n, n-1, q, ldq, work, b,
695 $ lda, work( 2*n+1 ), iinfo )
698 CALL zunm2r(
'R',
'C', n, n, n-1, z, ldq, work( n+1 ),
699 $ b, lda, work( 2*n+1 ), iinfo )
709 a( jr, jc ) = rmagn( kamagn( jtype ) )*
711 b( jr, jc ) = rmagn( kbmagn( jtype ) )*
719 IF( iinfo.NE.0 )
THEN
720 WRITE( nounit, fmt = 9999 )
'Generator', iinfo, n, jtype,
735 IF( isort.EQ.0 )
THEN
745 CALL zlacpy(
'Full', n, n, a, lda, s, lda )
746 CALL zlacpy(
'Full', n, n, b, lda, t, lda )
747 ntest = 1 + rsub + isort
748 result( 1+rsub+isort ) = ulpinv
749 CALL zgges(
'V',
'V', sort, zlctes, n, s, lda, t, lda,
750 $ sdim, alpha, beta, q, ldq, z, ldq, work,
751 $ lwork, rwork, bwork, iinfo )
752 IF( iinfo.NE.0 .AND. iinfo.NE.n+2 )
THEN
753 result( 1+rsub+isort ) = ulpinv
754 WRITE( nounit, fmt = 9999 )
'ZGGES', iinfo, n, jtype,
764 IF( isort.EQ.0 )
THEN
765 CALL zget51( 1, n, a, lda, s, lda, q, ldq, z, ldq,
766 $ work, rwork, result( 1 ) )
767 CALL zget51( 1, n, b, lda, t, lda, q, ldq, z, ldq,
768 $ work, rwork, result( 2 ) )
770 CALL zget54( n, a, lda, b, lda, s, lda, t, lda, q,
771 $ ldq, z, ldq, work, result( 2+rsub ) )
774 CALL zget51( 3, n, b, lda, t, lda, q, ldq, q, ldq, work,
775 $ rwork, result( 3+rsub ) )
776 CALL zget51( 3, n, b, lda, t, lda, z, ldq, z, ldq, work,
777 $ rwork, result( 4+rsub ) )
788 temp2 = ( abs1( alpha( j )-s( j, j ) ) /
789 $ max( safmin, abs1( alpha( j ) ), abs1( s( j,
790 $ j ) ) )+abs1( beta( j )-t( j, j ) ) /
791 $ max( safmin, abs1( beta( j ) ), abs1( t( j,
795 IF( s( j+1, j ).NE.zero )
THEN
797 result( 5+rsub ) = ulpinv
801 IF( s( j, j-1 ).NE.zero )
THEN
803 result( 5+rsub ) = ulpinv
806 temp1 = max( temp1, temp2 )
808 WRITE( nounit, fmt = 9998 )j, n, jtype, ioldsd
811 result( 6+rsub ) = temp1
813 IF( isort.GE.1 )
THEN
821 IF( zlctes( alpha( i ), beta( i ) ) )
822 $ knteig = knteig + 1
825 $ result( 13 ) = ulpinv
834 ntestt = ntestt + ntest
839 IF( result( jr ).GE.thresh )
THEN
844 IF( nerrs.EQ.0 )
THEN
845 WRITE( nounit, fmt = 9997 )
'ZGS'
849 WRITE( nounit, fmt = 9996 )
850 WRITE( nounit, fmt = 9995 )
851 WRITE( nounit, fmt = 9994 )
'Unitary'
855 WRITE( nounit, fmt = 9993 )
'unitary',
'''',
856 $
'transpose', (
'''', j = 1, 8 )
860 IF( result( jr ).LT.10000.0d0 )
THEN
861 WRITE( nounit, fmt = 9992 )n, jtype, ioldsd, jr,
864 WRITE( nounit, fmt = 9991 )n, jtype, ioldsd, jr,
875 CALL alasvm(
'ZGS', nounit, nerrs, ntestt, 0 )
881 9999
FORMAT(
' ZDRGES: ', a,
' returned INFO=', i6,
'.', / 9x,
'N=',
882 $ i6,
', JTYPE=', i6,
', ISEED=(', 4( i4,
',' ), i5,
')' )
884 9998
FORMAT(
' ZDRGES: S not in Schur form at eigenvalue ', i6,
'.',
885 $ / 9x,
'N=', i6,
', JTYPE=', i6,
', ISEED=(', 3( i5,
',' ),
888 9997
FORMAT( / 1x, a3,
' -- Complex Generalized Schur from problem ',
891 9996
FORMAT(
' Matrix types (see ZDRGES for details): ' )
893 9995
FORMAT(
' Special Matrices:', 23x,
894 $
'(J''=transposed Jordan block)',
895 $ /
' 1=(0,0) 2=(I,0) 3=(0,I) 4=(I,I) 5=(J'',J'') ',
896 $
'6=(diag(J'',I), diag(I,J''))', /
' Diagonal Matrices: ( ',
897 $
'D=diag(0,1,2,...) )', /
' 7=(D,I) 9=(large*D, small*I',
898 $
') 11=(large*I, small*D) 13=(large*D, large*I)', /
899 $
' 8=(I,D) 10=(small*D, large*I) 12=(small*I, large*D) ',
900 $
' 14=(small*D, small*I)', /
' 15=(D, reversed D)' )
901 9994
FORMAT(
' Matrices Rotated by Random ', a,
' Matrices U, V:',
902 $ /
' 16=Transposed Jordan Blocks 19=geometric ',
903 $
'alpha, beta=0,1', /
' 17=arithm. alpha&beta ',
904 $
' 20=arithmetic alpha, beta=0,1', /
' 18=clustered ',
905 $
'alpha, beta=0,1 21=random alpha, beta=0,1',
906 $ /
' Large & Small Matrices:', /
' 22=(large, small) ',
907 $
'23=(small,large) 24=(small,small) 25=(large,large)',
908 $ /
' 26=random O(1) matrices.' )
910 9993
FORMAT( /
' Tests performed: (S is Schur, T is triangular, ',
911 $
'Q and Z are ', a,
',', / 19x,
912 $
'l and r are the appropriate left and right', / 19x,
913 $
'eigenvectors, resp., a is alpha, b is beta, and', / 19x, a,
914 $
' means ', a,
'.)', /
' Without ordering: ',
915 $ /
' 1 = | A - Q S Z', a,
916 $
' | / ( |A| n ulp ) 2 = | B - Q T Z', a,
917 $
' | / ( |B| n ulp )', /
' 3 = | I - QQ', a,
918 $
' | / ( n ulp ) 4 = | I - ZZ', a,
919 $
' | / ( n ulp )', /
' 5 = A is in Schur form S',
920 $ /
' 6 = difference between (alpha,beta)',
921 $
' and diagonals of (S,T)', /
' With ordering: ',
922 $ /
' 7 = | (A,B) - Q (S,T) Z', a,
' | / ( |(A,B)| n ulp )',
923 $ /
' 8 = | I - QQ', a,
924 $
' | / ( n ulp ) 9 = | I - ZZ', a,
925 $
' | / ( n ulp )', /
' 10 = A is in Schur form S',
926 $ /
' 11 = difference between (alpha,beta) and diagonals',
927 $
' of (S,T)', /
' 12 = SDIM is the correct number of ',
928 $
'selected eigenvalues', / )
929 9992
FORMAT(
' Matrix order=', i5,
', type=', i2,
', seed=',
930 $ 4( i4,
',' ),
' result ', i2,
' is', 0p, f8.2 )
931 9991
FORMAT(
' Matrix order=', i5,
', type=', i2,
', seed=',
932 $ 4( i4,
',' ),
' result ', i2,
' is', 1p, d10.3 )
subroutine dlabad(SMALL, LARGE)
DLABAD
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
subroutine zlatm4(ITYPE, N, NZ1, NZ2, RSIGN, AMAGN, RCOND, TRIANG, IDIST, ISEED, A, LDA)
ZLATM4
subroutine zget51(ITYPE, N, A, LDA, B, LDB, U, LDU, V, LDV, WORK, RWORK, RESULT)
ZGET51
subroutine zdrges(NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH, NOUNIT, A, LDA, B, S, T, Q, LDQ, Z, ALPHA, BETA, WORK, LWORK, RWORK, RESULT, BWORK, INFO)
ZDRGES
subroutine zget54(N, A, LDA, B, LDB, S, LDS, T, LDT, U, LDU, V, LDV, WORK, RESULT)
ZGET54
subroutine zgges(JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, WORK, LWORK, RWORK, BWORK, INFO)
ZGGES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE m...
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
subroutine zunm2r(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO)
ZUNM2R multiplies a general matrix by the unitary matrix from a QR factorization determined by cgeqrf...