1
        2
        3
        4
        5
        6
        7
        8
        9
       10
       11
       12
       13
       14
       15
       16
       17
       18
       19
       20
       21
       22
       23
       24
       25
       26
       27
       28
       29
       30
       31
       32
       33
       34
       35
       36
       37
       38
       39
       40
       41
       42
       43
       44
       45
       46
       47
       48
       49
       50
       51
       52
       53
       54
       55
       56
       57
       58
       59
       60
       61
       62
       63
       64
       65
       66
       67
       68
       69
       70
       71
       72
       73
       74
       75
       76
       77
       78
       79
       80
       81
       82
       83
       84
       85
       86
       87
       88
       89
       90
       91
       92
       93
       94
       95
       96
       97
       98
       99
      100
      101
      102
      103
      104
      105
      106
      107
      108
      109
      110
      111
      112
      113
      114
      115
      116
      117
      118
      119
      120
      121
      122
      123
      124
      125
      126
      127
      128
      129
      130
      131
      132
      133
      134
      135
      136
      137
      138
      139
      140
      141
      142
      143
      144
      145
      146
      147
      148
      149
      150
      151
      152
      153
      154
      155
      156
      157
      158
      159
      160
      161
      162
      163
      164
      165
      166
      167
      168
      169
      170
      171
      172
      173
      174
      175
      176
      177
      178
      179
      180
      181
      182
      183
      184
      185
      186
      187
      188
      189
      190
      191
      192
      193
      194
      195
      196
      197
      198
      199
      200
      201
      202
      203
      204
      205
      206
      207
      208
      209
      210
      211
      212
      213
      214
      215
      216
      217
      218
      219
      220
      221
      222
      223
      224
      225
      226
      227
      228
      229
      230
      231
      232
      233
      234
      235
      236
      237
      238
      239
      240
      241
      242
      243
      244
      245
      246
      247
      248
      249
      250
      251
      252
      253
      254
      255
      256
      257
      258
      259
      260
      261
      262
      263
      264
      265
      266
      267
      268
      269
      270
      271
      272
      273
      274
      275
      276
      277
      278
      279
      280
      281
      282
      283
      284
      285
      286
      287
      288
      289
      290
      291
      292
      293
      294
      295
      296
      297
      298
      299
      300
      301
      302
      303
      304
      305
      306
      307
      308
      309
      310
      311
      312
      313
      314
      315
      316
      317
      318
      319
      320
      321
      322
      323
      324
      325
      326
      327
      328
      329
      330
      331
      332
      333
      334
      335
      336
      337
      338
      339
      340
      341
      342
      343
      344
      345
      346
      347
      348
      349
      350
      351
      352
      353
      354
      355
      356
      357
      358
      359
      360
      361
      362
      363
      364
      365
      366
      367
      368
      369
      370
      371
      372
      373
      374
      375
      376
      377
      378
      379
      380
      381
      382
      383
      384
      385
      386
      387
      388
      389
      390
      391
      392
      393
      394
      395
      396
      397
      398
      399
      400
      401
      402
      403
      404
      405
      406
      407
      408
      409
      410
      411
      412
      413
      414
      415
      416
      417
      418
      419
      420
      421
      422
      423
      424
      425
      426
      427
      428
      429
      430
      431
      432
      433
      434
      435
      436
      437
      438
      439
      440
      441
      442
      443
      444
      445
      446
      447
      448
      449
      450
      451
      452
      453
      454
      455
      456
      457
      458
      459
      460
      461
      462
      463
      464
      465
      466
      467
      468
      469
      470
      471
      472
      473
      474
      475
      476
      477
      478
      479
      480
      481
      482
      483
      484
      485
      486
      487
      488
      489
      490
      491
      492
      493
      494
      495
      496
      497
      498
      499
      500
      501
      502
      503
      504
      505
      506
      507
      508
      509
      510
      511
      512
      513
      514
      515
      516
      517
      518
      519
      520
      521
      522
      523
      524
      525
      526
      527
      528
      529
      530
      531
      532
      533
      534
      535
      536
      537
      538
      539
      540
      541
      542
      543
      544
      545
      546
      547
      548
      549
      550
      551
      552
      553
      554
      555
      556
      557
      558
      559
      560
      561
      562
      563
      564
      565
      566
      567
      568
      569
      570
      571
      572
      573
      574
      575
      576
      577
      578
      579
      580
      581
      582
      583
      584
      585
      586
      587
      588
      589
      590
      591
      592
      593
      594
      595
      596
      597
      598
      599
      600
      601
      602
      603
      604
      605
      606
      607
      608
      609
      610
      611
      612
      613
      614
      615
      616
      617
      618
      619
      620
      621
      622
      623
      624
      625
      626
      627
      628
      629
      630
      631
      632
      633
      634
      635
      636
      637
      638
      639
      640
      641
      642
      643
      644
      645
      646
      647
      648
      649
      650
      651
      652
      653
      654
      655
      656
      657
      658
      659
      660
      661
      662
      663
      664
      665
      666
      667
      668
      669
      670
      671
      672
      673
      674
      675
      676
      677
      678
      679
      680
      681
      682
      683
      684
      685
      686
      687
      688
      689
      690
      691
      692
      693
      694
      695
      696
      697
      698
      699
      700
      701
      702
      703
      704
      705
      706
      707
      708
      709
      710
      711
      712
      713
      714
      715
      716
      717
      718
      719
      720
      721
      722
      723
      724
      725
      726
      727
      728
      729
      730
      731
      732
      733
      734
      735
      736
      737
      738
      739
      740
      741
      742
      743
      744
      745
      746
      747
      748
      749
      750
      751
      752
      753
      754
      755
      756
      757
      758
      759
      760
      761
      762
      763
      764
      765
      766
      767
      768
      769
      770
      771
      772
      773
      774
      775
      776
      777
      778
      779
      780
      781
      782
      783
      784
      785
      786
      787
      788
      789
      790
      791
      792
      793
      794
      795
      796
      797
      798
      799
      800
      801
      802
      803
      804
      805
      806
      807
      808
      809
      810
      811
      812
      813
      814
      815
      816
      817
      818
      819
      820
      821
      822
      823
      824
      825
      826
      827
      828
      829
      830
      831
      832
      833
      834
      835
      836
      837
      838
      839
      840
      841
      842
      843
      844
      845
      846
      847
      848
      849
      850
      851
      852
      853
      854
      855
      856
      857
      858
      859
      860
      861
      862
      863
      864
      865
      866
      867
      868
      869
      870
      871
      872
      873
      874
      875
      876
      877
      878
      879
      880
      881
      882
      883
      884
      885
      886
      887
      888
      889
      890
      891
      892
      893
      894
      895
      896
      897
      898
      899
      900
      901
      902
      903
      904
      905
      906
      907
      908
      909
      910
      911
      912
      913
      914
      915
      916
      917
      918
      919
      920
      921
      922
      923
      924
      925
      926
      927
      928
      929
      930
      931
      932
      933
      934
      935
      936
      937
      938
      939
      940
      941
      942
      943
      944
      945
      946
      947
      948
      949
      950
      951
      952
      953
      954
      955
      956
      957
      958
      959
      960
      961
      962
      963
      964
      965
      966
      967
      968
      969
      970
      971
      972
      973
      974
      975
      976
      977
      978
      979
      980
      981
      982
      983
      984
      985
      986
      987
      988
      989
      990
      991
      992
      993
      994
      995
      996
      997
      998
      999
     1000
     1001
     1002
     1003
     1004
     1005
     1006
     1007
     1008
     1009
     1010
     1011
     1012
     1013
     1014
     1015
     1016
     1017
     1018
     1019
     1020
     1021
     1022
     1023
     1024
     1025
     1026
     1027
     1028
     1029
     1030
     1031
     1032
     1033
     1034
     1035
     1036
     1037
     1038
     1039
     1040
     1041
     1042
     1043
     1044
     1045
     1046
     1047
     1048
     1049
     1050
     1051
     1052
     1053
     1054
     1055
     1056
     1057
     1058
     1059
     1060
     1061
     1062
     1063
     1064
     1065
     1066
     1067
     1068
     1069
     1070
     1071
     1072
     1073
     1074
     1075
     1076
     1077
     1078
     1079
     1080
     1081
     1082
     1083
     1084
     1085
      SUBROUTINE CCHKGGNSIZESNNNTYPESDOTYPEISEEDTHRESH,
     $                   TSTDIFTHRSHNNOUNITALDABHTS1,
     $                   S2P1P2ULDUVQZALPHA1BETA1,
     $                   ALPHA3BETA3EVECTLEVECTRWORKLWORK,
     $                   RWORKLLWORKRESULTINFO )
