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
      SUBROUTINE DDRVGGNSIZESNNNTYPESDOTYPEISEEDTHRESH,
     $                   THRSHNNOUNITALDABSTS2T2Q,
     $                   LDQZALPHR1ALPHI1BETA1ALPHR2ALPHI2,
     $                   BETA2VLVRWORKLWORKRESULTINFO )
*
*  -- LAPACK test routine (version 3.1) --
*     Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd..
*     November 2006
*
*     .. Scalar Arguments ..
      INTEGER            INFOLDALDQLWORKNOUNITNSIZESNTYPES
      DOUBLE PRECISION   THRESHTHRSHN
*     ..
*     .. Array Arguments ..
      LOGICAL            DOTYPE* )
      INTEGER            ISEED4 ), NN* )
      DOUBLE PRECISION   ALDA* ), ALPHI1* ), ALPHI2* ),
     $                   ALPHR1* ), ALPHR2* ), BLDA* ),
     $                   BETA1* ), BETA2* ), QLDQ* ),
     $                   RESULT* ), SLDA* ), S2LDA* ),
     $                   TLDA* ), T2LDA* ), VLLDQ* ),
     $                   VRLDQ* ), WORK* ), ZLDQ* )
*     ..
*
*  Purpose
*  =======
*
*  DDRVGG  checks the nonsymmetric generalized eigenvalue driver
*  routines.
*                                T          T        T
*  DGEGS factors A and B as Q S Z  and Q T Z , where   means
*  transpose, T is upper triangular, S is in generalized Schur form
*  (block upper triangular, with 1x1 and 2x2 blocks on the diagonal,
*  the 2x2 blocks corresponding to complex conjugate pairs of
*  generalized eigenvalues), and Q and Z are orthogonal.  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
*
*  DGEGV computes the generalized eigenvalues (alpha(1),beta(1)), ...,
*  (alpha(n),beta(n)), the matrix L whose columns contain the
*  generalized left eigenvectors l, and the matrix R whose columns
*  contain the generalized right eigenvectors r for the pair (A,B).
*
*  When DDRVGG 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, 7
*  tests will be performed and compared with the threshhold THRESH:
*
*  Results from DGEGS:
*
*                   T
*  (1)   | A - Q S Z  | / ( |A| n ulp )
*
*                   T
*  (2)   | B - Q T Z  | / ( |B| n ulp )
*
*                T
*  (3)   | I - QQ  | / ( n ulp )
*
*                T
*  (4)   | I - ZZ  | / ( n ulp )
*
*  (5)   maximum over j of D(j)  where:
*
*  if alpha(j) is real:
*                      |alpha(j) - S(j,j)|        |beta(j) - T(j,j)|
*            D(j) = ------------------------ + -----------------------
*                   max(|alpha(j)|,|S(j,j)|)   max(|beta(j)|,|T(j,j)|)
*
*  if alpha(j) is complex:
*                                  | det( s S - w T ) |
*            D(j) = ---------------------------------------------------
*                   ulp max( s norm(S), |w| norm(T) )*norm( s S - w T )
*
*            and S and T are here the 2 x 2 diagonal blocks of S and T
*            corresponding to the j-th eigenvalue.
*
*  Results from DGEGV:
*
*  (6)   max over all left eigenvalue/-vector pairs (beta/alpha,l) of
*
*     | l**H * (beta A - alpha B) | / ( ulp max( |beta A|, |alpha B| ) )
*
*        where l**H is the conjugate tranpose of l.
*
*  (7)   max over all right eigenvalue/-vector pairs (beta/alpha,r) of
*
*        | (beta A - alpha B) r | / ( ulp max( |beta A|, |alpha B| ) )
*
*  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 diag( 0, 1,..., N-1 ) (a diagonal
*                        matrix with those diagonal entries.)
*  (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 is diag( 0, 0, 1, ..., N-3, 0 ) and
*                         D2 is diag( 0, N-3, N-4,..., 1, 0, 0 )
*            t   t
*  (16) Q ( J , J ) Z     where Q and Z are random orthogonal matrices.
*
*  (17) Q ( T1, T2 ) Z    where T1 and T2 are upper triangular matrices
*                         with random O(1) entries above the diagonal
*                         and diagonal entries diag(T1) =
*                         ( 0, 0, 1, ..., N-3, 0 ) and diag(T2) =
*                         ( 0, N-3, N-4,..., 1, 0, 0 )
*
*  (18) Q ( T1, T2 ) Z    diag(T1) = ( 0, 0, 1, 1, s, ..., s, 0 )
*                         diag(T2) = ( 0, 1, 0, 1,..., 1, 0 )
*                         s = machine precision.
*
*  (19) Q ( T1, T2 ) Z    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) Q ( T1, T2 ) Z    diag(T1)=( 0, 0, 1, 1, a, ..., a   =s, 0 )
*                         diag(T2) = ( 0, 1, 0, 1, ..., 1, 0, 0 )
*
*  (21) Q ( T1, T2 ) Z    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) Q ( big*T1, small*T2 ) Z    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
*                                   diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (23) Q ( small*T1, big*T2 ) Z    diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
*                                   diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (24) Q ( small*T1, small*T2 ) Z  diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
*                                   diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (25) Q ( big*T1, big*T2 ) Z      diag(T1) = ( 0, 0, 1, ..., N-3, 0 )
*                                   diag(T2) = ( 0, 1, ..., 1, 0, 0 )
*
*  (26) Q ( T1, T2 ) Z     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,
*          DDRVGG 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, DDRVGG
*          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 DDRVGG to continue the same random number
*          sequence.
*
*  THRESH  (input) DOUBLE PRECISION
*          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.
*
*  THRSHN  (input) DOUBLE PRECISION
*          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) DOUBLE PRECISION 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, S, T, S2, and T2.
*          It must be at least 1 and at least max( NN ).
*
*  B       (input/workspace) DOUBLE PRECISION 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.
*
*  S       (workspace) DOUBLE PRECISION array, dimension (LDA, max(NN))
*          The Schur form matrix computed from A by DGEGS.  On exit, S
*          contains the Schur form matrix corresponding to the matrix
*          in A.
*
*  T       (workspace) DOUBLE PRECISION array, dimension (LDA, max(NN))
*          The upper triangular matrix computed from B by DGEGS.
*
*  S2      (workspace) DOUBLE PRECISION array, dimension (LDA, max(NN))
*          The matrix computed from A by DGEGV.  This will be the
*          Schur form of some matrix related to A, but will not, in
*          general, be the same as S.
*
*  T2      (workspace) DOUBLE PRECISION array, dimension (LDA, max(NN))
*          The matrix computed from B by DGEGV.  This will be the
*          Schur form of some matrix related to B, but will not, in
*          general, be the same as T.
*
*  Q       (workspace) DOUBLE PRECISION array, dimension (LDQ, max(NN))
*          The (left) orthogonal matrix computed by DGEGS.
*
*  LDQ     (input) INTEGER
*          The leading dimension of Q, Z, VL, and VR.  It must
*          be at least 1 and at least max( NN ).
*
*  Z       (workspace) DOUBLE PRECISION array of
*                             dimension( LDQ, max(NN) )
*          The (right) orthogonal matrix computed by DGEGS.
*
*  ALPHR1  (workspace) DOUBLE PRECISION array, dimension (max(NN))
*  ALPHI1  (workspace) DOUBLE PRECISION array, dimension (max(NN))
*  BETA1   (workspace) DOUBLE PRECISION array, dimension (max(NN))
*
*          The generalized eigenvalues of (A,B) computed by DGEGS.
*          ( ALPHR1(k)+ALPHI1(k)*i ) / BETA1(k) is the k-th
*          generalized eigenvalue of the matrices in A and B.
*
*  ALPHR2  (workspace) DOUBLE PRECISION array, dimension (max(NN))
*  ALPHI2  (workspace) DOUBLE PRECISION array, dimension (max(NN))
*  BETA2   (workspace) DOUBLE PRECISION array, dimension (max(NN))
*
*          The generalized eigenvalues of (A,B) computed by DGEGV.
*          ( ALPHR2(k)+ALPHI2(k)*i ) / BETA2(k) is the k-th
*          generalized eigenvalue of the matrices in A and B.
*
*  VL      (workspace) DOUBLE PRECISION array, dimension (LDQ, max(NN))
*          The (block lower triangular) left eigenvector matrix for
*          the matrices in A and B.  (See DTGEVC for the format.)
*
*  VR      (workspace) DOUBLE PRECISION array, dimension (LDQ, max(NN))
*          The (block upper triangular) right eigenvector matrix for
*          the matrices in A and B.  (See DTGEVC for the format.)
*
*  WORK    (workspace) DOUBLE PRECISION array, dimension (LWORK)
*
*  LWORK   (input) INTEGER
*          The number of entries in WORK.  This must be at least
*          2*N + MAX( 6*N, N*(NB+1), (k+1)*(2*k+N+1) ), where
*          "k" is the sum of the blocksize and number-of-shifts for
*          DHGEQZ, and NB is the greatest of the blocksizes for
*          DGEQRF, DORMQR, and DORGQR.  (The blocksizes and the
*          number-of-shifts are retrieved through calls to ILAENV.)
*
*  RESULT  (output) DOUBLE PRECISION 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 ..
      DOUBLE PRECISION   ZEROONE
      PARAMETER          ( ZERO = 0.0D0ONE = 1.0D0 )
      INTEGER            MAXTYP
      PARAMETER          ( MAXTYP = 26 )
