1*> \brief \b CGET24 2* 3* =========== DOCUMENTATION =========== 4* 5* Online html documentation available at 6* http://www.netlib.org/lapack/explore-html/ 7* 8* Definition: 9* =========== 10* 11* SUBROUTINE CGET24( COMP, JTYPE, THRESH, ISEED, NOUNIT, N, A, LDA, 12* H, HT, W, WT, WTMP, VS, LDVS, VS1, RCDEIN, 13* RCDVIN, NSLCT, ISLCT, ISRT, RESULT, WORK, 14* LWORK, RWORK, BWORK, INFO ) 15* 16* .. Scalar Arguments .. 17* LOGICAL COMP 18* INTEGER INFO, ISRT, JTYPE, LDA, LDVS, LWORK, N, NOUNIT, 19* $ NSLCT 20* REAL RCDEIN, RCDVIN, THRESH 21* .. 22* .. Array Arguments .. 23* LOGICAL BWORK( * ) 24* INTEGER ISEED( 4 ), ISLCT( * ) 25* REAL RESULT( 17 ), RWORK( * ) 26* COMPLEX A( LDA, * ), H( LDA, * ), HT( LDA, * ), 27* $ VS( LDVS, * ), VS1( LDVS, * ), W( * ), 28* $ WORK( * ), WT( * ), WTMP( * ) 29* .. 30* 31* 32*> \par Purpose: 33* ============= 34*> 35*> \verbatim 36*> 37*> CGET24 checks the nonsymmetric eigenvalue (Schur form) problem 38*> expert driver CGEESX. 39*> 40*> If COMP = .FALSE., the first 13 of the following tests will be 41*> be performed on the input matrix A, and also tests 14 and 15 42*> if LWORK is sufficiently large. 43*> If COMP = .TRUE., all 17 test will be performed. 44*> 45*> (1) 0 if T is in Schur form, 1/ulp otherwise 46*> (no sorting of eigenvalues) 47*> 48*> (2) | A - VS T VS' | / ( n |A| ulp ) 49*> 50*> Here VS is the matrix of Schur eigenvectors, and T is in Schur 51*> form (no sorting of eigenvalues). 52*> 53*> (3) | I - VS VS' | / ( n ulp ) (no sorting of eigenvalues). 54*> 55*> (4) 0 if W are eigenvalues of T 56*> 1/ulp otherwise 57*> (no sorting of eigenvalues) 58*> 59*> (5) 0 if T(with VS) = T(without VS), 60*> 1/ulp otherwise 61*> (no sorting of eigenvalues) 62*> 63*> (6) 0 if eigenvalues(with VS) = eigenvalues(without VS), 64*> 1/ulp otherwise 65*> (no sorting of eigenvalues) 66*> 67*> (7) 0 if T is in Schur form, 1/ulp otherwise 68*> (with sorting of eigenvalues) 69*> 70*> (8) | A - VS T VS' | / ( n |A| ulp ) 71*> 72*> Here VS is the matrix of Schur eigenvectors, and T is in Schur 73*> form (with sorting of eigenvalues). 74*> 75*> (9) | I - VS VS' | / ( n ulp ) (with sorting of eigenvalues). 76*> 77*> (10) 0 if W are eigenvalues of T 78*> 1/ulp otherwise 79*> If workspace sufficient, also compare W with and 80*> without reciprocal condition numbers 81*> (with sorting of eigenvalues) 82*> 83*> (11) 0 if T(with VS) = T(without VS), 84*> 1/ulp otherwise 85*> If workspace sufficient, also compare T with and without 86*> reciprocal condition numbers 87*> (with sorting of eigenvalues) 88*> 89*> (12) 0 if eigenvalues(with VS) = eigenvalues(without VS), 90*> 1/ulp otherwise 91*> If workspace sufficient, also compare VS with and without 92*> reciprocal condition numbers 93*> (with sorting of eigenvalues) 94*> 95*> (13) if sorting worked and SDIM is the number of 96*> eigenvalues which were SELECTed 97*> If workspace sufficient, also compare SDIM with and 98*> without reciprocal condition numbers 99*> 100*> (14) if RCONDE the same no matter if VS and/or RCONDV computed 101*> 102*> (15) if RCONDV the same no matter if VS and/or RCONDE computed 103*> 104*> (16) |RCONDE - RCDEIN| / cond(RCONDE) 105*> 106*> RCONDE is the reciprocal average eigenvalue condition number 107*> computed by CGEESX and RCDEIN (the precomputed true value) 108*> is supplied as input. cond(RCONDE) is the condition number 109*> of RCONDE, and takes errors in computing RCONDE into account, 110*> so that the resulting quantity should be O(ULP). cond(RCONDE) 111*> is essentially given by norm(A)/RCONDV. 112*> 113*> (17) |RCONDV - RCDVIN| / cond(RCONDV) 114*> 115*> RCONDV is the reciprocal right invariant subspace condition 116*> number computed by CGEESX and RCDVIN (the precomputed true 117*> value) is supplied as input. cond(RCONDV) is the condition 118*> number of RCONDV, and takes errors in computing RCONDV into 119*> account, so that the resulting quantity should be O(ULP). 120*> cond(RCONDV) is essentially given by norm(A)/RCONDE. 121*> \endverbatim 122* 123* Arguments: 124* ========== 125* 126*> \param[in] COMP 127*> \verbatim 128*> COMP is LOGICAL 129*> COMP describes which input tests to perform: 130*> = .FALSE. if the computed condition numbers are not to 131*> be tested against RCDVIN and RCDEIN 132*> = .TRUE. if they are to be compared 133*> \endverbatim 134*> 135*> \param[in] JTYPE 136*> \verbatim 137*> JTYPE is INTEGER 138*> Type of input matrix. Used to label output if error occurs. 139*> \endverbatim 140*> 141*> \param[in] ISEED 142*> \verbatim 143*> ISEED is INTEGER array, dimension (4) 144*> If COMP = .FALSE., the random number generator seed 145*> used to produce matrix. 146*> If COMP = .TRUE., ISEED(1) = the number of the example. 147*> Used to label output if error occurs. 148*> \endverbatim 149*> 150*> \param[in] THRESH 151*> \verbatim 152*> THRESH is REAL 153*> A test will count as "failed" if the "error", computed as 154*> described above, exceeds THRESH. Note that the error 155*> is scaled to be O(1), so THRESH should be a reasonably 156*> small multiple of 1, e.g., 10 or 100. In particular, 157*> it should not depend on the precision (single vs. double) 158*> or the size of the matrix. It must be at least zero. 159*> \endverbatim 160*> 161*> \param[in] NOUNIT 162*> \verbatim 163*> NOUNIT is INTEGER 164*> The FORTRAN unit number for printing out error messages 165*> (e.g., if a routine returns INFO not equal to 0.) 166*> \endverbatim 167*> 168*> \param[in] N 169*> \verbatim 170*> N is INTEGER 171*> The dimension of A. N must be at least 0. 172*> \endverbatim 173*> 174*> \param[in,out] A 175*> \verbatim 176*> A is COMPLEX array, dimension (LDA, N) 177*> Used to hold the matrix whose eigenvalues are to be 178*> computed. 179*> \endverbatim 180*> 181*> \param[in] LDA 182*> \verbatim 183*> LDA is INTEGER 184*> The leading dimension of A, and H. LDA must be at 185*> least 1 and at least N. 186*> \endverbatim 187*> 188*> \param[out] H 189*> \verbatim 190*> H is COMPLEX array, dimension (LDA, N) 191*> Another copy of the test matrix A, modified by CGEESX. 192*> \endverbatim 193*> 194*> \param[out] HT 195*> \verbatim 196*> HT is COMPLEX array, dimension (LDA, N) 197*> Yet another copy of the test matrix A, modified by CGEESX. 198*> \endverbatim 199*> 200*> \param[out] W 201*> \verbatim 202*> W is COMPLEX array, dimension (N) 203*> The computed eigenvalues of A. 204*> \endverbatim 205*> 206*> \param[out] WT 207*> \verbatim 208*> WT is COMPLEX array, dimension (N) 209*> Like W, this array contains the eigenvalues of A, 210*> but those computed when CGEESX only computes a partial 211*> eigendecomposition, i.e. not Schur vectors 212*> \endverbatim 213*> 214*> \param[out] WTMP 215*> \verbatim 216*> WTMP is COMPLEX array, dimension (N) 217*> Like W, this array contains the eigenvalues of A, 218*> but sorted by increasing real or imaginary part. 219*> \endverbatim 220*> 221*> \param[out] VS 222*> \verbatim 223*> VS is COMPLEX array, dimension (LDVS, N) 224*> VS holds the computed Schur vectors. 225*> \endverbatim 226*> 227*> \param[in] LDVS 228*> \verbatim 229*> LDVS is INTEGER 230*> Leading dimension of VS. Must be at least max(1, N). 231*> \endverbatim 232*> 233*> \param[out] VS1 234*> \verbatim 235*> VS1 is COMPLEX array, dimension (LDVS, N) 236*> VS1 holds another copy of the computed Schur vectors. 237*> \endverbatim 238*> 239*> \param[in] RCDEIN 240*> \verbatim 241*> RCDEIN is REAL 242*> When COMP = .TRUE. RCDEIN holds the precomputed reciprocal 243*> condition number for the average of selected eigenvalues. 244*> \endverbatim 245*> 246*> \param[in] RCDVIN 247*> \verbatim 248*> RCDVIN is REAL 249*> When COMP = .TRUE. RCDVIN holds the precomputed reciprocal 250*> condition number for the selected right invariant subspace. 251*> \endverbatim 252*> 253*> \param[in] NSLCT 254*> \verbatim 255*> NSLCT is INTEGER 256*> When COMP = .TRUE. the number of selected eigenvalues 257*> corresponding to the precomputed values RCDEIN and RCDVIN. 258*> \endverbatim 259*> 260*> \param[in] ISLCT 261*> \verbatim 262*> ISLCT is INTEGER array, dimension (NSLCT) 263*> When COMP = .TRUE. ISLCT selects the eigenvalues of the 264*> input matrix corresponding to the precomputed values RCDEIN 265*> and RCDVIN. For I=1, ... ,NSLCT, if ISLCT(I) = J, then the 266*> eigenvalue with the J-th largest real or imaginary part is 267*> selected. The real part is used if ISRT = 0, and the 268*> imaginary part if ISRT = 1. 269*> Not referenced if COMP = .FALSE. 270*> \endverbatim 271*> 272*> \param[in] ISRT 273*> \verbatim 274*> ISRT is INTEGER 275*> When COMP = .TRUE., ISRT describes how ISLCT is used to 276*> choose a subset of the spectrum. 277*> Not referenced if COMP = .FALSE. 278*> \endverbatim 279*> 280*> \param[out] RESULT 281*> \verbatim 282*> RESULT is REAL array, dimension (17) 283*> The values computed by the 17 tests described above. 284*> The values are currently limited to 1/ulp, to avoid 285*> overflow. 286*> \endverbatim 287*> 288*> \param[out] WORK 289*> \verbatim 290*> WORK is COMPLEX array, dimension (2*N*N) 291*> \endverbatim 292*> 293*> \param[in] LWORK 294*> \verbatim 295*> LWORK is INTEGER 296*> The number of entries in WORK to be passed to CGEESX. This 297*> must be at least 2*N, and N*(N+1)/2 if tests 14--16 are to 298*> be performed. 299*> \endverbatim 300*> 301*> \param[out] RWORK 302*> \verbatim 303*> RWORK is REAL array, dimension (N) 304*> \endverbatim 305*> 306*> \param[out] BWORK 307*> \verbatim 308*> BWORK is LOGICAL array, dimension (N) 309*> \endverbatim 310*> 311*> \param[out] INFO 312*> \verbatim 313*> INFO is INTEGER 314*> If 0, successful exit. 315*> If <0, input parameter -INFO had an incorrect value. 316*> If >0, CGEESX returned an error code, the absolute 317*> value of which is returned. 318*> \endverbatim 319* 320* Authors: 321* ======== 322* 323*> \author Univ. of Tennessee 324*> \author Univ. of California Berkeley 325*> \author Univ. of Colorado Denver 326*> \author NAG Ltd. 327* 328*> \date November 2011 329* 330*> \ingroup complex_eig 331* 332* ===================================================================== 333 SUBROUTINE CGET24( COMP, JTYPE, THRESH, ISEED, NOUNIT, N, A, LDA, 334 $ H, HT, W, WT, WTMP, VS, LDVS, VS1, RCDEIN, 335 $ RCDVIN, NSLCT, ISLCT, ISRT, RESULT, WORK, 336 $ LWORK, RWORK, BWORK, INFO ) 337* 338* -- LAPACK test routine (version 3.4.0) -- 339* -- LAPACK is a software package provided by Univ. of Tennessee, -- 340* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 341* November 2011 342* 343* .. Scalar Arguments .. 344 LOGICAL COMP 345 INTEGER INFO, ISRT, JTYPE, LDA, LDVS, LWORK, N, NOUNIT, 346 $ NSLCT 347 REAL RCDEIN, RCDVIN, THRESH 348* .. 349* .. Array Arguments .. 350 LOGICAL BWORK( * ) 351 INTEGER ISEED( 4 ), ISLCT( * ) 352 REAL RESULT( 17 ), RWORK( * ) 353 COMPLEX A( LDA, * ), H( LDA, * ), HT( LDA, * ), 354 $ VS( LDVS, * ), VS1( LDVS, * ), W( * ), 355 $ WORK( * ), WT( * ), WTMP( * ) 356* .. 357* 358* ===================================================================== 359* 360* .. Parameters .. 361 COMPLEX CZERO, CONE 362 PARAMETER ( CZERO = ( 0.0E+0, 0.0E+0 ), 363 $ CONE = ( 1.0E+0, 0.0E+0 ) ) 364 REAL ZERO, ONE 365 PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 ) 366 REAL EPSIN 367 PARAMETER ( EPSIN = 5.9605E-8 ) 368* .. 369* .. Local Scalars .. 370 CHARACTER SORT 371 INTEGER I, IINFO, ISORT, ITMP, J, KMIN, KNTEIG, RSUB, 372 $ SDIM, SDIM1 373 REAL ANORM, EPS, RCNDE1, RCNDV1, RCONDE, RCONDV, 374 $ SMLNUM, TOL, TOLIN, ULP, ULPINV, V, VRICMP, 375 $ VRIMIN, WNORM 376 COMPLEX CTMP 377* .. 378* .. Local Arrays .. 379 INTEGER IPNT( 20 ) 380* .. 381* .. External Functions .. 382 LOGICAL CSLECT 383 REAL CLANGE, SLAMCH 384 EXTERNAL CSLECT, CLANGE, SLAMCH 385* .. 386* .. External Subroutines .. 387 EXTERNAL CCOPY, CGEESX, CGEMM, CLACPY, CUNT01, XERBLA 388* .. 389* .. Intrinsic Functions .. 390 INTRINSIC ABS, AIMAG, MAX, MIN, REAL 391* .. 392* .. Arrays in Common .. 393 LOGICAL SELVAL( 20 ) 394 REAL SELWI( 20 ), SELWR( 20 ) 395* .. 396* .. Scalars in Common .. 397 INTEGER SELDIM, SELOPT 398* .. 399* .. Common blocks .. 400 COMMON / SSLCT / SELOPT, SELDIM, SELVAL, SELWR, SELWI 401* .. 402* .. Executable Statements .. 403* 404* Check for errors 405* 406 INFO = 0 407 IF( THRESH.LT.ZERO ) THEN 408 INFO = -3 409 ELSE IF( NOUNIT.LE.0 ) THEN 410 INFO = -5 411 ELSE IF( N.LT.0 ) THEN 412 INFO = -6 413 ELSE IF( LDA.LT.1 .OR. LDA.LT.N ) THEN 414 INFO = -8 415 ELSE IF( LDVS.LT.1 .OR. LDVS.LT.N ) THEN 416 INFO = -15 417 ELSE IF( LWORK.LT.2*N ) THEN 418 INFO = -24 419 END IF 420* 421 IF( INFO.NE.0 ) THEN 422 CALL XERBLA( 'CGET24', -INFO ) 423 RETURN 424 END IF 425* 426* Quick return if nothing to do 427* 428 DO 10 I = 1, 17 429 RESULT( I ) = -ONE 430 10 CONTINUE 431* 432 IF( N.EQ.0 ) 433 $ RETURN 434* 435* Important constants 436* 437 SMLNUM = SLAMCH( 'Safe minimum' ) 438 ULP = SLAMCH( 'Precision' ) 439 ULPINV = ONE / ULP 440* 441* Perform tests (1)-(13) 442* 443 SELOPT = 0 444 DO 90 ISORT = 0, 1 445 IF( ISORT.EQ.0 ) THEN 446 SORT = 'N' 447 RSUB = 0 448 ELSE 449 SORT = 'S' 450 RSUB = 6 451 END IF 452* 453* Compute Schur form and Schur vectors, and test them 454* 455 CALL CLACPY( 'F', N, N, A, LDA, H, LDA ) 456 CALL CGEESX( 'V', SORT, CSLECT, 'N', N, H, LDA, SDIM, W, VS, 457 $ LDVS, RCONDE, RCONDV, WORK, LWORK, RWORK, BWORK, 458 $ IINFO ) 459 IF( IINFO.NE.0 ) THEN 460 RESULT( 1+RSUB ) = ULPINV 461 IF( JTYPE.NE.22 ) THEN 462 WRITE( NOUNIT, FMT = 9998 )'CGEESX1', IINFO, N, JTYPE, 463 $ ISEED 464 ELSE 465 WRITE( NOUNIT, FMT = 9999 )'CGEESX1', IINFO, N, 466 $ ISEED( 1 ) 467 END IF 468 INFO = ABS( IINFO ) 469 RETURN 470 END IF 471 IF( ISORT.EQ.0 ) THEN 472 CALL CCOPY( N, W, 1, WTMP, 1 ) 473 END IF 474* 475* Do Test (1) or Test (7) 476* 477 RESULT( 1+RSUB ) = ZERO 478 DO 30 J = 1, N - 1 479 DO 20 I = J + 1, N 480 IF( H( I, J ).NE.CZERO ) 481 $ RESULT( 1+RSUB ) = ULPINV 482 20 CONTINUE 483 30 CONTINUE 484* 485* Test (2) or (8): Compute norm(A - Q*H*Q') / (norm(A) * N * ULP) 486* 487* Copy A to VS1, used as workspace 488* 489 CALL CLACPY( ' ', N, N, A, LDA, VS1, LDVS ) 490* 491* Compute Q*H and store in HT. 492* 493 CALL CGEMM( 'No transpose', 'No transpose', N, N, N, CONE, VS, 494 $ LDVS, H, LDA, CZERO, HT, LDA ) 495* 496* Compute A - Q*H*Q' 497* 498 CALL CGEMM( 'No transpose', 'Conjugate transpose', N, N, N, 499 $ -CONE, HT, LDA, VS, LDVS, CONE, VS1, LDVS ) 500* 501 ANORM = MAX( CLANGE( '1', N, N, A, LDA, RWORK ), SMLNUM ) 502 WNORM = CLANGE( '1', N, N, VS1, LDVS, RWORK ) 503* 504 IF( ANORM.GT.WNORM ) THEN 505 RESULT( 2+RSUB ) = ( WNORM / ANORM ) / ( N*ULP ) 506 ELSE 507 IF( ANORM.LT.ONE ) THEN 508 RESULT( 2+RSUB ) = ( MIN( WNORM, N*ANORM ) / ANORM ) / 509 $ ( N*ULP ) 510 ELSE 511 RESULT( 2+RSUB ) = MIN( WNORM / ANORM, REAL( N ) ) / 512 $ ( N*ULP ) 513 END IF 514 END IF 515* 516* Test (3) or (9): Compute norm( I - Q'*Q ) / ( N * ULP ) 517* 518 CALL CUNT01( 'Columns', N, N, VS, LDVS, WORK, LWORK, RWORK, 519 $ RESULT( 3+RSUB ) ) 520* 521* Do Test (4) or Test (10) 522* 523 RESULT( 4+RSUB ) = ZERO 524 DO 40 I = 1, N 525 IF( H( I, I ).NE.W( I ) ) 526 $ RESULT( 4+RSUB ) = ULPINV 527 40 CONTINUE 528* 529* Do Test (5) or Test (11) 530* 531 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 532 CALL CGEESX( 'N', SORT, CSLECT, 'N', N, HT, LDA, SDIM, WT, VS, 533 $ LDVS, RCONDE, RCONDV, WORK, LWORK, RWORK, BWORK, 534 $ IINFO ) 535 IF( IINFO.NE.0 ) THEN 536 RESULT( 5+RSUB ) = ULPINV 537 IF( JTYPE.NE.22 ) THEN 538 WRITE( NOUNIT, FMT = 9998 )'CGEESX2', IINFO, N, JTYPE, 539 $ ISEED 540 ELSE 541 WRITE( NOUNIT, FMT = 9999 )'CGEESX2', IINFO, N, 542 $ ISEED( 1 ) 543 END IF 544 INFO = ABS( IINFO ) 545 GO TO 220 546 END IF 547* 548 RESULT( 5+RSUB ) = ZERO 549 DO 60 J = 1, N 550 DO 50 I = 1, N 551 IF( H( I, J ).NE.HT( I, J ) ) 552 $ RESULT( 5+RSUB ) = ULPINV 553 50 CONTINUE 554 60 CONTINUE 555* 556* Do Test (6) or Test (12) 557* 558 RESULT( 6+RSUB ) = ZERO 559 DO 70 I = 1, N 560 IF( W( I ).NE.WT( I ) ) 561 $ RESULT( 6+RSUB ) = ULPINV 562 70 CONTINUE 563* 564* Do Test (13) 565* 566 IF( ISORT.EQ.1 ) THEN 567 RESULT( 13 ) = ZERO 568 KNTEIG = 0 569 DO 80 I = 1, N 570 IF( CSLECT( W( I ) ) ) 571 $ KNTEIG = KNTEIG + 1 572 IF( I.LT.N ) THEN 573 IF( CSLECT( W( I+1 ) ) .AND. 574 $ ( .NOT.CSLECT( W( I ) ) ) )RESULT( 13 ) = ULPINV 575 END IF 576 80 CONTINUE 577 IF( SDIM.NE.KNTEIG ) 578 $ RESULT( 13 ) = ULPINV 579 END IF 580* 581 90 CONTINUE 582* 583* If there is enough workspace, perform tests (14) and (15) 584* as well as (10) through (13) 585* 586 IF( LWORK.GE.( N*( N+1 ) ) / 2 ) THEN 587* 588* Compute both RCONDE and RCONDV with VS 589* 590 SORT = 'S' 591 RESULT( 14 ) = ZERO 592 RESULT( 15 ) = ZERO 593 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 594 CALL CGEESX( 'V', SORT, CSLECT, 'B', N, HT, LDA, SDIM1, WT, 595 $ VS1, LDVS, RCONDE, RCONDV, WORK, LWORK, RWORK, 596 $ BWORK, IINFO ) 597 IF( IINFO.NE.0 ) THEN 598 RESULT( 14 ) = ULPINV 599 RESULT( 15 ) = ULPINV 600 IF( JTYPE.NE.22 ) THEN 601 WRITE( NOUNIT, FMT = 9998 )'CGEESX3', IINFO, N, JTYPE, 602 $ ISEED 603 ELSE 604 WRITE( NOUNIT, FMT = 9999 )'CGEESX3', IINFO, N, 605 $ ISEED( 1 ) 606 END IF 607 INFO = ABS( IINFO ) 608 GO TO 220 609 END IF 610* 611* Perform tests (10), (11), (12), and (13) 612* 613 DO 110 I = 1, N 614 IF( W( I ).NE.WT( I ) ) 615 $ RESULT( 10 ) = ULPINV 616 DO 100 J = 1, N 617 IF( H( I, J ).NE.HT( I, J ) ) 618 $ RESULT( 11 ) = ULPINV 619 IF( VS( I, J ).NE.VS1( I, J ) ) 620 $ RESULT( 12 ) = ULPINV 621 100 CONTINUE 622 110 CONTINUE 623 IF( SDIM.NE.SDIM1 ) 624 $ RESULT( 13 ) = ULPINV 625* 626* Compute both RCONDE and RCONDV without VS, and compare 627* 628 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 629 CALL CGEESX( 'N', SORT, CSLECT, 'B', N, HT, LDA, SDIM1, WT, 630 $ VS1, LDVS, RCNDE1, RCNDV1, WORK, LWORK, RWORK, 631 $ BWORK, IINFO ) 632 IF( IINFO.NE.0 ) THEN 633 RESULT( 14 ) = ULPINV 634 RESULT( 15 ) = ULPINV 635 IF( JTYPE.NE.22 ) THEN 636 WRITE( NOUNIT, FMT = 9998 )'CGEESX4', IINFO, N, JTYPE, 637 $ ISEED 638 ELSE 639 WRITE( NOUNIT, FMT = 9999 )'CGEESX4', IINFO, N, 640 $ ISEED( 1 ) 641 END IF 642 INFO = ABS( IINFO ) 643 GO TO 220 644 END IF 645* 646* Perform tests (14) and (15) 647* 648 IF( RCNDE1.NE.RCONDE ) 649 $ RESULT( 14 ) = ULPINV 650 IF( RCNDV1.NE.RCONDV ) 651 $ RESULT( 15 ) = ULPINV 652* 653* Perform tests (10), (11), (12), and (13) 654* 655 DO 130 I = 1, N 656 IF( W( I ).NE.WT( I ) ) 657 $ RESULT( 10 ) = ULPINV 658 DO 120 J = 1, N 659 IF( H( I, J ).NE.HT( I, J ) ) 660 $ RESULT( 11 ) = ULPINV 661 IF( VS( I, J ).NE.VS1( I, J ) ) 662 $ RESULT( 12 ) = ULPINV 663 120 CONTINUE 664 130 CONTINUE 665 IF( SDIM.NE.SDIM1 ) 666 $ RESULT( 13 ) = ULPINV 667* 668* Compute RCONDE with VS, and compare 669* 670 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 671 CALL CGEESX( 'V', SORT, CSLECT, 'E', N, HT, LDA, SDIM1, WT, 672 $ VS1, LDVS, RCNDE1, RCNDV1, WORK, LWORK, RWORK, 673 $ BWORK, IINFO ) 674 IF( IINFO.NE.0 ) THEN 675 RESULT( 14 ) = ULPINV 676 IF( JTYPE.NE.22 ) THEN 677 WRITE( NOUNIT, FMT = 9998 )'CGEESX5', IINFO, N, JTYPE, 678 $ ISEED 679 ELSE 680 WRITE( NOUNIT, FMT = 9999 )'CGEESX5', IINFO, N, 681 $ ISEED( 1 ) 682 END IF 683 INFO = ABS( IINFO ) 684 GO TO 220 685 END IF 686* 687* Perform test (14) 688* 689 IF( RCNDE1.NE.RCONDE ) 690 $ RESULT( 14 ) = ULPINV 691* 692* Perform tests (10), (11), (12), and (13) 693* 694 DO 150 I = 1, N 695 IF( W( I ).NE.WT( I ) ) 696 $ RESULT( 10 ) = ULPINV 697 DO 140 J = 1, N 698 IF( H( I, J ).NE.HT( I, J ) ) 699 $ RESULT( 11 ) = ULPINV 700 IF( VS( I, J ).NE.VS1( I, J ) ) 701 $ RESULT( 12 ) = ULPINV 702 140 CONTINUE 703 150 CONTINUE 704 IF( SDIM.NE.SDIM1 ) 705 $ RESULT( 13 ) = ULPINV 706* 707* Compute RCONDE without VS, and compare 708* 709 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 710 CALL CGEESX( 'N', SORT, CSLECT, 'E', N, HT, LDA, SDIM1, WT, 711 $ VS1, LDVS, RCNDE1, RCNDV1, WORK, LWORK, RWORK, 712 $ BWORK, IINFO ) 713 IF( IINFO.NE.0 ) THEN 714 RESULT( 14 ) = ULPINV 715 IF( JTYPE.NE.22 ) THEN 716 WRITE( NOUNIT, FMT = 9998 )'CGEESX6', IINFO, N, JTYPE, 717 $ ISEED 718 ELSE 719 WRITE( NOUNIT, FMT = 9999 )'CGEESX6', IINFO, N, 720 $ ISEED( 1 ) 721 END IF 722 INFO = ABS( IINFO ) 723 GO TO 220 724 END IF 725* 726* Perform test (14) 727* 728 IF( RCNDE1.NE.RCONDE ) 729 $ RESULT( 14 ) = ULPINV 730* 731* Perform tests (10), (11), (12), and (13) 732* 733 DO 170 I = 1, N 734 IF( W( I ).NE.WT( I ) ) 735 $ RESULT( 10 ) = ULPINV 736 DO 160 J = 1, N 737 IF( H( I, J ).NE.HT( I, J ) ) 738 $ RESULT( 11 ) = ULPINV 739 IF( VS( I, J ).NE.VS1( I, J ) ) 740 $ RESULT( 12 ) = ULPINV 741 160 CONTINUE 742 170 CONTINUE 743 IF( SDIM.NE.SDIM1 ) 744 $ RESULT( 13 ) = ULPINV 745* 746* Compute RCONDV with VS, and compare 747* 748 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 749 CALL CGEESX( 'V', SORT, CSLECT, 'V', N, HT, LDA, SDIM1, WT, 750 $ VS1, LDVS, RCNDE1, RCNDV1, WORK, LWORK, RWORK, 751 $ BWORK, IINFO ) 752 IF( IINFO.NE.0 ) THEN 753 RESULT( 15 ) = ULPINV 754 IF( JTYPE.NE.22 ) THEN 755 WRITE( NOUNIT, FMT = 9998 )'CGEESX7', IINFO, N, JTYPE, 756 $ ISEED 757 ELSE 758 WRITE( NOUNIT, FMT = 9999 )'CGEESX7', IINFO, N, 759 $ ISEED( 1 ) 760 END IF 761 INFO = ABS( IINFO ) 762 GO TO 220 763 END IF 764* 765* Perform test (15) 766* 767 IF( RCNDV1.NE.RCONDV ) 768 $ RESULT( 15 ) = ULPINV 769* 770* Perform tests (10), (11), (12), and (13) 771* 772 DO 190 I = 1, N 773 IF( W( I ).NE.WT( I ) ) 774 $ RESULT( 10 ) = ULPINV 775 DO 180 J = 1, N 776 IF( H( I, J ).NE.HT( I, J ) ) 777 $ RESULT( 11 ) = ULPINV 778 IF( VS( I, J ).NE.VS1( I, J ) ) 779 $ RESULT( 12 ) = ULPINV 780 180 CONTINUE 781 190 CONTINUE 782 IF( SDIM.NE.SDIM1 ) 783 $ RESULT( 13 ) = ULPINV 784* 785* Compute RCONDV without VS, and compare 786* 787 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 788 CALL CGEESX( 'N', SORT, CSLECT, 'V', N, HT, LDA, SDIM1, WT, 789 $ VS1, LDVS, RCNDE1, RCNDV1, WORK, LWORK, RWORK, 790 $ BWORK, IINFO ) 791 IF( IINFO.NE.0 ) THEN 792 RESULT( 15 ) = ULPINV 793 IF( JTYPE.NE.22 ) THEN 794 WRITE( NOUNIT, FMT = 9998 )'CGEESX8', IINFO, N, JTYPE, 795 $ ISEED 796 ELSE 797 WRITE( NOUNIT, FMT = 9999 )'CGEESX8', IINFO, N, 798 $ ISEED( 1 ) 799 END IF 800 INFO = ABS( IINFO ) 801 GO TO 220 802 END IF 803* 804* Perform test (15) 805* 806 IF( RCNDV1.NE.RCONDV ) 807 $ RESULT( 15 ) = ULPINV 808* 809* Perform tests (10), (11), (12), and (13) 810* 811 DO 210 I = 1, N 812 IF( W( I ).NE.WT( I ) ) 813 $ RESULT( 10 ) = ULPINV 814 DO 200 J = 1, N 815 IF( H( I, J ).NE.HT( I, J ) ) 816 $ RESULT( 11 ) = ULPINV 817 IF( VS( I, J ).NE.VS1( I, J ) ) 818 $ RESULT( 12 ) = ULPINV 819 200 CONTINUE 820 210 CONTINUE 821 IF( SDIM.NE.SDIM1 ) 822 $ RESULT( 13 ) = ULPINV 823* 824 END IF 825* 826 220 CONTINUE 827* 828* If there are precomputed reciprocal condition numbers, compare 829* computed values with them. 830* 831 IF( COMP ) THEN 832* 833* First set up SELOPT, SELDIM, SELVAL, SELWR and SELWI so that 834* the logical function CSLECT selects the eigenvalues specified 835* by NSLCT, ISLCT and ISRT. 836* 837 SELDIM = N 838 SELOPT = 1 839 EPS = MAX( ULP, EPSIN ) 840 DO 230 I = 1, N 841 IPNT( I ) = I 842 SELVAL( I ) = .FALSE. 843 SELWR( I ) = REAL( WTMP( I ) ) 844 SELWI( I ) = AIMAG( WTMP( I ) ) 845 230 CONTINUE 846 DO 250 I = 1, N - 1 847 KMIN = I 848 IF( ISRT.EQ.0 ) THEN 849 VRIMIN = REAL( WTMP( I ) ) 850 ELSE 851 VRIMIN = AIMAG( WTMP( I ) ) 852 END IF 853 DO 240 J = I + 1, N 854 IF( ISRT.EQ.0 ) THEN 855 VRICMP = REAL( WTMP( J ) ) 856 ELSE 857 VRICMP = AIMAG( WTMP( J ) ) 858 END IF 859 IF( VRICMP.LT.VRIMIN ) THEN 860 KMIN = J 861 VRIMIN = VRICMP 862 END IF 863 240 CONTINUE 864 CTMP = WTMP( KMIN ) 865 WTMP( KMIN ) = WTMP( I ) 866 WTMP( I ) = CTMP 867 ITMP = IPNT( I ) 868 IPNT( I ) = IPNT( KMIN ) 869 IPNT( KMIN ) = ITMP 870 250 CONTINUE 871 DO 260 I = 1, NSLCT 872 SELVAL( IPNT( ISLCT( I ) ) ) = .TRUE. 873 260 CONTINUE 874* 875* Compute condition numbers 876* 877 CALL CLACPY( 'F', N, N, A, LDA, HT, LDA ) 878 CALL CGEESX( 'N', 'S', CSLECT, 'B', N, HT, LDA, SDIM1, WT, VS1, 879 $ LDVS, RCONDE, RCONDV, WORK, LWORK, RWORK, BWORK, 880 $ IINFO ) 881 IF( IINFO.NE.0 ) THEN 882 RESULT( 16 ) = ULPINV 883 RESULT( 17 ) = ULPINV 884 WRITE( NOUNIT, FMT = 9999 )'CGEESX9', IINFO, N, ISEED( 1 ) 885 INFO = ABS( IINFO ) 886 GO TO 270 887 END IF 888* 889* Compare condition number for average of selected eigenvalues 890* taking its condition number into account 891* 892 ANORM = CLANGE( '1', N, N, A, LDA, RWORK ) 893 V = MAX( REAL( N )*EPS*ANORM, SMLNUM ) 894 IF( ANORM.EQ.ZERO ) 895 $ V = ONE 896 IF( V.GT.RCONDV ) THEN 897 TOL = ONE 898 ELSE 899 TOL = V / RCONDV 900 END IF 901 IF( V.GT.RCDVIN ) THEN 902 TOLIN = ONE 903 ELSE 904 TOLIN = V / RCDVIN 905 END IF 906 TOL = MAX( TOL, SMLNUM / EPS ) 907 TOLIN = MAX( TOLIN, SMLNUM / EPS ) 908 IF( EPS*( RCDEIN-TOLIN ).GT.RCONDE+TOL ) THEN 909 RESULT( 16 ) = ULPINV 910 ELSE IF( RCDEIN-TOLIN.GT.RCONDE+TOL ) THEN 911 RESULT( 16 ) = ( RCDEIN-TOLIN ) / ( RCONDE+TOL ) 912 ELSE IF( RCDEIN+TOLIN.LT.EPS*( RCONDE-TOL ) ) THEN 913 RESULT( 16 ) = ULPINV 914 ELSE IF( RCDEIN+TOLIN.LT.RCONDE-TOL ) THEN 915 RESULT( 16 ) = ( RCONDE-TOL ) / ( RCDEIN+TOLIN ) 916 ELSE 917 RESULT( 16 ) = ONE 918 END IF 919* 920* Compare condition numbers for right invariant subspace 921* taking its condition number into account 922* 923 IF( V.GT.RCONDV*RCONDE ) THEN 924 TOL = RCONDV 925 ELSE 926 TOL = V / RCONDE 927 END IF 928 IF( V.GT.RCDVIN*RCDEIN ) THEN 929 TOLIN = RCDVIN 930 ELSE 931 TOLIN = V / RCDEIN 932 END IF 933 TOL = MAX( TOL, SMLNUM / EPS ) 934 TOLIN = MAX( TOLIN, SMLNUM / EPS ) 935 IF( EPS*( RCDVIN-TOLIN ).GT.RCONDV+TOL ) THEN 936 RESULT( 17 ) = ULPINV 937 ELSE IF( RCDVIN-TOLIN.GT.RCONDV+TOL ) THEN 938 RESULT( 17 ) = ( RCDVIN-TOLIN ) / ( RCONDV+TOL ) 939 ELSE IF( RCDVIN+TOLIN.LT.EPS*( RCONDV-TOL ) ) THEN 940 RESULT( 17 ) = ULPINV 941 ELSE IF( RCDVIN+TOLIN.LT.RCONDV-TOL ) THEN 942 RESULT( 17 ) = ( RCONDV-TOL ) / ( RCDVIN+TOLIN ) 943 ELSE 944 RESULT( 17 ) = ONE 945 END IF 946* 947 270 CONTINUE 948* 949 END IF 950* 951 9999 FORMAT( ' CGET24: ', A, ' returned INFO=', I6, '.', / 9X, 'N=', 952 $ I6, ', INPUT EXAMPLE NUMBER = ', I4 ) 953 9998 FORMAT( ' CGET24: ', A, ' returned INFO=', I6, '.', / 9X, 'N=', 954 $ I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5, ')' ) 955* 956 RETURN 957* 958* End of CGET24 959* 960 END 961