*
*  -- LAPACK test routine (version 3.1) --
*     Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
*     November 2006
*
*     .. Scalar Arguments ..
      LOGICAL            TSTDIF
      INTEGER            INFOLDALDULWORKNOUNITNSIZESNTYPES
      REAL               THRESHTHRSHN
*     ..
*     .. Array Arguments ..
      LOGICAL            DOTYPE* ), LLWORK* )
      INTEGER            ISEED4 ), NN* )
      REAL               RESULT15 ), RWORK* )
      COMPLEX            ALDA* ), ALPHA1* ), ALPHA3* ),
     $                   BLDA* ), BETA1* ), BETA3* ),
     $                   EVECTLLDU* ), EVECTRLDU* ),
     $                   HLDA* ), P1LDA* ), P2LDA* ),
     $                   QLDU* ), S1LDA* ), S2LDA* ),
     $                   TLDA* ), ULDU* ), VLDU* ),
     $                   WORK* ), ZLDU* )
*     ..
*
*  Purpose
*  =======
*
*  CCHKGG  checks the nonsymmetric generalized eigenvalue problem
*  routines.
*                                 H          H        H
*  CGGHRD factors A and B as U H V  and U T V , where   means conjugate
*  transpose, H is hessenberg, T is triangular and U and V are unitary.
*
*                                  H          H
*  CHGEQZ factors H and T as  Q S Z  and Q P Z , where P and S are upper
*  triangular and Q and Z are unitary.  It also computes the generalized
*  eigenvalues (alpha(1),beta(1)),...,(alpha(n),beta(n)), where
*  alpha(j)=S(j,j) and beta(j)=P(j,j) -- thus, w(j) = alpha(j)/beta(j)
*  is a root of the generalized eigenvalue problem
*
*      det( A - w(j) B ) = 0
*
*  and m(j) = beta(j)/alpha(j) is a root of the essentially equivalent
*  problem
*
*      det( m(j) A - B ) = 0
*
*  CTGEVC computes the matrix L of left eigenvectors and the matrix R
*  of right eigenvectors for the matrix pair ( S, P ).  In the
*  description below,  l and r are left and right eigenvectors
*  corresponding to the generalized eigenvalues (alpha,beta).
*
*  When CCHKGG is called, a number of matrix "sizes" ("n's") and a
*  number of matrix "types" are specified.  For each size ("n")
*  and each type of matrix, one matrix will be generated and used
*  to test the nonsymmetric eigenroutines.  For each matrix, 13
*  tests will be performed.  The first twelve "test ratios" should be
*  small -- O(1).  They will be compared with the threshhold THRESH:
*
*                   H
*  (1)   | A - U H V  | / ( |A| n ulp )
*
*                   H
*  (2)   | B - U T V  | / ( |B| n ulp )
*
*                H
*  (3)   | I - UU  | / ( n ulp )
*
*                H
*  (4)   | I - VV  | / ( n ulp )
*
*                   H
*  (5)   | H - Q S Z  | / ( |H| n ulp )
*
*                   H
*  (6)   | T - Q P Z  | / ( |T| n ulp )
*
*                H
*  (7)   | I - QQ  | / ( n ulp )
*
*                H
*  (8)   | I - ZZ  | / ( n ulp )
*
*  (9)   max over all left eigenvalue/-vector pairs (beta/alpha,l) of
*                            H
*        | (beta A - alpha B) l | / ( ulp max( |beta A|, |alpha B| ) )
*
*  (10)  max over all left eigenvalue/-vector pairs (beta/alpha,l') of
*                            H
*        | (beta H - alpha T) l' | / ( ulp max( |beta H|, |alpha T| ) )
*
*        where the eigenvectors l' are the result of passing Q to
*        STGEVC and back transforming (JOB='B').
*
*  (11)  max over all right eigenvalue/-vector pairs (beta/alpha,r) of
*
*        | (beta A - alpha B) r | / ( ulp max( |beta A|, |alpha B| ) )
*
*  (12)  max over all right eigenvalue/-vector pairs (beta/alpha,r') of
*
*        | (beta H - alpha T) r' | / ( ulp max( |beta H|, |alpha T| ) )
*
*        where the eigenvectors r' are the result of passing Z to
*        STGEVC and back transforming (JOB='B').
*
*  The last three test ratios will usually be small, but there is no
*  mathematical requirement that they be so.  They are therefore
*  compared with THRESH only if TSTDIF is .TRUE.
*
*  (13)  | S(Q,Z computed) - S(Q,Z not computed) | / ( |S| ulp )
*
*  (14)  | P(Q,Z computed) - P(Q,Z not computed) | / ( |P| ulp )
*
*  (15)  max( |alpha(Q,Z computed) - alpha(Q,Z not computed)|/|S| ,
*             |beta(Q,Z computed) - beta(Q,Z not computed)|/|P| ) / ulp
*
*  In addition, the normalization of L and R are checked, and compared
*  with the threshhold THRSHN.
*
*  Test Matrices
*  ---- --------
*
*  The sizes of the test matrices are specified by an array
*  NN(1:NSIZES); the value of each element NN(j) specifies one size.
*  The "types" are specified by a logical array DOTYPE( 1:NTYPES ); if
*  DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
*  Currently, the list of possible types is:
*
*  (1)  ( 0, 0 )         (a pair of zero matrices)
*
*  (2)  ( I, 0 )         (an identity and a zero matrix)
*
*  (3)  ( 0, I )         (an identity and a zero matrix)
*
*  (4)  ( I, I )         (a pair of identity matrices)
*
*          t   t
*  (5)  ( J , J  )       (a pair of transposed Jordan blocks)
*
*                                      t                ( I   0  )
*  (6)  ( X, Y )         where  X = ( J   0  )  and Y = (      t )
*                                   ( 0   I  )          ( 0   J  )
*                        and I is a k x k identity and J a (k+1)x(k+1)
*                        Jordan block; k=(N-1)/2
*
*  (7)  ( D, I )         where D is P*D1, P is a random unitary diagonal
*                        matrix (i.e., with random magnitude 1 entries
*                        on the diagonal), and D1=diag( 0, 1,..., N-1 )
*                        (i.e., a diagonal matrix with D1(1,1)=0,
*                        D1(2,2)=1, ..., D1(N,N)=N-1.)
*  (8)  ( I, D )
*
*  (9)  ( big*D, small*I ) where "big" is near overflow and small=1/big
*
*  (10) ( small*D, big*I )
*
*  (11) ( big*I, small*D )
*
*  (12) ( small*I, big*D )
*
*  (13) ( big*D, big*I )
*
*  (14) ( small*D, small*I )
*
*  (15) ( D1, D2 )        where D1=P*diag( 0, 0, 1, ..., N-3, 0 ) and
*                         D2=Q*diag( 0, N-3, N-4,..., 1, 0, 0 ), and
*                         P and Q are random unitary diagonal matrices.
*            t   t
*  (16) U ( J , J ) V     where U and V are random unitary matrices.
*
*  (17) U ( T1, T2 ) V    where T1 and T2 are upper triangular matrices
*                         with random O(1) entries above the diagonal
*                         and diagonal entries diag(T1) =
*                         P*( 0, 0, 1, ..., N-3, 0 ) and diag(T2) =
*                         Q*( 0, N-3, N-4,..., 1, 0, 0 )
*
*  (18) U ( T1, T2 ) V    diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
*                         diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
*                         s = machine precision.
*
*  (19) U ( T1, T2 ) V    diag(T1)=( 0,0,1,1, 1-d, ..., 1-(N-5)*d=s, 0 )
*                         diag(T2) = ( 0, 1, 0, 1, ..., 1, 0 )
*
*                                                         N-5
*  (20) U ( T1, T2 ) V    diag(T1)=( 0, 0, 1, 1, a, ..., a   =s, 0 )
*                         diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
*
*  (21) U ( T1, T2 ) V    diag(T1)=( 0, 0, 1, r1, r2, ..., r(N-4), 0 )
*                         diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
*                         where r1,..., r(N-4) are random.
*
*  (22) U ( big*T1, small*T2 ) V   diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
*                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (23) U ( small*T1, big*T2 ) V   diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
*                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (24) U ( small*T1, small*T2 ) V diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
*                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (25) U ( big*T1, big*T2 ) V     diag(T1) = P*( 0, 0, 1, ..., N-3, 0 )
*                                  diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (26) U ( T1, T2 ) V     where T1 and T2 are random upper-triangular
*                          matrices.
*
*  Arguments
*  =========
*
*  NSIZES  (input) INTEGER
*          The number of sizes of matrices to use.  If it is zero,
*          CCHKGG does nothing.  It must be at least zero.
*
*  NN      (input) INTEGER array, dimension (NSIZES)
*          An array containing the sizes to be used for the matrices.
*          Zero values will be skipped.  The values must be at least
*          zero.
*
*  NTYPES  (input) INTEGER
*          The number of elements in DOTYPE.   If it is zero, CCHKGG
*          does nothing.  It must be at least zero.  If it is MAXTYP+1
*          and NSIZES is 1, then an additional type, MAXTYP+1 is
*          defined, which is to use whatever matrix is in A.  This
*          is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
*          DOTYPE(MAXTYP+1) is .TRUE. .
*
*  DOTYPE  (input) LOGICAL array, dimension (NTYPES)
*          If DOTYPE(j) is .TRUE., then for each size in NN a
*          matrix of that size and of type j will be generated.
*          If NTYPES is smaller than the maximum number of types
*          defined (PARAMETER MAXTYP), then types NTYPES+1 through
*          MAXTYP will not be generated.  If NTYPES is larger
*          than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
*          will be ignored.
*
*  ISEED   (input/output) INTEGER array, dimension (4)
*          On entry ISEED specifies the seed of the random number
*          generator. The array elements should be between 0 and 4095;
*          if not they will be reduced mod 4096.  Also, ISEED(4) must
*          be odd.  The random number generator uses a linear
*          congruential sequence limited to small integers, and so
*          should produce machine independent random numbers. The
*          values of ISEED are changed on exit, and can be used in the
*          next call to CCHKGG to continue the same random number
*          sequence.
*
*  THRESH  (input) REAL
*          A test will count as "failed" if the "error", computed as
*          described above, exceeds THRESH.  Note that the error
*          is scaled to be O(1), so THRESH should be a reasonably
*          small multiple of 1, e.g., 10 or 100.  In particular,
*          it should not depend on the precision (single vs. double)
*          or the size of the matrix.  It must be at least zero.
*
*  TSTDIF  (input) LOGICAL
*          Specifies whether test ratios 13-15 will be computed and
*          compared with THRESH.
*          = .FALSE.: Only test ratios 1-12 will be computed and tested.
*                     Ratios 13-15 will be set to zero.
*          = .TRUE.:  All the test ratios 1-15 will be computed and
*                     tested.
*
*  THRSHN  (input) REAL
*          Threshhold for reporting eigenvector normalization error.
*          If the normalization of any eigenvector differs from 1 by
*          more than THRSHN*ulp, then a special error message will be
*          printed.  (This is handled separately from the other tests,
*          since only a compiler or programming error should cause an
*          error message, at least if THRSHN is at least 5--10.)
*
*  NOUNIT  (input) INTEGER
*          The FORTRAN unit number for printing out error messages
*          (e.g., if a routine returns IINFO not equal to 0.)
*
*  A       (input/workspace) COMPLEX array, dimension (LDA, max(NN))
*          Used to hold the original A matrix.  Used as input only
*          if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
*          DOTYPE(MAXTYP+1)=.TRUE.
*
*  LDA     (input) INTEGER
*          The leading dimension of A, B, H, T, S1, P1, S2, and P2.
*          It must be at least 1 and at least max( NN ).
*
*  B       (input/workspace) COMPLEX array, dimension (LDA, max(NN))
*          Used to hold the original B matrix.  Used as input only
*          if NTYPES=MAXTYP+1, DOTYPE(1:MAXTYP)=.FALSE., and
*          DOTYPE(MAXTYP+1)=.TRUE.
*
*  H       (workspace) COMPLEX array, dimension (LDA, max(NN))
*          The upper Hessenberg matrix computed from A by CGGHRD.
*
*  T       (workspace) COMPLEX array, dimension (LDA, max(NN))
*          The upper triangular matrix computed from B by CGGHRD.
*
*  S1      (workspace) COMPLEX array, dimension (LDA, max(NN))
*          The Schur (upper triangular) matrix computed from H by CHGEQZ
*          when Q and Z are also computed.
*
*  S2      (workspace) COMPLEX array, dimension (LDA, max(NN))
*          The Schur (upper triangular) matrix computed from H by CHGEQZ
*          when Q and Z are not computed.
*
*  P1      (workspace) COMPLEX array, dimension (LDA, max(NN))
*          The upper triangular matrix computed from T by CHGEQZ
*          when Q and Z are also computed.
*
*  P2      (workspace) COMPLEX array, dimension (LDA, max(NN))
*          The upper triangular matrix computed from T by CHGEQZ
*          when Q and Z are not computed.
*
*  U       (workspace) COMPLEX array, dimension (LDU, max(NN))
*          The (left) unitary matrix computed by CGGHRD.
*
*  LDU     (input) INTEGER
*          The leading dimension of U, V, Q, Z, EVECTL, and EVECTR.  It
*          must be at least 1 and at least max( NN ).
*
*  V       (workspace) COMPLEX array, dimension (LDU, max(NN))
*          The (right) unitary matrix computed by CGGHRD.
*
*  Q       (workspace) COMPLEX array, dimension (LDU, max(NN))
*          The (left) unitary matrix computed by CHGEQZ.
*
*  Z       (workspace) COMPLEX array, dimension (LDU, max(NN))
*          The (left) unitary matrix computed by CHGEQZ.
*
*  ALPHA1  (workspace) COMPLEX array, dimension (max(NN))
*  BETA1   (workspace) COMPLEX array, dimension (max(NN))
*          The generalized eigenvalues of (A,B) computed by CHGEQZ
*          when Q, Z, and the full Schur matrices are computed.
*
*  ALPHA3  (workspace) COMPLEX array, dimension (max(NN))
*  BETA3   (workspace) COMPLEX array, dimension (max(NN))
*          The generalized eigenvalues of (A,B) computed by CHGEQZ
*          when neither Q, Z, nor the Schur matrices are computed.
*
*  EVECTL  (workspace) COMPLEX array, dimension (LDU, max(NN))
*          The (lower triangular) left eigenvector matrix for the
*          matrices in S1 and P1.
*
*  EVECTR  (workspace) COMPLEX array, dimension (LDU, max(NN))
*          The (upper triangular) right eigenvector matrix for the
*          matrices in S1 and P1.
*
*  WORK    (workspace) COMPLEX array, dimension (LWORK)
*
*  LWORK   (input) INTEGER
*          The number of entries in WORK.  This must be at least
*          max( 4*N, 2 * N**2, 1 ), for all N=NN(j).
*
*  RWORK   (workspace) REAL array, dimension (2*max(NN))
*
*  LLWORK  (workspace) LOGICAL array, dimension (max(NN))
*
*  RESULT  (output) REAL array, dimension (15)
*          The values computed by the tests described above.
*          The values are currently limited to 1/ulp, to avoid
*          overflow.
*
*  INFO    (output) INTEGER
*          = 0:  successful exit.
*          < 0:  if INFO = -i, the i-th argument had an illegal value.
*          > 0:  A routine returned an error code.  INFO is the
*                absolute value of the INFO value returned.
*
*  =====================================================================
*
*     .. Parameters ..
      REAL               ZEROONE
      PARAMETER          ( ZERO = 0.0E+0ONE = 1.0E+0 )
      COMPLEX            CZEROCONE
      PARAMETER          ( CZERO = ( 0.0E+00.0E+0 ),
     $                   CONE = ( 1.0E+00.0E+0 ) )
      INTEGER            MAXTYP
      PARAMETER          ( MAXTYP = 26 )
*     ..
*     .. Local Scalars ..
      LOGICAL            BADNN
      INTEGER            I1IADDIINFOINJJCJRJSIZEJTYPE,
     $                   LWKOPTMTYPESNN1NERRSNMATSNMAX,
     $                   NTESTNTESTT
      REAL               ANORMBNORMSAFMAXSAFMINTEMP1TEMP2,
     $                   ULPULPINV
      COMPLEX            CTEMP
*     ..
*     .. Local Arrays ..
      LOGICAL            LASIGNMAXTYP ), LBSIGNMAXTYP )
      INTEGER            IOLDSD4 ), KADD6 ), KAMAGNMAXTYP ),
     $                   KATYPEMAXTYP ), KAZEROMAXTYP ),
     $                   KBMAGNMAXTYP ), KBTYPEMAXTYP ),
     $                   KBZEROMAXTYP ), KCLASSMAXTYP ),
     $                   KTRIANMAXTYP ), KZ16 ), KZ26 )
      REAL               DUMMA4 ), RMAGN03 )
      COMPLEX            CDUMMA4 )