*     ..
*     .. Local Scalars ..
      LOGICAL            BADNNILABAD
      INTEGER            I1IADDIINFOINJJCJRJSIZEJTYPE,
     $                   LWKOPTMTYPESNN1NBNBZNERRSNMATS,
     $                   NMAXNSNTESTNTESTT
      DOUBLE PRECISION   SAFMAXSAFMINTEMP1TEMP2ULPULPINV
*     ..
*     .. Local Arrays ..
      INTEGER            IASIGNMAXTYP ), IBSIGNMAXTYP ),
     $                   IOLDSD4 ), KADD6 ), KAMAGNMAXTYP ),
     $                   KATYPEMAXTYP ), KAZEROMAXTYP ),
     $                   KBMAGNMAXTYP ), KBTYPEMAXTYP ),
     $                   KBZEROMAXTYP ), KCLASSMAXTYP ),
     $                   KTRIANMAXTYP ), KZ16 ), KZ26 )
      DOUBLE PRECISION   DUMMA4 ), RMAGN03 )
*     ..
*     .. External Functions ..
      INTEGER            ILAENV
      DOUBLE PRECISION   DLAMCHDLARND
      EXTERNAL           ILAENVDLAMCHDLARND
*     ..
*     .. External Subroutines ..
      EXTERNAL           ALASVMDGEGSDGEGVDGET51DGET52DGET53,
     $                   DLABADDLACPYDLARFGDLASETDLATM4DORM2R,
     $                   XERBLA
