LAPACK 3.3.0
|
00001 SUBROUTINE ZDRGES( NSIZES, NN, NTYPES, DOTYPE, ISEED, THRESH, 00002 $ NOUNIT, A, LDA, B, S, T, Q, LDQ, Z, ALPHA, 00003 $ BETA, WORK, LWORK, RWORK, RESULT, BWORK, INFO ) 00004 * 00005 * -- LAPACK test routine (version 3.1.1) -- 00006 * Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. 00007 * February 2007 00008 * 00009 * .. Scalar Arguments .. 00010 INTEGER INFO, LDA, LDQ, LWORK, NOUNIT, NSIZES, NTYPES 00011 DOUBLE PRECISION THRESH 00012 * .. 00013 * .. Array Arguments .. 00014 LOGICAL BWORK( * ), DOTYPE( * ) 00015 INTEGER ISEED( 4 ), NN( * ) 00016 DOUBLE PRECISION RESULT( 13 ), RWORK( * ) 00017 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDA, * ), 00018 $ BETA( * ), Q( LDQ, * ), S( LDA, * ), 00019 $ T( LDA, * ), WORK( * ), Z( LDQ, * ) 00020 * .. 00021 * 00022 * Purpose 00023 * ======= 00024 * 00025 * ZDRGES checks the nonsymmetric generalized eigenvalue (Schur form) 00026 * problem driver ZGGES. 00027 * 00028 * ZGGES factors A and B as Q*S*Z' and Q*T*Z' , where ' means conjugate 00029 * transpose, S and T are upper triangular (i.e., in generalized Schur 00030 * form), and Q and Z are unitary. It also computes the generalized 00031 * eigenvalues (alpha(j),beta(j)), j=1,...,n. Thus, 00032 * w(j) = alpha(j)/beta(j) is a root of the characteristic equation 00033 * 00034 * det( A - w(j) B ) = 0 00035 * 00036 * Optionally it also reorder the eigenvalues so that a selected 00037 * cluster of eigenvalues appears in the leading diagonal block of the 00038 * Schur forms. 00039 * 00040 * When ZDRGES is called, a number of matrix "sizes" ("N's") and a 00041 * number of matrix "TYPES" are specified. For each size ("N") 00042 * and each TYPE of matrix, a pair of matrices (A, B) will be generated 00043 * and used for testing. For each matrix pair, the following 13 tests 00044 * will be performed and compared with the threshhold THRESH except 00045 * the tests (5), (11) and (13). 00046 * 00047 * 00048 * (1) | A - Q S Z' | / ( |A| n ulp ) (no sorting of eigenvalues) 00049 * 00050 * 00051 * (2) | B - Q T Z' | / ( |B| n ulp ) (no sorting of eigenvalues) 00052 * 00053 * 00054 * (3) | I - QQ' | / ( n ulp ) (no sorting of eigenvalues) 00055 * 00056 * 00057 * (4) | I - ZZ' | / ( n ulp ) (no sorting of eigenvalues) 00058 * 00059 * (5) if A is in Schur form (i.e. triangular form) (no sorting of 00060 * eigenvalues) 00061 * 00062 * (6) if eigenvalues = diagonal elements of the Schur form (S, T), 00063 * i.e., test the maximum over j of D(j) where: 00064 * 00065 * |alpha(j) - S(j,j)| |beta(j) - T(j,j)| 00066 * D(j) = ------------------------ + ----------------------- 00067 * max(|alpha(j)|,|S(j,j)|) max(|beta(j)|,|T(j,j)|) 00068 * 00069 * (no sorting of eigenvalues) 00070 * 00071 * (7) | (A,B) - Q (S,T) Z' | / ( |(A,B)| n ulp ) 00072 * (with sorting of eigenvalues). 00073 * 00074 * (8) | I - QQ' | / ( n ulp ) (with sorting of eigenvalues). 00075 * 00076 * (9) | I - ZZ' | / ( n ulp ) (with sorting of eigenvalues). 00077 * 00078 * (10) if A is in Schur form (i.e. quasi-triangular form) 00079 * (with sorting of eigenvalues). 00080 * 00081 * (11) if eigenvalues = diagonal elements of the Schur form (S, T), 00082 * i.e. test the maximum over j of D(j) where: 00083 * 00084 * |alpha(j) - S(j,j)| |beta(j) - T(j,j)| 00085 * D(j) = ------------------------ + ----------------------- 00086 * max(|alpha(j)|,|S(j,j)|) max(|beta(j)|,|T(j,j)|) 00087 * 00088 * (with sorting of eigenvalues). 00089 * 00090 * (12) if sorting worked and SDIM is the number of eigenvalues 00091 * which were CELECTed. 00092 * 00093 * Test Matrices 00094 * ============= 00095 * 00096 * The sizes of the test matrices are specified by an array 00097 * NN(1:NSIZES); the value of each element NN(j) specifies one size. 00098 * The "types" are specified by a logical array DOTYPE( 1:NTYPES ); if 00099 * DOTYPE(j) is .TRUE., then matrix type "j" will be generated. 00100 * Currently, the list of possible types is: 00101 * 00102 * (1) ( 0, 0 ) (a pair of zero matrices) 00103 * 00104 * (2) ( I, 0 ) (an identity and a zero matrix) 00105 * 00106 * (3) ( 0, I ) (an identity and a zero matrix) 00107 * 00108 * (4) ( I, I ) (a pair of identity matrices) 00109 * 00110 * t t 00111 * (5) ( J , J ) (a pair of transposed Jordan blocks) 00112 * 00113 * t ( I 0 ) 00114 * (6) ( X, Y ) where X = ( J 0 ) and Y = ( t ) 00115 * ( 0 I ) ( 0 J ) 00116 * and I is a k x k identity and J a (k+1)x(k+1) 00117 * Jordan block; k=(N-1)/2 00118 * 00119 * (7) ( D, I ) where D is diag( 0, 1,..., N-1 ) (a diagonal 00120 * matrix with those diagonal entries.) 00121 * (8) ( I, D ) 00122 * 00123 * (9) ( big*D, small*I ) where "big" is near overflow and small=1/big 00124 * 00125 * (10) ( small*D, big*I ) 00126 * 00127 * (11) ( big*I, small*D ) 00128 * 00129 * (12) ( small*I, big*D ) 00130 * 00131 * (13) ( big*D, big*I ) 00132 * 00133 * (14) ( small*D, small*I ) 00134 * 00135 * (15) ( D1, D2 ) where D1 is diag( 0, 0, 1, ..., N-3, 0 ) and 00136 * D2 is diag( 0, N-3, N-4,..., 1, 0, 0 ) 00137 * t t 00138 * (16) Q ( J , J ) Z where Q and Z are random orthogonal matrices. 00139 * 00140 * (17) Q ( T1, T2 ) Z where T1 and T2 are upper triangular matrices 00141 * with random O(1) entries above the diagonal 00142 * and diagonal entries diag(T1) = 00143 * ( 0, 0, 1, ..., N-3, 0 ) and diag(T2) = 00144 * ( 0, N-3, N-4,..., 1, 0, 0 ) 00145 * 00146 * (18) Q ( T1, T2 ) Z diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 ) 00147 * diag(T2) = ( 0, 1, 0, 1,..., 1, 0 ) 00148 * s = machine precision. 00149 * 00150 * (19) Q ( T1, T2 ) Z diag(T1)=( 0,0,1,1, 1-d, ..., 1-(N-5)*d=s, 0 ) 00151 * diag(T2) = ( 0, 1, 0, 1, ..., 1, 0 ) 00152 * 00153 * N-5 00154 * (20) Q ( T1, T2 ) Z diag(T1)=( 0, 0, 1, 1, a, ..., a =s, 0 ) 00155 * diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 ) 00156 * 00157 * (21) Q ( T1, T2 ) Z diag(T1)=( 0, 0, 1, r1, r2, ..., r(N-4), 0 ) 00158 * diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 ) 00159 * where r1,..., r(N-4) are random. 00160 * 00161 * (22) Q ( big*T1, small*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 ) 00162 * diag(T2) = ( 0, 1, ..., 1, 0, 0 ) 00163 * 00164 * (23) Q ( small*T1, big*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 ) 00165 * diag(T2) = ( 0, 1, ..., 1, 0, 0 ) 00166 * 00167 * (24) Q ( small*T1, small*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 ) 00168 * diag(T2) = ( 0, 1, ..., 1, 0, 0 ) 00169 * 00170 * (25) Q ( big*T1, big*T2 ) Z diag(T1) = ( 0, 0, 1, ..., N-3, 0 ) 00171 * diag(T2) = ( 0, 1, ..., 1, 0, 0 ) 00172 * 00173 * (26) Q ( T1, T2 ) Z where T1 and T2 are random upper-triangular 00174 * matrices. 00175 * 00176 * 00177 * Arguments 00178 * ========= 00179 * 00180 * NSIZES (input) INTEGER 00181 * The number of sizes of matrices to use. If it is zero, 00182 * DDRGES does nothing. NSIZES >= 0. 00183 * 00184 * NN (input) INTEGER array, dimension (NSIZES) 00185 * An array containing the sizes to be used for the matrices. 00186 * Zero values will be skipped. NN >= 0. 00187 * 00188 * NTYPES (input) INTEGER 00189 * The number of elements in DOTYPE. If it is zero, DDRGES 00190 * does nothing. It must be at least zero. If it is MAXTYP+1 00191 * and NSIZES is 1, then an additional type, MAXTYP+1 is 00192 * defined, which is to use whatever matrix is in A on input. 00193 * This is only useful if DOTYPE(1:MAXTYP) is .FALSE. and 00194 * DOTYPE(MAXTYP+1) is .TRUE. . 00195 * 00196 * DOTYPE (input) LOGICAL array, dimension (NTYPES) 00197 * If DOTYPE(j) is .TRUE., then for each size in NN a 00198 * matrix of that size and of type j will be generated. 00199 * If NTYPES is smaller than the maximum number of types 00200 * defined (PARAMETER MAXTYP), then types NTYPES+1 through 00201 * MAXTYP will not be generated. If NTYPES is larger 00202 * than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES) 00203 * will be ignored. 00204 * 00205 * ISEED (input/output) INTEGER array, dimension (4) 00206 * On entry ISEED specifies the seed of the random number 00207 * generator. The array elements should be between 0 and 4095; 00208 * if not they will be reduced mod 4096. Also, ISEED(4) must 00209 * be odd. The random number generator uses a linear 00210 * congruential sequence limited to small integers, and so 00211 * should produce machine independent random numbers. The 00212 * values of ISEED are changed on exit, and can be used in the 00213 * next call to DDRGES to continue the same random number 00214 * sequence. 00215 * 00216 * THRESH (input) DOUBLE PRECISION 00217 * A test will count as "failed" if the "error", computed as 00218 * described above, exceeds THRESH. Note that the error is 00219 * scaled to be O(1), so THRESH should be a reasonably small 00220 * multiple of 1, e.g., 10 or 100. In particular, it should 00221 * not depend on the precision (single vs. double) or the size 00222 * of the matrix. THRESH >= 0. 00223 * 00224 * NOUNIT (input) INTEGER 00225 * The FORTRAN unit number for printing out error messages 00226 * (e.g., if a routine returns IINFO not equal to 0.) 00227 * 00228 * A (input/workspace) COMPLEX*16 array, dimension(LDA, max(NN)) 00229 * Used to hold the original A matrix. Used as input only 00230 * if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and 00231 * DOTYPE(MAXTYP+1)=.TRUE. 00232 * 00233 * LDA (input) INTEGER 00234 * The leading dimension of A, B, S, and T. 00235 * It must be at least 1 and at least max( NN ). 00236 * 00237 * B (input/workspace) COMPLEX*16 array, dimension(LDA, max(NN)) 00238 * Used to hold the original B matrix. Used as input only 00239 * if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and 00240 * DOTYPE(MAXTYP+1)=.TRUE. 00241 * 00242 * S (workspace) COMPLEX*16 array, dimension (LDA, max(NN)) 00243 * The Schur form matrix computed from A by ZGGES. On exit, S 00244 * contains the Schur form matrix corresponding to the matrix 00245 * in A. 00246 * 00247 * T (workspace) COMPLEX*16 array, dimension (LDA, max(NN)) 00248 * The upper triangular matrix computed from B by ZGGES. 00249 * 00250 * Q (workspace) COMPLEX*16 array, dimension (LDQ, max(NN)) 00251 * The (left) orthogonal matrix computed by ZGGES. 00252 * 00253 * LDQ (input) INTEGER 00254 * The leading dimension of Q and Z. It must 00255 * be at least 1 and at least max( NN ). 00256 * 00257 * Z (workspace) COMPLEX*16 array, dimension( LDQ, max(NN) ) 00258 * The (right) orthogonal matrix computed by ZGGES. 00259 * 00260 * ALPHA (workspace) COMPLEX*16 array, dimension (max(NN)) 00261 * BETA (workspace) COMPLEX*16 array, dimension (max(NN)) 00262 * The generalized eigenvalues of (A,B) computed by ZGGES. 00263 * ALPHA(k) / BETA(k) is the k-th generalized eigenvalue of A 00264 * and B. 00265 * 00266 * WORK (workspace) COMPLEX*16 array, dimension (LWORK) 00267 * 00268 * LWORK (input) INTEGER 00269 * The dimension of the array WORK. LWORK >= 3*N*N. 00270 * 00271 * RWORK (workspace) DOUBLE PRECISION array, dimension ( 8*N ) 00272 * Real workspace. 00273 * 00274 * RESULT (output) DOUBLE PRECISION array, dimension (15) 00275 * The values computed by the tests described above. 00276 * The values are currently limited to 1/ulp, to avoid overflow. 00277 * 00278 * BWORK (workspace) LOGICAL array, dimension (N) 00279 * 00280 * INFO (output) INTEGER 00281 * = 0: successful exit 00282 * < 0: if INFO = -i, the i-th argument had an illegal value. 00283 * > 0: A routine returned an error code. INFO is the 00284 * absolute value of the INFO value returned. 00285 * 00286 * ===================================================================== 00287 * 00288 * .. Parameters .. 00289 DOUBLE PRECISION ZERO, ONE 00290 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 00291 COMPLEX*16 CZERO, CONE 00292 PARAMETER ( CZERO = ( 0.0D+0, 0.0D+0 ), 00293 $ CONE = ( 1.0D+0, 0.0D+0 ) ) 00294 INTEGER MAXTYP 00295 PARAMETER ( MAXTYP = 26 ) 00296 * .. 00297 * .. Local Scalars .. 00298 LOGICAL BADNN, ILABAD 00299 CHARACTER SORT 00300 INTEGER I, IADD, IINFO, IN, ISORT, J, JC, JR, JSIZE, 00301 $ JTYPE, KNTEIG, MAXWRK, MINWRK, MTYPES, N, N1, 00302 $ NB, NERRS, NMATS, NMAX, NTEST, NTESTT, RSUB, 00303 $ SDIM 00304 DOUBLE PRECISION SAFMAX, SAFMIN, TEMP1, TEMP2, ULP, ULPINV 00305 COMPLEX*16 CTEMP, X 00306 * .. 00307 * .. Local Arrays .. 00308 LOGICAL LASIGN( MAXTYP ), LBSIGN( MAXTYP ) 00309 INTEGER IOLDSD( 4 ), KADD( 6 ), KAMAGN( MAXTYP ), 00310 $ KATYPE( MAXTYP ), KAZERO( MAXTYP ), 00311 $ KBMAGN( MAXTYP ), KBTYPE( MAXTYP ), 00312 $ KBZERO( MAXTYP ), KCLASS( MAXTYP ), 00313 $ KTRIAN( MAXTYP ), KZ1( 6 ), KZ2( 6 ) 00314 DOUBLE PRECISION RMAGN( 0: 3 ) 00315 * .. 00316 * .. External Functions .. 00317 LOGICAL ZLCTES 00318 INTEGER ILAENV 00319 DOUBLE PRECISION DLAMCH 00320 COMPLEX*16 ZLARND 00321 EXTERNAL ZLCTES, ILAENV, DLAMCH, ZLARND 00322 * .. 00323 * .. External Subroutines .. 00324 EXTERNAL ALASVM, DLABAD, XERBLA, ZGET51, ZGET54, ZGGES, 00325 $ ZLACPY, ZLARFG, ZLASET, ZLATM4, ZUNM2R 00326 * .. 00327 * .. Intrinsic Functions .. 00328 INTRINSIC ABS, DBLE, DCONJG, DIMAG, MAX, MIN, SIGN 00329 * .. 00330 * .. Statement Functions .. 00331 DOUBLE PRECISION ABS1 00332 * .. 00333 * .. Statement Function definitions .. 00334 ABS1( X ) = ABS( DBLE( X ) ) + ABS( DIMAG( X ) ) 00335 * .. 00336 * .. Data statements .. 00337 DATA KCLASS / 15*1, 10*2, 1*3 / 00338 DATA KZ1 / 0, 1, 2, 1, 3, 3 / 00339 DATA KZ2 / 0, 0, 1, 2, 1, 1 / 00340 DATA KADD / 0, 0, 0, 0, 3, 2 / 00341 DATA KATYPE / 0, 1, 0, 1, 2, 3, 4, 1, 4, 4, 1, 1, 4, 00342 $ 4, 4, 2, 4, 5, 8, 7, 9, 4*4, 0 / 00343 DATA KBTYPE / 0, 0, 1, 1, 2, -3, 1, 4, 1, 1, 4, 4, 00344 $ 1, 1, -4, 2, -4, 8*8, 0 / 00345 DATA KAZERO / 6*1, 2, 1, 2*2, 2*1, 2*2, 3, 1, 3, 00346 $ 4*5, 4*3, 1 / 00347 DATA KBZERO / 6*1, 1, 2, 2*1, 2*2, 2*1, 4, 1, 4, 00348 $ 4*6, 4*4, 1 / 00349 DATA KAMAGN / 8*1, 2, 3, 2, 3, 2, 3, 7*1, 2, 3, 3, 00350 $ 2, 1 / 00351 DATA KBMAGN / 8*1, 3, 2, 3, 2, 2, 3, 7*1, 3, 2, 3, 00352 $ 2, 1 / 00353 DATA KTRIAN / 16*0, 10*1 / 00354 DATA LASIGN / 6*.FALSE., .TRUE., .FALSE., 2*.TRUE., 00355 $ 2*.FALSE., 3*.TRUE., .FALSE., .TRUE., 00356 $ 3*.FALSE., 5*.TRUE., .FALSE. / 00357 DATA LBSIGN / 7*.FALSE., .TRUE., 2*.FALSE., 00358 $ 2*.TRUE., 2*.FALSE., .TRUE., .FALSE., .TRUE., 00359 $ 9*.FALSE. / 00360 * .. 00361 * .. Executable Statements .. 00362 * 00363 * Check for errors 00364 * 00365 INFO = 0 00366 * 00367 BADNN = .FALSE. 00368 NMAX = 1 00369 DO 10 J = 1, NSIZES 00370 NMAX = MAX( NMAX, NN( J ) ) 00371 IF( NN( J ).LT.0 ) 00372 $ BADNN = .TRUE. 00373 10 CONTINUE 00374 * 00375 IF( NSIZES.LT.0 ) THEN 00376 INFO = -1 00377 ELSE IF( BADNN ) THEN 00378 INFO = -2 00379 ELSE IF( NTYPES.LT.0 ) THEN 00380 INFO = -3 00381 ELSE IF( THRESH.LT.ZERO ) THEN 00382 INFO = -6 00383 ELSE IF( LDA.LE.1 .OR. LDA.LT.NMAX ) THEN 00384 INFO = -9 00385 ELSE IF( LDQ.LE.1 .OR. LDQ.LT.NMAX ) THEN 00386 INFO = -14 00387 END IF 00388 * 00389 * Compute workspace 00390 * (Note: Comments in the code beginning "Workspace:" describe the 00391 * minimal amount of workspace needed at that point in the code, 00392 * as well as the preferred amount for good performance. 00393 * NB refers to the optimal block size for the immediately 00394 * following subroutine, as returned by ILAENV. 00395 * 00396 MINWRK = 1 00397 IF( INFO.EQ.0 .AND. LWORK.GE.1 ) THEN 00398 MINWRK = 3*NMAX*NMAX 00399 NB = MAX( 1, ILAENV( 1, 'ZGEQRF', ' ', NMAX, NMAX, -1, -1 ), 00400 $ ILAENV( 1, 'ZUNMQR', 'LC', NMAX, NMAX, NMAX, -1 ), 00401 $ ILAENV( 1, 'ZUNGQR', ' ', NMAX, NMAX, NMAX, -1 ) ) 00402 MAXWRK = MAX( NMAX+NMAX*NB, 3*NMAX*NMAX ) 00403 WORK( 1 ) = MAXWRK 00404 END IF 00405 * 00406 IF( LWORK.LT.MINWRK ) 00407 $ INFO = -19 00408 * 00409 IF( INFO.NE.0 ) THEN 00410 CALL XERBLA( 'ZDRGES', -INFO ) 00411 RETURN 00412 END IF 00413 * 00414 * Quick return if possible 00415 * 00416 IF( NSIZES.EQ.0 .OR. NTYPES.EQ.0 ) 00417 $ RETURN 00418 * 00419 ULP = DLAMCH( 'Precision' ) 00420 SAFMIN = DLAMCH( 'Safe minimum' ) 00421 SAFMIN = SAFMIN / ULP 00422 SAFMAX = ONE / SAFMIN 00423 CALL DLABAD( SAFMIN, SAFMAX ) 00424 ULPINV = ONE / ULP 00425 * 00426 * The values RMAGN(2:3) depend on N, see below. 00427 * 00428 RMAGN( 0 ) = ZERO 00429 RMAGN( 1 ) = ONE 00430 * 00431 * Loop over matrix sizes 00432 * 00433 NTESTT = 0 00434 NERRS = 0 00435 NMATS = 0 00436 * 00437 DO 190 JSIZE = 1, NSIZES 00438 N = NN( JSIZE ) 00439 N1 = MAX( 1, N ) 00440 RMAGN( 2 ) = SAFMAX*ULP / DBLE( N1 ) 00441 RMAGN( 3 ) = SAFMIN*ULPINV*DBLE( N1 ) 00442 * 00443 IF( NSIZES.NE.1 ) THEN 00444 MTYPES = MIN( MAXTYP, NTYPES ) 00445 ELSE 00446 MTYPES = MIN( MAXTYP+1, NTYPES ) 00447 END IF 00448 * 00449 * Loop over matrix types 00450 * 00451 DO 180 JTYPE = 1, MTYPES 00452 IF( .NOT.DOTYPE( JTYPE ) ) 00453 $ GO TO 180 00454 NMATS = NMATS + 1 00455 NTEST = 0 00456 * 00457 * Save ISEED in case of an error. 00458 * 00459 DO 20 J = 1, 4 00460 IOLDSD( J ) = ISEED( J ) 00461 20 CONTINUE 00462 * 00463 * Initialize RESULT 00464 * 00465 DO 30 J = 1, 13 00466 RESULT( J ) = ZERO 00467 30 CONTINUE 00468 * 00469 * Generate test matrices A and B 00470 * 00471 * Description of control parameters: 00472 * 00473 * KZLASS: =1 means w/o rotation, =2 means w/ rotation, 00474 * =3 means random. 00475 * KATYPE: the "type" to be passed to ZLATM4 for computing A. 00476 * KAZERO: the pattern of zeros on the diagonal for A: 00477 * =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ), 00478 * =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ), 00479 * =6: ( 0, 1, 0, xxx, 0 ). (xxx means a string of 00480 * non-zero entries.) 00481 * KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1), 00482 * =2: large, =3: small. 00483 * LASIGN: .TRUE. if the diagonal elements of A are to be 00484 * multiplied by a random magnitude 1 number. 00485 * KBTYPE, KBZERO, KBMAGN, LBSIGN: the same, but for B. 00486 * KTRIAN: =0: don't fill in the upper triangle, =1: do. 00487 * KZ1, KZ2, KADD: used to implement KAZERO and KBZERO. 00488 * RMAGN: used to implement KAMAGN and KBMAGN. 00489 * 00490 IF( MTYPES.GT.MAXTYP ) 00491 $ GO TO 110 00492 IINFO = 0 00493 IF( KCLASS( JTYPE ).LT.3 ) THEN 00494 * 00495 * Generate A (w/o rotation) 00496 * 00497 IF( ABS( KATYPE( JTYPE ) ).EQ.3 ) THEN 00498 IN = 2*( ( N-1 ) / 2 ) + 1 00499 IF( IN.NE.N ) 00500 $ CALL ZLASET( 'Full', N, N, CZERO, CZERO, A, LDA ) 00501 ELSE 00502 IN = N 00503 END IF 00504 CALL ZLATM4( KATYPE( JTYPE ), IN, KZ1( KAZERO( JTYPE ) ), 00505 $ KZ2( KAZERO( JTYPE ) ), LASIGN( JTYPE ), 00506 $ RMAGN( KAMAGN( JTYPE ) ), ULP, 00507 $ RMAGN( KTRIAN( JTYPE )*KAMAGN( JTYPE ) ), 2, 00508 $ ISEED, A, LDA ) 00509 IADD = KADD( KAZERO( JTYPE ) ) 00510 IF( IADD.GT.0 .AND. IADD.LE.N ) 00511 $ A( IADD, IADD ) = RMAGN( KAMAGN( JTYPE ) ) 00512 * 00513 * Generate B (w/o rotation) 00514 * 00515 IF( ABS( KBTYPE( JTYPE ) ).EQ.3 ) THEN 00516 IN = 2*( ( N-1 ) / 2 ) + 1 00517 IF( IN.NE.N ) 00518 $ CALL ZLASET( 'Full', N, N, CZERO, CZERO, B, LDA ) 00519 ELSE 00520 IN = N 00521 END IF 00522 CALL ZLATM4( KBTYPE( JTYPE ), IN, KZ1( KBZERO( JTYPE ) ), 00523 $ KZ2( KBZERO( JTYPE ) ), LBSIGN( JTYPE ), 00524 $ RMAGN( KBMAGN( JTYPE ) ), ONE, 00525 $ RMAGN( KTRIAN( JTYPE )*KBMAGN( JTYPE ) ), 2, 00526 $ ISEED, B, LDA ) 00527 IADD = KADD( KBZERO( JTYPE ) ) 00528 IF( IADD.NE.0 .AND. IADD.LE.N ) 00529 $ B( IADD, IADD ) = RMAGN( KBMAGN( JTYPE ) ) 00530 * 00531 IF( KCLASS( JTYPE ).EQ.2 .AND. N.GT.0 ) THEN 00532 * 00533 * Include rotations 00534 * 00535 * Generate Q, Z as Householder transformations times 00536 * a diagonal matrix. 00537 * 00538 DO 50 JC = 1, N - 1 00539 DO 40 JR = JC, N 00540 Q( JR, JC ) = ZLARND( 3, ISEED ) 00541 Z( JR, JC ) = ZLARND( 3, ISEED ) 00542 40 CONTINUE 00543 CALL ZLARFG( N+1-JC, Q( JC, JC ), Q( JC+1, JC ), 1, 00544 $ WORK( JC ) ) 00545 WORK( 2*N+JC ) = SIGN( ONE, DBLE( Q( JC, JC ) ) ) 00546 Q( JC, JC ) = CONE 00547 CALL ZLARFG( N+1-JC, Z( JC, JC ), Z( JC+1, JC ), 1, 00548 $ WORK( N+JC ) ) 00549 WORK( 3*N+JC ) = SIGN( ONE, DBLE( Z( JC, JC ) ) ) 00550 Z( JC, JC ) = CONE 00551 50 CONTINUE 00552 CTEMP = ZLARND( 3, ISEED ) 00553 Q( N, N ) = CONE 00554 WORK( N ) = CZERO 00555 WORK( 3*N ) = CTEMP / ABS( CTEMP ) 00556 CTEMP = ZLARND( 3, ISEED ) 00557 Z( N, N ) = CONE 00558 WORK( 2*N ) = CZERO 00559 WORK( 4*N ) = CTEMP / ABS( CTEMP ) 00560 * 00561 * Apply the diagonal matrices 00562 * 00563 DO 70 JC = 1, N 00564 DO 60 JR = 1, N 00565 A( JR, JC ) = WORK( 2*N+JR )* 00566 $ DCONJG( WORK( 3*N+JC ) )* 00567 $ A( JR, JC ) 00568 B( JR, JC ) = WORK( 2*N+JR )* 00569 $ DCONJG( WORK( 3*N+JC ) )* 00570 $ B( JR, JC ) 00571 60 CONTINUE 00572 70 CONTINUE 00573 CALL ZUNM2R( 'L', 'N', N, N, N-1, Q, LDQ, WORK, A, 00574 $ LDA, WORK( 2*N+1 ), IINFO ) 00575 IF( IINFO.NE.0 ) 00576 $ GO TO 100 00577 CALL ZUNM2R( 'R', 'C', N, N, N-1, Z, LDQ, WORK( N+1 ), 00578 $ A, LDA, WORK( 2*N+1 ), IINFO ) 00579 IF( IINFO.NE.0 ) 00580 $ GO TO 100 00581 CALL ZUNM2R( 'L', 'N', N, N, N-1, Q, LDQ, WORK, B, 00582 $ LDA, WORK( 2*N+1 ), IINFO ) 00583 IF( IINFO.NE.0 ) 00584 $ GO TO 100 00585 CALL ZUNM2R( 'R', 'C', N, N, N-1, Z, LDQ, WORK( N+1 ), 00586 $ B, LDA, WORK( 2*N+1 ), IINFO ) 00587 IF( IINFO.NE.0 ) 00588 $ GO TO 100 00589 END IF 00590 ELSE 00591 * 00592 * Random matrices 00593 * 00594 DO 90 JC = 1, N 00595 DO 80 JR = 1, N 00596 A( JR, JC ) = RMAGN( KAMAGN( JTYPE ) )* 00597 $ ZLARND( 4, ISEED ) 00598 B( JR, JC ) = RMAGN( KBMAGN( JTYPE ) )* 00599 $ ZLARND( 4, ISEED ) 00600 80 CONTINUE 00601 90 CONTINUE 00602 END IF 00603 * 00604 100 CONTINUE 00605 * 00606 IF( IINFO.NE.0 ) THEN 00607 WRITE( NOUNIT, FMT = 9999 )'Generator', IINFO, N, JTYPE, 00608 $ IOLDSD 00609 INFO = ABS( IINFO ) 00610 RETURN 00611 END IF 00612 * 00613 110 CONTINUE 00614 * 00615 DO 120 I = 1, 13 00616 RESULT( I ) = -ONE 00617 120 CONTINUE 00618 * 00619 * Test with and without sorting of eigenvalues 00620 * 00621 DO 150 ISORT = 0, 1 00622 IF( ISORT.EQ.0 ) THEN 00623 SORT = 'N' 00624 RSUB = 0 00625 ELSE 00626 SORT = 'S' 00627 RSUB = 5 00628 END IF 00629 * 00630 * Call ZGGES to compute H, T, Q, Z, alpha, and beta. 00631 * 00632 CALL ZLACPY( 'Full', N, N, A, LDA, S, LDA ) 00633 CALL ZLACPY( 'Full', N, N, B, LDA, T, LDA ) 00634 NTEST = 1 + RSUB + ISORT 00635 RESULT( 1+RSUB+ISORT ) = ULPINV 00636 CALL ZGGES( 'V', 'V', SORT, ZLCTES, N, S, LDA, T, LDA, 00637 $ SDIM, ALPHA, BETA, Q, LDQ, Z, LDQ, WORK, 00638 $ LWORK, RWORK, BWORK, IINFO ) 00639 IF( IINFO.NE.0 .AND. IINFO.NE.N+2 ) THEN 00640 RESULT( 1+RSUB+ISORT ) = ULPINV 00641 WRITE( NOUNIT, FMT = 9999 )'ZGGES', IINFO, N, JTYPE, 00642 $ IOLDSD 00643 INFO = ABS( IINFO ) 00644 GO TO 160 00645 END IF 00646 * 00647 NTEST = 4 + RSUB 00648 * 00649 * Do tests 1--4 (or tests 7--9 when reordering ) 00650 * 00651 IF( ISORT.EQ.0 ) THEN 00652 CALL ZGET51( 1, N, A, LDA, S, LDA, Q, LDQ, Z, LDQ, 00653 $ WORK, RWORK, RESULT( 1 ) ) 00654 CALL ZGET51( 1, N, B, LDA, T, LDA, Q, LDQ, Z, LDQ, 00655 $ WORK, RWORK, RESULT( 2 ) ) 00656 ELSE 00657 CALL ZGET54( N, A, LDA, B, LDA, S, LDA, T, LDA, Q, 00658 $ LDQ, Z, LDQ, WORK, RESULT( 2+RSUB ) ) 00659 END IF 00660 * 00661 CALL ZGET51( 3, N, B, LDA, T, LDA, Q, LDQ, Q, LDQ, WORK, 00662 $ RWORK, RESULT( 3+RSUB ) ) 00663 CALL ZGET51( 3, N, B, LDA, T, LDA, Z, LDQ, Z, LDQ, WORK, 00664 $ RWORK, RESULT( 4+RSUB ) ) 00665 * 00666 * Do test 5 and 6 (or Tests 10 and 11 when reordering): 00667 * check Schur form of A and compare eigenvalues with 00668 * diagonals. 00669 * 00670 NTEST = 6 + RSUB 00671 TEMP1 = ZERO 00672 * 00673 DO 130 J = 1, N 00674 ILABAD = .FALSE. 00675 TEMP2 = ( ABS1( ALPHA( J )-S( J, J ) ) / 00676 $ MAX( SAFMIN, ABS1( ALPHA( J ) ), ABS1( S( J, 00677 $ J ) ) )+ABS1( BETA( J )-T( J, J ) ) / 00678 $ MAX( SAFMIN, ABS1( BETA( J ) ), ABS1( T( J, 00679 $ J ) ) ) ) / ULP 00680 * 00681 IF( J.LT.N ) THEN 00682 IF( S( J+1, J ).NE.ZERO ) THEN 00683 ILABAD = .TRUE. 00684 RESULT( 5+RSUB ) = ULPINV 00685 END IF 00686 END IF 00687 IF( J.GT.1 ) THEN 00688 IF( S( J, J-1 ).NE.ZERO ) THEN 00689 ILABAD = .TRUE. 00690 RESULT( 5+RSUB ) = ULPINV 00691 END IF 00692 END IF 00693 TEMP1 = MAX( TEMP1, TEMP2 ) 00694 IF( ILABAD ) THEN 00695 WRITE( NOUNIT, FMT = 9998 )J, N, JTYPE, IOLDSD 00696 END IF 00697 130 CONTINUE 00698 RESULT( 6+RSUB ) = TEMP1 00699 * 00700 IF( ISORT.GE.1 ) THEN 00701 * 00702 * Do test 12 00703 * 00704 NTEST = 12 00705 RESULT( 12 ) = ZERO 00706 KNTEIG = 0 00707 DO 140 I = 1, N 00708 IF( ZLCTES( ALPHA( I ), BETA( I ) ) ) 00709 $ KNTEIG = KNTEIG + 1 00710 140 CONTINUE 00711 IF( SDIM.NE.KNTEIG ) 00712 $ RESULT( 13 ) = ULPINV 00713 END IF 00714 * 00715 150 CONTINUE 00716 * 00717 * End of Loop -- Check for RESULT(j) > THRESH 00718 * 00719 160 CONTINUE 00720 * 00721 NTESTT = NTESTT + NTEST 00722 * 00723 * Print out tests which fail. 00724 * 00725 DO 170 JR = 1, NTEST 00726 IF( RESULT( JR ).GE.THRESH ) THEN 00727 * 00728 * If this is the first test to fail, 00729 * print a header to the data file. 00730 * 00731 IF( NERRS.EQ.0 ) THEN 00732 WRITE( NOUNIT, FMT = 9997 )'ZGS' 00733 * 00734 * Matrix types 00735 * 00736 WRITE( NOUNIT, FMT = 9996 ) 00737 WRITE( NOUNIT, FMT = 9995 ) 00738 WRITE( NOUNIT, FMT = 9994 )'Unitary' 00739 * 00740 * Tests performed 00741 * 00742 WRITE( NOUNIT, FMT = 9993 )'unitary', '''', 00743 $ 'transpose', ( '''', J = 1, 8 ) 00744 * 00745 END IF 00746 NERRS = NERRS + 1 00747 IF( RESULT( JR ).LT.10000.0D0 ) THEN 00748 WRITE( NOUNIT, FMT = 9992 )N, JTYPE, IOLDSD, JR, 00749 $ RESULT( JR ) 00750 ELSE 00751 WRITE( NOUNIT, FMT = 9991 )N, JTYPE, IOLDSD, JR, 00752 $ RESULT( JR ) 00753 END IF 00754 END IF 00755 170 CONTINUE 00756 * 00757 180 CONTINUE 00758 190 CONTINUE 00759 * 00760 * Summary 00761 * 00762 CALL ALASVM( 'ZGS', NOUNIT, NERRS, NTESTT, 0 ) 00763 * 00764 WORK( 1 ) = MAXWRK 00765 * 00766 RETURN 00767 * 00768 9999 FORMAT( ' ZDRGES: ', A, ' returned INFO=', I6, '.', / 9X, 'N=', 00769 $ I6, ', JTYPE=', I6, ', ISEED=(', 4( I4, ',' ), I5, ')' ) 00770 * 00771 9998 FORMAT( ' ZDRGES: S not in Schur form at eigenvalue ', I6, '.', 00772 $ / 9X, 'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), 00773 $ I5, ')' ) 00774 * 00775 9997 FORMAT( / 1X, A3, ' -- Complex Generalized Schur from problem ', 00776 $ 'driver' ) 00777 * 00778 9996 FORMAT( ' Matrix types (see ZDRGES for details): ' ) 00779 * 00780 9995 FORMAT( ' Special Matrices:', 23X, 00781 $ '(J''=transposed Jordan block)', 00782 $ / ' 1=(0,0) 2=(I,0) 3=(0,I) 4=(I,I) 5=(J'',J'') ', 00783 $ '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices: ( ', 00784 $ 'D=diag(0,1,2,...) )', / ' 7=(D,I) 9=(large*D, small*I', 00785 $ ') 11=(large*I, small*D) 13=(large*D, large*I)', / 00786 $ ' 8=(I,D) 10=(small*D, large*I) 12=(small*I, large*D) ', 00787 $ ' 14=(small*D, small*I)', / ' 15=(D, reversed D)' ) 00788 9994 FORMAT( ' Matrices Rotated by Random ', A, ' Matrices U, V:', 00789 $ / ' 16=Transposed Jordan Blocks 19=geometric ', 00790 $ 'alpha, beta=0,1', / ' 17=arithm. alpha&beta ', 00791 $ ' 20=arithmetic alpha, beta=0,1', / ' 18=clustered ', 00792 $ 'alpha, beta=0,1 21=random alpha, beta=0,1', 00793 $ / ' Large & Small Matrices:', / ' 22=(large, small) ', 00794 $ '23=(small,large) 24=(small,small) 25=(large,large)', 00795 $ / ' 26=random O(1) matrices.' ) 00796 * 00797 9993 FORMAT( / ' Tests performed: (S is Schur, T is triangular, ', 00798 $ 'Q and Z are ', A, ',', / 19X, 00799 $ 'l and r are the appropriate left and right', / 19X, 00800 $ 'eigenvectors, resp., a is alpha, b is beta, and', / 19X, A, 00801 $ ' means ', A, '.)', / ' Without ordering: ', 00802 $ / ' 1 = | A - Q S Z', A, 00803 $ ' | / ( |A| n ulp ) 2 = | B - Q T Z', A, 00804 $ ' | / ( |B| n ulp )', / ' 3 = | I - QQ', A, 00805 $ ' | / ( n ulp ) 4 = | I - ZZ', A, 00806 $ ' | / ( n ulp )', / ' 5 = A is in Schur form S', 00807 $ / ' 6 = difference between (alpha,beta)', 00808 $ ' and diagonals of (S,T)', / ' With ordering: ', 00809 $ / ' 7 = | (A,B) - Q (S,T) Z', A, ' | / ( |(A,B)| n ulp )', 00810 $ / ' 8 = | I - QQ', A, 00811 $ ' | / ( n ulp ) 9 = | I - ZZ', A, 00812 $ ' | / ( n ulp )', / ' 10 = A is in Schur form S', 00813 $ / ' 11 = difference between (alpha,beta) and diagonals', 00814 $ ' of (S,T)', / ' 12 = SDIM is the correct number of ', 00815 $ 'selected eigenvalues', / ) 00816 9992 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=', 00817 $ 4( I4, ',' ), ' result ', I2, ' is', 0P, F8.2 ) 00818 9991 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=', 00819 $ 4( I4, ',' ), ' result ', I2, ' is', 1P, D10.3 ) 00820 * 00821 * End of ZDRGES 00822 * 00823 END