*     ..
*     .. External Functions ..
      REAL               CLANGESLAMCH
      COMPLEX            CLARND
      EXTERNAL           CLANGESLAMCHCLARND
*     ..
*     .. External Subroutines ..
      EXTERNAL           CGEQR2CGET51CGET52CGGHRDCHGEQZCLACPY,
     $                   CLARFGCLASETCLATM4CTGEVCCUNM2RSLABAD,
     $                   SLASUMXERBLA
*     ..
*     .. Intrinsic Functions ..
      INTRINSIC          ABSCONJGMAXMINREALSIGN
*     ..
*     .. Data statements ..
      DATA               KCLASS / 15*110*21*3 /
      DATA               KZ1 / 012133 /
      DATA               KZ2 / 001211 /
      DATA               KADD / 000032 /
      DATA               KATYPE / 0101234144114,
     $                   442458794*40 /
      DATA               KBTYPE / 00112-3141144,
     $                   11-42-48*80 /
      DATA               KAZERO / 6*1212*22*12*2313,
     $                   4*54*31 /
      DATA               KBZERO / 6*1122*12*22*1414,
     $                   4*64*41 /
      DATA               KAMAGN / 8*12323237*1233,
     $                   21 /
      DATA               KBMAGN / 8*13232237*1323,
     $                   21 /
      DATA               KTRIAN / 16*010*1 /
      DATA               LASIGN / 6*.FALSE..TRUE..FALSE.2*.TRUE.,
     $                   2*.FALSE.3*.TRUE..FALSE..TRUE.,
     $                   3*.FALSE.5*.TRUE..FALSE. /
      DATA               LBSIGN / 7*.FALSE..TRUE.2*.FALSE.,
     $                   2*.TRUE.2*.FALSE..TRUE..FALSE..TRUE.,
     $                   9*.FALSE. /