*     ..
*     .. Intrinsic Functions ..
      INTRINSIC          ABSDBLEMAXMINSIGN
*     ..
*     .. 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               IASIGN / 6*0202*22*03*2023*0,
     $                   5*20 /
      DATA               IBSIGN / 7*022*02*22*02029*0 /
*     ..
*     .. Executable Statements ..
*
*     Check for errors
*
      INFO = 0
*
      BADNN = .FALSE.
      NMAX = 1
      DO 10 J = 1NSIZES
         NMAX = MAXNMAXNNJ ) )
         IFNNJ ).LT.0 )
     $      BADNN = .TRUE.
   10 CONTINUE
*
*     Maximum blocksize and shift -- we assume that blocksize and number
*     of shifts are monotone increasing functions of N.
*
      NB = MAX1ILAENV1'DGEQRF'' 'NMAXNMAX-1-1 ),
     $     ILAENV1'DORMQR''LT'NMAXNMAXNMAX-1 ),
     $     ILAENV1'DORGQR'' 'NMAXNMAXNMAX-1 ) )
      NBZ = ILAENV1'DHGEQZ''SII'NMAX1NMAX0 )
      NS = ILAENV4'DHGEQZ''SII'NMAX1NMAX0 )
      I1 = NBZ + NS
      LWKOPT = 2*NMAX + MAX6*NMAXNMAX*NB+1 ),
     $         ( 2*I1+NMAX+1 )*I1+1 ) )