*     ..
*     .. Executable Statements ..
*
*     Check for errors
*
      INFO = 0
*
      BADNN = .FALSE.
      NMAX = 1
      DO 10 J = 1NSIZES
         NMAX = MAXNMAXNNJ ) )
         IFNNJ ).LT.0 )
     $      BADNN = .TRUE.
   10 CONTINUE
*
      LWKOPT = MAX2*NMAX*NMAX4*NMAX1 )
*
*     Check for errors
*
      IFNSIZES.LT.0 ) THEN
         INFO = -1
      ELSE IFBADNN ) THEN
         INFO = -2
      ELSE IFNTYPES.LT.0 ) THEN
         INFO = -3
      ELSE IFTHRESH.LT.ZERO ) THEN
         INFO = -6
      ELSE IFLDA.LE.1 .OR. LDA.LT.NMAX ) THEN
         INFO = -10
      ELSE IFLDU.LE.1 .OR. LDU.LT.NMAX ) THEN
         INFO = -19
      ELSE IFLWKOPT.GT.LWORK ) THEN
         INFO = -30
      END IF
*
      IFINFO.NE.0 ) THEN
         CALL XERBLA'CCHKGG'-INFO )
         RETURN
      END IF
*
*     Quick return if possible
*
      IFNSIZES.EQ.0 .OR. NTYPES.EQ.0 )
     $   RETURN
*
      SAFMIN = SLAMCH'Safe minimum' )
      ULP = SLAMCH'Epsilon' )*SLAMCH'Base' )
      SAFMIN = SAFMIN / ULP
      SAFMAX = ONE / SAFMIN
      CALL SLABADSAFMINSAFMAX )
      ULPINV = ONE / ULP
*
*     The values RMAGN(2:3) depend on N, see below.
*
      RMAGN0 ) = ZERO
      RMAGN1 ) = ONE
*
*     Loop over sizes, types
*
      NTESTT = 0
      NERRS = 0
      NMATS = 0
*
      DO 240 JSIZE = 1NSIZES
         N = NNJSIZE )
         N1 = MAX1N )
         RMAGN2 ) = SAFMAX*ULP / REALN1 )
         RMAGN3 ) = SAFMIN*ULPINV*N1
*
         IFNSIZES.NE.1 ) THEN
            MTYPES = MINMAXTYPNTYPES )
         ELSE
            MTYPES = MINMAXTYP+1NTYPES )
         END IF
*
         DO 230 JTYPE = 1MTYPES
            IF.NOT.DOTYPEJTYPE ) )
     $         GO TO 230
            NMATS = NMATS + 1
            NTEST = 0
*
*           Save ISEED in case of an error.
*
            DO 20 J = 14
               IOLDSDJ ) = ISEEDJ )
   20       CONTINUE
*
*           Initialize RESULT
*
            DO 30 J = 115
               RESULTJ ) = ZERO
   30       CONTINUE
*
*           Compute A and B
*
*           Description of control parameters:
*
*           KCLASS: =1 means w/o rotation, =2 means w/ rotation,
*                   =3 means random.
*           KATYPE: the "type" to be passed to CLATM4 for computing A.
*           KAZERO: the pattern of zeros on the diagonal for A:
*                   =1: ( xxx ), =2: (0, xxx ) =3: ( 0, 0, xxx, 0 ),
*                   =4: ( 0, xxx, 0, 0 ), =5: ( 0, 0, 1, xxx, 0 ),
*                   =6: ( 0, 1, 0, xxx, 0 ).  (xxx means a string of
*                   non-zero entries.)
*           KAMAGN: the magnitude of the matrix: =0: zero, =1: O(1),
*                   =2: large, =3: small.
*           LASIGN: .TRUE. if the diagonal elements of A are to be
*                   multiplied by a random magnitude 1 number.
*           KBTYPE, KBZERO, KBMAGN, LBSIGN: the same, but for B.
*           KTRIAN: =0: don't fill in the upper triangle, =1: do.
*           KZ1, KZ2, KADD: used to implement KAZERO and KBZERO.
*           RMAGN:  used to implement KAMAGN and KBMAGN.
*
            IFMTYPES.GT.MAXTYP )
     $         GO TO 110
            IINFO = 0
            IFKCLASSJTYPE ).LT.3 ) THEN