*
*     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 IFLDQ.LE.1 .OR. LDQ.LT.NMAX ) THEN
         INFO = -19
      ELSE IFLWKOPT.GT.LWORK ) THEN
         INFO = -30
      END IF
*
      IFINFO.NE.0 ) THEN
         CALL XERBLA'DDRVGG'-INFO )
         RETURN
      END IF
*
*     Quick return if possible
*
      IFNSIZES.EQ.0 .OR. NTYPES.EQ.0 )
     $   RETURN
*
      SAFMIN = DLAMCH'Safe minimum' )
      ULP = DLAMCH'Epsilon' )*DLAMCH'Base' )
      SAFMIN = SAFMIN / ULP
      SAFMAX = ONE / SAFMIN
      CALL DLABADSAFMINSAFMAX )
      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 170 JSIZE = 1NSIZES
         N = NNJSIZE )
         N1 = MAX1N )
         RMAGN2 ) = SAFMAX*ULP / DBLEN1 )
         RMAGN3 ) = SAFMIN*ULPINV*N1
*
         IFNSIZES.NE.1 ) THEN
            MTYPES = MINMAXTYPNTYPES )
         ELSE
            MTYPES = MINMAXTYP+1NTYPES )
         END IF
*
         DO 160 JTYPE = 1MTYPES
            IF.NOT.DOTYPEJTYPE ) )
     $         GO TO 160
            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:
*
*           KZLASS: =1 means w/o rotation, =2 means w/ rotation,
*                   =3 means random.
*           KATYPE: the "type" to be passed to DLATM4 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.
*           IASIGN: 1 if the diagonal elements of A are to be
*                   multiplied by a random magnitude 1 number, =2 if
*                   randomly chosen diagonal blocks are to be rotated
*                   to form 2x2 blocks.
*           KBTYPE, KBZERO, KBMAGN, IBSIGN: 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 DLASET'Full'NNZEROZEROALDA )
               ELSE
                  IN = N
               END IF
               CALL DLATM4KATYPEJTYPE ), INKZ1KAZEROJTYPE ) ),
     $                      KZ2KAZEROJTYPE ) ), IASIGNJTYPE ),
     $                      RMAGNKAMAGNJTYPE ) ), ULP,
     $                      RMAGNKTRIANJTYPE )*KAMAGNJTYPE ) ), 2,
     $                      ISEEDALDA )
               IADD = KADDKAZEROJTYPE ) )
               IFIADD.GT.0 .AND. IADD.LE.N )
     $            AIADDIADD ) = ONE
*
*              Generate B (w/o rotation)
*
               IFABSKBTYPEJTYPE ) ).EQ.3 ) THEN
                  IN = 2*( ( N-1 ) / 2 ) + 1
                  IFIN.NE.N )
     $               CALL DLASET'Full'NNZEROZEROBLDA )
               ELSE
                  IN = N
               END IF
               CALL DLATM4KBTYPEJTYPE ), INKZ1KBZEROJTYPE ) ),
     $                      KZ2KBZEROJTYPE ) ), IBSIGNJTYPE ),
     $                      RMAGNKBMAGNJTYPE ) ), ONE,
     $                      RMAGNKTRIANJTYPE )*KBMAGNJTYPE ) ), 2,
     $                      ISEEDBLDA )
               IADD = KADDKBZEROJTYPE ) )
               IFIADD.NE.0 .AND. IADD.LE.N )
     $            BIADDIADD ) = ONE