*
*              Generate A (w/o rotation)
*
               IFABSKATYPEJTYPE ) ).EQ.3 ) THEN
                  IN = 2*( ( N-1 ) / 2 ) + 1
                  IFIN.NE.N )
     $               CALL CLASET'Full'NNCZEROCZEROALDA )
               ELSE
                  IN = N
               END IF
               CALL CLATM4KATYPEJTYPE ), INKZ1KAZEROJTYPE ) ),
     $                      KZ2KAZEROJTYPE ) ), LASIGNJTYPE ),
     $                      RMAGNKAMAGNJTYPE ) ), ULP,
     $                      RMAGNKTRIANJTYPE )*KAMAGNJTYPE ) ), 4,
     $                      ISEEDALDA )
               IADD = KADDKAZEROJTYPE ) )
               IFIADD.GT.0 .AND. IADD.LE.N )
     $            AIADDIADD ) = RMAGNKAMAGNJTYPE ) )
*
*              Generate B (w/o rotation)
*
               IFABSKBTYPEJTYPE ) ).EQ.3 ) THEN
                  IN = 2*( ( N-1 ) / 2 ) + 1
                  IFIN.NE.N )
     $               CALL CLASET'Full'NNCZEROCZEROBLDA )
               ELSE
                  IN = N
               END IF
               CALL CLATM4KBTYPEJTYPE ), INKZ1KBZEROJTYPE ) ),
     $                      KZ2KBZEROJTYPE ) ), LBSIGNJTYPE ),
     $                      RMAGNKBMAGNJTYPE ) ), ONE,
     $                      RMAGNKTRIANJTYPE )*KBMAGNJTYPE ) ), 4,
     $                      ISEEDBLDA )
               IADD = KADDKBZEROJTYPE ) )
               IFIADD.NE.0 )
     $            BIADDIADD ) = RMAGNKBMAGNJTYPE ) )
*
               IFKCLASSJTYPE ).EQ.2 .AND. N.GT.0 ) THEN
*
*                 Include rotations
*
*                 Generate U, V as Householder transformations times a
*                 diagonal matrix.  (Note that CLARFG makes U(j,j) and
*                 V(j,j) real.)
*
                  DO 50 JC = 1N - 1
                     DO 40 JR = JCN
                        UJRJC ) = CLARND3ISEED )
                        VJRJC ) = CLARND3ISEED )
   40                CONTINUE
                     CALL CLARFGN+1-JCUJCJC ), UJC+1JC ), 1,
     $                            WORKJC ) )
                     WORK2*N+JC ) = SIGNONEREALUJCJC ) ) )
                     UJCJC ) = CONE
                     CALL CLARFGN+1-JCVJCJC ), VJC+1JC ), 1,
     $                            WORKN+JC ) )
                     WORK3*N+JC ) = SIGNONEREALVJCJC ) ) )
                     VJCJC ) = CONE
   50             CONTINUE
                  CTEMP = CLARND3ISEED )
                  UNN ) = CONE
                  WORKN ) = CZERO
                  WORK3*N ) = CTEMP / ABSCTEMP )
                  CTEMP = CLARND3ISEED )
                  VNN ) = CONE
                  WORK2*N ) = CZERO
                  WORK4*N ) = CTEMP / ABSCTEMP )
*
*                 Apply the diagonal matrices
*
                  DO 70 JC = 1N
                     DO 60 JR = 1N
                        AJRJC ) = WORK2*N+JR )*
     $                                CONJGWORK3*N+JC ) )*
     $                                AJRJC )
                        BJRJC ) = WORK2*N+JR )*
     $                                CONJGWORK3*N+JC ) )*
     $                                BJRJC )
   60                CONTINUE
   70             CONTINUE
                  CALL CUNM2R'L''N'NNN-1ULDUWORKA,
     $                         LDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
                  CALL CUNM2R'R''C'NNN-1VLDUWORKN+1 ),
     $                         ALDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
                  CALL CUNM2R'L''N'NNN-1ULDUWORKB,
     $                         LDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
                  CALL CUNM2R'R''C'NNN-1VLDUWORKN+1 ),
     $                         BLDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
               END IF
            ELSE
*
*              Random matrices
*
               DO 90 JC = 1N
                  DO 80 JR = 1N
                     AJRJC ) = RMAGNKAMAGNJTYPE ) )*
     $                             CLARND4ISEED )
                     BJRJC ) = RMAGNKBMAGNJTYPE ) )*
     $                             CLARND4ISEED )
   80             CONTINUE
   90          CONTINUE
            END IF
*
            ANORM = CLANGE'1'NNALDARWORK )
            BNORM = CLANGE'1'NNBLDARWORK )
*
  100       CONTINUE
*
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'Generator', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               RETURN
            END IF
*
  110       CONTINUE
*
*           Call CGEQR2, CUNM2R, and CGGHRD to compute H, T, U, and V
*
            CALL CLACPY' 'NNALDAHLDA )
            CALL CLACPY' 'NNBLDATLDA )
            NTEST = 1
            RESULT1 ) = ULPINV
*
            CALL CGEQR2NNTLDAWORKWORKN+1 ), IINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CGEQR2', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CUNM2R'L''C'NNNTLDAWORKHLDA,
     $                   WORKN+1 ), IINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CUNM2R', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CLASET'Full'NNCZEROCONEULDU )
            CALL CUNM2R'R''N'NNNTLDAWORKULDU,
     $                   WORKN+1 ), IINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CUNM2R', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CGGHRD'V''I'N1NHLDATLDAULDUV,
     $                   LDUIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CGGHRD', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
            NTEST = 4
*
*           Do tests 1--4
*
            CALL CGET511NALDAHLDAULDUVLDUWORK,
     $                   RWORKRESULT1 ) )
            CALL CGET511NBLDATLDAULDUVLDUWORK,
     $                   RWORKRESULT2 ) )
            CALL CGET513NBLDATLDAULDUULDUWORK,
     $                   RWORKRESULT3 ) )
            CALL CGET513NBLDATLDAVLDUVLDUWORK,
     $                   RWORKRESULT4 ) )
*
*           Call CHGEQZ to compute S1, P1, S2, P2, Q, and Z, do tests.
*
*           Compute T1 and UZ
*
*           Eigenvalues only
*
            CALL CLACPY' 'NNHLDAS2LDA )
            CALL CLACPY' 'NNTLDAP2LDA )
            NTEST = 5
            RESULT5 ) = ULPINV
*
            CALL CHGEQZ'E''N''N'N1NS2LDAP2LDA,
     $                   ALPHA3BETA3QLDUZLDUWORKLWORK,
     $                   RWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CHGEQZ(E)', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