*
               IFKCLASSJTYPE ).EQ.2 .AND. N.GT.0 ) THEN
*
*                 Include rotations
*
*                 Generate Q, Z as Householder transformations times
*                 a diagonal matrix.
*
                  DO 50 JC = 1N - 1
                     DO 40 JR = JCN
                        QJRJC ) = DLARND3ISEED )
                        ZJRJC ) = DLARND3ISEED )
   40                CONTINUE
                     CALL DLARFGN+1-JCQJCJC ), QJC+1JC ), 1,
     $                            WORKJC ) )
                     WORK2*N+JC ) = SIGNONEQJCJC ) )
                     QJCJC ) = ONE
                     CALL DLARFGN+1-JCZJCJC ), ZJC+1JC ), 1,
     $                            WORKN+JC ) )
                     WORK3*N+JC ) = SIGNONEZJCJC ) )
                     ZJCJC ) = ONE
   50             CONTINUE
                  QNN ) = ONE
                  WORKN ) = ZERO
                  WORK3*N ) = SIGNONEDLARND2ISEED ) )
                  ZNN ) = ONE
                  WORK2*N ) = ZERO
                  WORK4*N ) = SIGNONEDLARND2ISEED ) )
*
*                 Apply the diagonal matrices
*
                  DO 70 JC = 1N
                     DO 60 JR = 1N
                        AJRJC ) = WORK2*N+JR )*WORK3*N+JC )*
     $                                AJRJC )
                        BJRJC ) = WORK2*N+JR )*WORK3*N+JC )*
     $                                BJRJC )
   60                CONTINUE
   70             CONTINUE
                  CALL DORM2R'L''N'NNN-1QLDQWORKA,
     $                         LDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
                  CALL DORM2R'R''T'NNN-1ZLDQWORKN+1 ),
     $                         ALDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
                  CALL DORM2R'L''N'NNN-1QLDQWORKB,
     $                         LDAWORK2*N+1 ), IINFO )
                  IFIINFO.NE.0 )
     $               GO TO 100
                  CALL DORM2R'R''T'NNN-1ZLDQWORKN+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 ) )*
     $                             DLARND2ISEED )
                     BJRJC ) = RMAGNKBMAGNJTYPE ) )*
     $                             DLARND2ISEED )
   80             CONTINUE
   90          CONTINUE
            END IF
*
  100       CONTINUE
*
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'Generator', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               RETURN
            END IF
*
  110       CONTINUE
*
*           Call DGEGS to compute H, T, Q, Z, alpha, and beta.
*
            CALL DLACPY' 'NNALDASLDA )
            CALL DLACPY' 'NNBLDATLDA )
            NTEST = 1
            RESULT1 ) = ULPINV
*
            CALL DGEGS'V''V'NSLDATLDAALPHR1ALPHI1,
     $                  BETA1QLDQZLDQWORKLWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'DGEGS', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 140
            END IF
*
            NTEST = 4
*
*           Do tests 1--4
*
            CALL DGET511NALDASLDAQLDQZLDQWORK,
     $                   RESULT1 ) )
            CALL DGET511NBLDATLDAQLDQZLDQWORK,
     $                   RESULT2 ) )
            CALL DGET513NBLDATLDAQLDQQLDQWORK,
     $                   RESULT3 ) )
            CALL DGET513NBLDATLDAZLDQZLDQWORK,
     $                   RESULT4 ) )
*
*           Do test 5: compare eigenvalues with diagonals.
*           Also check Schur form of A.
*
            TEMP1 = ZERO
*
            DO 120 J = 1N
               ILABAD = .FALSE.
               IFALPHI1J ).EQ.ZERO ) THEN
                  TEMP2 = ( ABSALPHR1J )-SJJ ) ) /
     $                    MAXSAFMINABSALPHR1J ) ), ABSSJ,
     $                    J ) ) )+ABSBETA1J )-TJJ ) ) /
     $                    MAXSAFMINABSBETA1J ) ), ABSTJ,
     $                    J ) ) ) ) / ULP
                  IFJ.LT.N ) THEN
                     IFSJ+1J ).NE.ZERO )
     $                  ILABAD = .TRUE.
                  END IF
                  IFJ.GT.1 ) THEN
                     IFSJJ-1 ).NE.ZERO )
     $                  ILABAD = .TRUE.
                  END IF
               ELSE
                  IFALPHI1J ).GT.ZERO ) THEN
                     I1 = J
                  ELSE
                     I1 = J - 1
                  END IF
                  IFI1.LE.0 .OR. I1.GE.N ) THEN
                     ILABAD = .TRUE.
                  ELSE IFI1.LT.N-1 ) THEN
                     IFSI1+2I1+1 ).NE.ZERO )
     $                  ILABAD = .TRUE.
                  ELSE IFI1.GT.1 ) THEN
                     IFSI1I1-1 ).NE.ZERO )
     $                  ILABAD = .TRUE.
                  END IF
                  IF.NOT.ILABAD ) THEN
                     CALL DGET53SI1I1 ), LDATI1I1 ), LDA,
     $                            BETA1J ), ALPHR1J ), ALPHI1J ),
     $                            TEMP2IINFO )
                     IFIINFO.GE.3 ) THEN
                        WRITENOUNIT, FMT = 9997 )IINFO, J, N, JTYPE,
     $                     IOLDSD
                        INFO = ABSIINFO )
                     END IF
                  ELSE
                     TEMP2 = ULPINV
                  END IF
               END IF
               TEMP1 = MAXTEMP1TEMP2 )
               IFILABAD ) THEN
                  WRITENOUNIT, FMT = 9996 )J, N, JTYPE, IOLDSD
               END IF
  120       CONTINUE
            RESULT5 ) = TEMP1
*
*           Call DGEGV to compute S2, T2, VL, and VR, do tests.
*
*           Eigenvalues and Eigenvectors
*
            CALL DLACPY' 'NNALDAS2LDA )
            CALL DLACPY' 'NNBLDAT2LDA )
            NTEST = 6
            RESULT6 ) = ULPINV
*
            CALL DGEGV'V''V'NS2LDAT2LDAALPHR2ALPHI2,
     $                  BETA2VLLDQVRLDQWORKLWORKIINFO )
            IFIINFO.NE.0 ) THEN
               WRITENOUNIT, FMT = 9999 )'DGEGV', IINFO, N, JTYPE,
     $            IOLDSD
               INFO = ABSIINFO )
               GO TO 140
            END IF
*
            NTEST = 7
*
*           Do Tests 6 and 7
*
            CALL DGET52.TRUE.NALDABLDAVLLDQALPHR2,
     $                   ALPHI2BETA2WORKDUMMA1 ) )
            RESULT6 ) = DUMMA1 )
            IFDUMMA2 ).GT.THRSHN ) THEN
               WRITENOUNIT, FMT = 9998 )'Left''DGEGV', DUMMA( 2 ),
     $            N, JTYPE, IOLDSD
            END IF
*
            CALL DGET52.FALSE.NALDABLDAVRLDQALPHR2,
     $                   ALPHI2BETA2WORKDUMMA1 ) )
            RESULT7 ) = DUMMA1 )
            IFDUMMA2 ).GT.THRESH ) THEN
               WRITENOUNIT, FMT = 9998 )'Right''DGEGV', DUMMA( 2 ),
     $            N, JTYPE, IOLDSD
            END IF
*
*           Check form of Complex eigenvalues.
*
            DO 130 J = 1N
               ILABAD = .FALSE.
               IFALPHI2J ).GT.ZERO ) THEN
                  IFJ.EQ.N ) THEN
                     ILABAD = .TRUE.
                  ELSE IFALPHI2J+1 ).GE.ZERO ) THEN
                     ILABAD = .TRUE.
                  END IF
               ELSE IFALPHI2J ).LT.ZERO ) THEN
                  IFJ.EQ.1 ) THEN
                     ILABAD = .TRUE.
                  ELSE IFALPHI2J-1 ).LE.ZERO ) THEN
                     ILABAD = .TRUE.
                  END IF
               END IF
               IFILABAD ) THEN
                  WRITENOUNIT, FMT = 9996 )J, N, JTYPE, IOLDSD
               END IF
  130       CONTINUE