*           Eigenvalues and Full Schur Form
*
            CALL CLACPY' 'NNHLDAS2LDA )
            CALL CLACPY' 'NNTLDAP2LDA )
*
            CALL CHGEQZ'S''N''N'N1NS2LDAP2LDA,
     $                   ALPHA1BETA1QLDUZLDUWORKLWORK,
     $                   RWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CHGEQZ(S)', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
*           Eigenvalues, Schur Form, and Schur Vectors
*
            CALL CLACPY' 'NNHLDAS1LDA )
            CALL CLACPY' 'NNTLDAP1LDA )
*
            CALL CHGEQZ'S''I''I'N1NS1LDAP1LDA,
     $                   ALPHA1BETA1QLDUZLDUWORKLWORK,
     $                   RWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CHGEQZ(V)', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            NTEST = 8
*
*           Do Tests 5--8
*
            CALL CGET511NHLDAS1LDAQLDUZLDUWORK,
     $                   RWORKRESULT5 ) )
            CALL CGET511NTLDAP1LDAQLDUZLDUWORK,
     $                   RWORKRESULT6 ) )
            CALL CGET513NTLDAP1LDAQLDUQLDUWORK,
     $                   RWORKRESULT7 ) )
            CALL CGET513NTLDAP1LDAZLDUZLDUWORK,
     $                   RWORKRESULT8 ) )
*
*           Compute the Left and Right Eigenvectors of (S1,P1)
*
*           9: Compute the left eigenvector Matrix without
*              back transforming:
*
            NTEST = 9
            RESULT9 ) = ULPINV
*
*           To test "SELECT" option, compute half of the eigenvectors
*           in one call, and half in another
*
            I1 = N / 2
            DO 120 J = 1I1
               LLWORKJ ) = .TRUE.
  120       CONTINUE
            DO 130 J = I1 + 1N
               LLWORKJ ) = .FALSE.
  130       CONTINUE
*
            CALL CTGEVC'L''S'LLWORKNS1LDAP1LDAEVECTL,
     $                   LDUCDUMMALDUNINWORKRWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CTGEVC(L,S1)', IINFO, N,
     $            JTYPE, IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            I1 = IN
            DO 140 J = 1I1
               LLWORKJ ) = .FALSE.
  140       CONTINUE
            DO 150 J = I1 + 1N
               LLWORKJ ) = .TRUE.
  150       CONTINUE
*
            CALL CTGEVC'L''S'LLWORKNS1LDAP1LDA,
     $                   EVECTL1I1+1 ), LDUCDUMMALDUNIN,
     $                   WORKRWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CTGEVC(L,S2)', IINFO, N,
     $            JTYPE, IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CGET52.TRUE.NS1LDAP1LDAEVECTLLDU,
     $                   ALPHA1BETA1WORKRWORKDUMMA1 ) )
            RESULT9 ) = DUMMA1 )
            IFDUMMA2 ).GT.THRSHN ) THEN
               WRITENOUNIT, FMT = 9998 )'Left''CTGEVC(HOWMNY=S)',
     $            DUMMA( 2 ), N, JTYPE, IOLDSD
            END IF
*
*           10: Compute the left eigenvector Matrix with
*               back transforming:
*
            NTEST = 10
            RESULT10 ) = ULPINV
            CALL CLACPY'F'NNQLDUEVECTLLDU )
            CALL CTGEVC'L''B'LLWORKNS1LDAP1LDAEVECTL,
     $                   LDUCDUMMALDUNINWORKRWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CTGEVC(L,B)', IINFO, N,
     $            JTYPE, IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CGET52.TRUE.NHLDATLDAEVECTLLDUALPHA1,
     $                   BETA1WORKRWORKDUMMA1 ) )
            RESULT10 ) = DUMMA1 )
            IFDUMMA2 ).GT.THRSHN ) THEN
               WRITENOUNIT, FMT = 9998 )'Left''CTGEVC(HOWMNY=B)',
     $            DUMMA( 2 ), N, JTYPE, IOLDSD
            END IF
*
*           11: Compute the right eigenvector Matrix without
*               back transforming:
*
            NTEST = 11
            RESULT11 ) = ULPINV
*
*           To test "SELECT" option, compute half of the eigenvectors
*           in one call, and half in another
*
            I1 = N / 2
            DO 160 J = 1I1
               LLWORKJ ) = .TRUE.
  160       CONTINUE
            DO 170 J = I1 + 1N
               LLWORKJ ) = .FALSE.
  170       CONTINUE
*
            CALL CTGEVC'R''S'LLWORKNS1LDAP1LDACDUMMA,
     $                   LDUEVECTRLDUNINWORKRWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CTGEVC(R,S1)', IINFO, N,
     $            JTYPE, IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            I1 = IN
            DO 180 J = 1I1
               LLWORKJ ) = .FALSE.
  180       CONTINUE
            DO 190 J = I1 + 1N
               LLWORKJ ) = .TRUE.
  190       CONTINUE
*
            CALL CTGEVC'R''S'LLWORKNS1LDAP1LDACDUMMA,
     $                   LDUEVECTR1I1+1 ), LDUNINWORK,
     $                   RWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CTGEVC(R,S2)', IINFO, N,
     $            JTYPE, IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CGET52.FALSE.NS1LDAP1LDAEVECTRLDU,
     $                   ALPHA1BETA1WORKRWORKDUMMA1 ) )
            RESULT11 ) = DUMMA1 )
            IFDUMMA2 ).GT.THRESH ) THEN
               WRITENOUNIT, FMT = 9998 )'Right''CTGEVC(HOWMNY=S)',
     $            DUMMA( 2 ), N, JTYPE, IOLDSD
            END IF
*
*           12: Compute the right eigenvector Matrix with
*               back transforming:
*
            NTEST = 12
            RESULT12 ) = ULPINV
            CALL CLACPY'F'NNZLDUEVECTRLDU )
            CALL CTGEVC'R''B'LLWORKNS1LDAP1LDACDUMMA,
     $                   LDUEVECTRLDUNINWORKRWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'CTGEVC(R,B)', IINFO, N,
     $            JTYPE, IOLDSD
               INFO = ABSIINFO )
               GO TO 210
            END IF