*
*           End of Loop -- Check for RESULT(j) > THRESH
*
  140       CONTINUE
*
            NTESTT = NTESTT + NTEST
*
*           Print out tests which fail.
*
            DO 150 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 = 9995 )'DGG'
*
*                    Matrix types
*
                     WRITENOUNIT, FMT = 9994 )
                     WRITENOUNIT, FMT = 9993 )
                     WRITENOUNIT, FMT = 9992 )'Orthogonal'
*
*                    Tests performed
*
                     WRITENOUNIT, FMT = 9991 )'orthogonal''''',
     $                  'transpose', ( '''', J = 15 )
*
                  END IF
                  NERRS = NERRS + 1
                  IFRESULTJR ).LT.10000.0D0 ) THEN
                     WRITENOUNIT, FMT = 9990 )N, JTYPE, IOLDSD, JR,
     $                  RESULT( JR )
                  ELSE
                     WRITENOUNIT, FMT = 9989 )N, JTYPE, IOLDSD, JR,
     $                  RESULT( JR )
                  END IF
               END IF
  150       CONTINUE
*
  160    CONTINUE
  170 CONTINUE
*
*     Summary
*
      CALL ALASVM'DGG'NOUNITNERRSNTESTT0 )
      RETURN
*
 9999 FORMAT( ' DDRVGG: ', A, ' returned INFO=', I6, '.', / 9X, 'N=',
     $      I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5, ')' )
*
 9998 FORMAT( ' DDRVGG: ', A, ' Eigenvectors from ', A, ' incorrectly ',
     $      'normalized.', / ' Bits of error=', 0P, G10.3, ',', 9X,
     $      'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ), I5,
     $      ')' )
*
 9997 FORMAT( ' DDRVGG: DGET53 returned INFO=', I1, ' for eigenvalue ',
     $      I6, '.', / 9X, 'N=', I6, ', JTYPE=', I6, ', ISEED=(',
     $      3( I5, ',' ), I5, ')' )
*
 9996 FORMAT( ' DDRVGG: S not in Schur form at eigenvalue ', I6, '.',
     $      / 9X, 'N=', I6, ', JTYPE=', I6, ', ISEED=(', 3( I5, ',' ),
     $      I5, ')' )
*
 9995 FORMAT( / 1X, A3, ' -- Real Generalized eigenvalue problem driver'
     $       )
*
 9994 FORMAT( ' Matrix types (see DDRVGG for details): ' )
*
 9993 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)' )
 9992 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.' )
*
 9991 FORMAT( / ' Tests performed:  (S is Schur, T is triangular, ',
     $      'Q and Z are ', A, ',', / 20X,
     $      'l and r are the appropriate left and right', / 19X,
     $      'eigenvectors, resp., a is alpha, b is beta, and', / 19X, A,
     $      ' means ', A, '.)', / ' 1 = | A - Q S Z', A,
     $      ' | / ( |A| n ulp )      2 = | B - Q T Z', A,
     $      ' | / ( |B| n ulp )', / ' 3 = | I - QQ', A,
     $      ' | / ( n ulp )             4 = | I - ZZ', A,
     $      ' | / ( n ulp )', /
     $      ' 5 = difference between (alpha,beta) and diagonals of',
     $      ' (S,T)', / ' 6 = max | ( b A - a B )', A,
     $      ' l | / const.   7 = max | ( b A - a B ) r | / const.',
     $      / 1X )
 9990 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
     $      4( I4, ',' ), ' result ', I3, ' is', 0P, F8.2 )
 9989 FORMAT( ' Matrix order=', I5, ', type=', I2, ', seed=',
     $      4( I4, ',' ), ' result ', I3, ' is', 1P, D10.3 )
*
*     End of DDRVGG
*
      END