*
            CALL CGET52.FALSE.NHLDATLDAEVECTRLDU,
     $                   ALPHA1BETA1WORKRWORKDUMMA1 ) )
            RESULT12 ) = DUMMA1 )
            IFDUMMA2 ).GT.THRESH ) THEN
               WRITENOUNIT, FMT = 9998 )'Right''CTGEVC(HOWMNY=B)',
     $            DUMMA( 2 ), N, JTYPE, IOLDSD
            END IF
*
*           Tests 13--15 are done only on request
*
            IFTSTDIF ) THEN
*
*              Do Tests 13--14
*
               CALL CGET512NS1LDAS2LDAQLDUZLDU,
     $                      WORKRWORKRESULT13 ) )
               CALL CGET512NP1LDAP2LDAQLDUZLDU,
     $                      WORKRWORKRESULT14 ) )
*
*              Do Test 15
*
               TEMP1 = ZERO
               TEMP2 = ZERO
               DO 200 J = 1N
                  TEMP1 = MAXTEMP1ABSALPHA1J )-ALPHA3J ) ) )
                  TEMP2 = MAXTEMP2ABSBETA1J )-BETA3J ) ) )
  200          CONTINUE
*
               TEMP1 = TEMP1 / MAXSAFMINULP*MAXTEMP1ANORM ) )
               TEMP2 = TEMP2 / MAXSAFMINULP*MAXTEMP2BNORM ) )
               RESULT15 ) = MAXTEMP1TEMP2 )
               NTEST = 15
            ELSE
               RESULT13 ) = ZERO
               RESULT14 ) = ZERO
               RESULT15 ) = ZERO
               NTEST = 12
            END IF
*
*           End of Loop -- Check for RESULT(j) > THRESH
*
  210       CONTINUE
*
            NTESTT = NTESTT + NTEST
*
*           Print out tests which fail.
*
            DO 220 JR = 1NTEST
               IFRESULTJR ).GE.THRESH ) THEN
*
*                 If this is the first test to fail,
*                 print a header to the data file.
*
                  IFNERRS.EQ.0 ) THEN
                     WRITENOUNIT, FMT = 9997 )'CGG'
*
*                    Matrix types
*
                     WRITENOUNIT, FMT = 9996 )
                     WRITENOUNIT, FMT = 9995 )
                     WRITENOUNIT, FMT = 9994 )'Unitary'
*
*                    Tests performed
*
                     WRITENOUNIT, FMT = 9993 )'unitary''*',
     $                  'conjugate transpose', ( '*', J = 110 )
*
                  END IF
                  NERRS = NERRS + 1
                  IFRESULTJR ).LT.10000.0 ) THEN
                     WRITENOUNIT, FMT = 9992 )N, JTYPE, IOLDSD, JR,
     $                  RESULT( JR )
                  ELSE
                     WRITENOUNIT, FMT = 9991 )N, JTYPE, IOLDSD, JR,
     $                  RESULT( JR )
                  END IF
               END IF
  220       CONTINUE
*
  230    CONTINUE
  240 CONTINUE
*
*     Summary
*
      CALL SLASUM'CGG'NOUNITNERRSNTESTT )
      RETURN
*
 9999 FORMAT( ' CCHKGG: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
     $      I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5, ')' )
*
 9998 FORMAT( ' CCHKGG: ', A, ' Eigenvectors from ', A, ' incorrectly ',
     $      'normalized.', / ' Bits of error=', 0P, G10.3, ',', 9X,
     $      'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5,
     $      ')' )
*
 9997 FORMAT( 1X, A3, ' -- Complex Generalized eigenvalue problem' )
*
 9996 FORMAT( ' Matrix types (see CCHKGG for details): ' )
*
 9995 FORMAT( ' Special Matrices:', 23X,
     $      '(J''=transposed Jordan block)',
     $      / '   1=(0,0)  2=(I,0)  3=(0,I)  4=(I,I)  5=(J'',J'')  ',
     $      '6=(diag(J'',I), diag(I,J''))', / ' Diagonal Matrices:  ( ',
     $      'D=diag(0,1,2,...) )', / '   7=(D,I)   9=(large*D, small*I',
     $      ')  11=(large*I, small*D)  13=(large*D, large*I)', /
     $      '   8=(I,D)  10=(small*D, large*I)  12=(small*I, large*D) ',
     $      ' 14=(small*D, small*I)', / '  15=(D, reversed D)' )
 9994 FORMAT( ' Matrices Rotated by Random ', A, ' Matrices U, V:',
     $      / '  16=Transposed Jordan Blocks             19=geometric ',
     $      'alpha, beta=0,1', / '  17=arithm. alpha&beta             ',
     $      '      20=arithmetic alpha, beta=0,1', / '  18=clustered ',
     $      'alpha, beta=0,1            21=random alpha, beta=0,1',
     $      / ' Large & Small Matrices:', / '  22=(large, small)   ',
     $      '23=(small,large)    24=(small,small)    25=(large,large)',
     $      / '  26=random O(1) matrices.' )
*
 9993 FORMAT( / ' Tests performed:   (H is Hessenberg, S is Schur, B, ',
     $      'T, P are triangular,', / 20X, 'U, V, Q, and Z are ', A,
     $      ', l and r are the', / 20X,
     $      'appropriate left and right eigenvectors, resp., a is',
     $      / 20X, 'alpha, b is beta, and ', A, ' means ', A, '.)',
     $      / ' 1 = | A - U H V', A,
     $      ' | / ( |A| n ulp )      2 = | B - U T V', A,
     $      ' | / ( |B| n ulp )', / ' 3 = | I - UU', A,
     $      ' | / ( n ulp )             4 = | I - VV', A,
     $      ' | / ( n ulp )', / ' 5 = | H - Q S Z', A,
     $      ' | / ( |H| n ulp )', 6X, '6 = | T - Q P Z', A,
     $      ' | / ( |T| n ulp )', / ' 7 = | I - QQ', A,
     $      ' | / ( n ulp )             8 = | I - ZZ', A,
     $      ' | / ( n ulp )', / ' 9 = max | ( b S - a P )', A,
     $      ' l | / const.  10 = max | ( b H - a T )', A,
     $      ' l | / const.', /
     $      ' 11= max | ( b S - a P ) r | / const.   12 = max | ( b H',
     $      ' - a T ) r | / const.', / 1X )
*
 9992 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
     $      4( I4, ',' ), ' result ', I2, ' is', 0P, F8.2 )
 9991 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
     $      4( I4, ',' ), ' result ', I2, ' is', 1P, E10.3 )
*
*     End of CCHKGG
*
      END