1*> \brief \b ZCHKHE_ROOK
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 ZCHKHE_ROOK( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
12*                               THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
13*                               XACT, WORK, RWORK, IWORK, NOUT )
14*
15*       .. Scalar Arguments ..
16*       LOGICAL            TSTERR
17*       INTEGER            NMAX, NN, NNB, NNS, NOUT
18*       DOUBLE PRECISION   THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
23*       DOUBLE PRECISION   RWORK( * )
24*       COMPLEX*16         A( * ), AFAC( * ), AINV( * ), B( * ),
25*      $                   WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> ZCHKHE_ROOK tests ZHETRF_ROOK, -TRI_ROOK, -TRS_ROOK,
35*> and -CON_ROOK.
36*> \endverbatim
37*
38*  Arguments:
39*  ==========
40*
41*> \param[in] DOTYPE
42*> \verbatim
43*>          DOTYPE is LOGICAL array, dimension (NTYPES)
44*>          The matrix types to be used for testing.  Matrices of type j
45*>          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
46*>          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
47*> \endverbatim
48*>
49*> \param[in] NN
50*> \verbatim
51*>          NN is INTEGER
52*>          The number of values of N contained in the vector NVAL.
53*> \endverbatim
54*>
55*> \param[in] NVAL
56*> \verbatim
57*>          NVAL is INTEGER array, dimension (NN)
58*>          The values of the matrix dimension N.
59*> \endverbatim
60*>
61*> \param[in] NNB
62*> \verbatim
63*>          NNB is INTEGER
64*>          The number of values of NB contained in the vector NBVAL.
65*> \endverbatim
66*>
67*> \param[in] NBVAL
68*> \verbatim
69*>          NBVAL is INTEGER array, dimension (NBVAL)
70*>          The values of the blocksize NB.
71*> \endverbatim
72*>
73*> \param[in] NNS
74*> \verbatim
75*>          NNS is INTEGER
76*>          The number of values of NRHS contained in the vector NSVAL.
77*> \endverbatim
78*>
79*> \param[in] NSVAL
80*> \verbatim
81*>          NSVAL is INTEGER array, dimension (NNS)
82*>          The values of the number of right hand sides NRHS.
83*> \endverbatim
84*>
85*> \param[in] THRESH
86*> \verbatim
87*>          THRESH is DOUBLE PRECISION
88*>          The threshold value for the test ratios.  A result is
89*>          included in the output file if RESULT >= THRESH.  To have
90*>          every test ratio printed, use THRESH = 0.
91*> \endverbatim
92*>
93*> \param[in] TSTERR
94*> \verbatim
95*>          TSTERR is LOGICAL
96*>          Flag that indicates whether error exits are to be tested.
97*> \endverbatim
98*>
99*> \param[in] NMAX
100*> \verbatim
101*>          NMAX is INTEGER
102*>          The maximum value permitted for N, used in dimensioning the
103*>          work arrays.
104*> \endverbatim
105*>
106*> \param[out] A
107*> \verbatim
108*>          A is CCOMPLEX*16 array, dimension (NMAX*NMAX)
109*> \endverbatim
110*>
111*> \param[out] AFAC
112*> \verbatim
113*>          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
114*> \endverbatim
115*>
116*> \param[out] AINV
117*> \verbatim
118*>          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
119*> \endverbatim
120*>
121*> \param[out] B
122*> \verbatim
123*>          B is CCOMPLEX*16 array, dimension (NMAX*NSMAX)
124*>          where NSMAX is the largest entry in NSVAL.
125*> \endverbatim
126*>
127*> \param[out] X
128*> \verbatim
129*>          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
130*> \endverbatim
131*>
132*> \param[out] XACT
133*> \verbatim
134*>          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
135*> \endverbatim
136*>
137*> \param[out] WORK
138*> \verbatim
139*>          WORK is COMPLEX*16 array, dimension (NMAX*max(3,NSMAX))
140*> \endverbatim
141*>
142*> \param[out] RWORK
143*> \verbatim
144*>          RWORK is DOUBLE PRECISION array, dimension (max(NMAX,2*NSMAX)
145*> \endverbatim
146*>
147*> \param[out] IWORK
148*> \verbatim
149*>          IWORK is INTEGER array, dimension (2*NMAX)
150*> \endverbatim
151*>
152*> \param[in] NOUT
153*> \verbatim
154*>          NOUT is INTEGER
155*>          The unit number for output.
156*> \endverbatim
157*
158*  Authors:
159*  ========
160*
161*> \author Univ. of Tennessee
162*> \author Univ. of California Berkeley
163*> \author Univ. of Colorado Denver
164*> \author NAG Ltd.
165*
166*> \date November 2013
167*
168*> \ingroup complex16_lin
169*
170*  =====================================================================
171      SUBROUTINE ZCHKHE_ROOK( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
172     $                        THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
173     $                        XACT, WORK, RWORK, IWORK, NOUT )
174*
175*  -- LAPACK test routine (version 3.5.0) --
176*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
177*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178*     November 2013
179*
180*     .. Scalar Arguments ..
181      LOGICAL            TSTERR
182      INTEGER            NMAX, NN, NNB, NNS, NOUT
183      DOUBLE PRECISION   THRESH
184*     ..
185*     .. Array Arguments ..
186      LOGICAL            DOTYPE( * )
187      INTEGER            IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
188      DOUBLE PRECISION   RWORK( * )
189      COMPLEX*16         A( * ), AFAC( * ), AINV( * ), B( * ),
190     $                   WORK( * ), X( * ), XACT( * )
191*     ..
192*
193*  =====================================================================
194*
195*     .. Parameters ..
196      DOUBLE PRECISION   ZERO, ONE
197      PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
198      DOUBLE PRECISION   ONEHALF
199      PARAMETER          ( ONEHALF = 0.5D+0 )
200      DOUBLE PRECISION   EIGHT, SEVTEN
201      PARAMETER          ( EIGHT = 8.0D+0, SEVTEN = 17.0D+0 )
202      COMPLEX*16         CZERO
203      PARAMETER          ( CZERO = ( 0.0D+0, 0.0D+0 ) )
204      INTEGER            NTYPES
205      PARAMETER          ( NTYPES = 10 )
206      INTEGER            NTESTS
207      PARAMETER          ( NTESTS = 7 )
208*     ..
209*     .. Local Scalars ..
210      LOGICAL            TRFCON, ZEROT
211      CHARACTER          DIST, TYPE, UPLO, XTYPE
212      CHARACTER*3        PATH, MATPATH
213      INTEGER            I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
214     $                   ITEMP, ITEMP2, IUPLO, IZERO, J, K, KL, KU, LDA,
215     $                   LWORK, MODE, N, NB, NERRS, NFAIL, NIMAT, NRHS,
216     $                   NRUN, NT
217      DOUBLE PRECISION   ALPHA, ANORM, CNDNUM, CONST, LAM_MAX, LAM_MIN,
218     $                   RCOND, RCONDC, DTEMP
219*     ..
220*     .. Local Arrays ..
221      CHARACTER          UPLOS( 2 )
222      INTEGER            ISEED( 4 ), ISEEDY( 4 ), IDUMMY( 1 )
223      DOUBLE PRECISION   RESULT( NTESTS )
224      COMPLEX*16         CDUMMY( 1 )
225*     ..
226*     .. External Functions ..
227      DOUBLE PRECISION   ZLANGE, ZLANHE, DGET06
228      EXTERNAL           ZLANGE, ZLANHE, DGET06
229*     ..
230*     .. External Subroutines ..
231      EXTERNAL           ALAERH, ALAHD, ALASUM, ZERRHE, ZHEEVX, ZGET04,
232     $                   ZLACPY, ZLARHS, ZLATB4, ZLATMS, ZPOT02,
233     $                   ZPOT03, ZHECON_ROOK, ZHET01_ROOK, ZHETRF_ROOK,
234     $                   ZHETRI_ROOK, ZHETRS_ROOK, XLAENV
235*     ..
236*     .. Intrinsic Functions ..
237      INTRINSIC          ABS, MAX, MIN, SQRT
238*     ..
239*     .. Scalars in Common ..
240      LOGICAL            LERR, OK
241      CHARACTER*32       SRNAMT
242      INTEGER            INFOT, NUNIT
243*     ..
244*     .. Common blocks ..
245      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
246      COMMON             / SRNAMC / SRNAMT
247*     ..
248*     .. Data statements ..
249      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
250      DATA               UPLOS / 'U', 'L' /
251*     ..
252*     .. Executable Statements ..
253*
254*     Initialize constants and the random number seed.
255*
256      ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT
257*
258*     Test path
259*
260      PATH( 1: 1 ) = 'Zomplex precision'
261      PATH( 2: 3 ) = 'HR'
262*
263*     Path to generate matrices
264*
265      MATPATH( 1: 1 ) = 'Zomplex precision'
266      MATPATH( 2: 3 ) = 'HE'
267*
268      NRUN = 0
269      NFAIL = 0
270      NERRS = 0
271      DO 10 I = 1, 4
272         ISEED( I ) = ISEEDY( I )
273   10 CONTINUE
274*
275*     Test the error exits
276*
277      IF( TSTERR )
278     $   CALL ZERRHE( PATH, NOUT )
279      INFOT = 0
280*
281*     Set the minimum block size for which the block routine should
282*     be used, which will be later returned by ILAENV
283*
284      CALL XLAENV( 2, 2 )
285*
286*     Do for each value of N in NVAL
287*
288      DO 270 IN = 1, NN
289         N = NVAL( IN )
290         LDA = MAX( N, 1 )
291         XTYPE = 'N'
292         NIMAT = NTYPES
293         IF( N.LE.0 )
294     $      NIMAT = 1
295*
296         IZERO = 0
297*
298*        Do for each value of matrix type IMAT
299*
300         DO 260 IMAT = 1, NIMAT
301*
302*           Do the tests only if DOTYPE( IMAT ) is true.
303*
304            IF( .NOT.DOTYPE( IMAT ) )
305     $         GO TO 260
306*
307*           Skip types 3, 4, 5, or 6 if the matrix size is too small.
308*
309            ZEROT = IMAT.GE.3 .AND. IMAT.LE.6
310            IF( ZEROT .AND. N.LT.IMAT-2 )
311     $         GO TO 260
312*
313*           Do first for UPLO = 'U', then for UPLO = 'L'
314*
315            DO 250 IUPLO = 1, 2
316               UPLO = UPLOS( IUPLO )
317*
318*                 Begin generate the test matrix A.
319*
320*                 Set up parameters with ZLATB4 for the matrix generator
321*                 based on the type of matrix to be generated.
322*
323                  CALL ZLATB4( MATPATH, IMAT, N, N, TYPE, KL, KU, ANORM,
324     $                         MODE, CNDNUM, DIST )
325*
326*                 Generate a matrix with ZLATMS.
327*
328                  SRNAMT = 'ZLATMS'
329                  CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
330     $                         CNDNUM, ANORM, KL, KU, UPLO, A, LDA,
331     $                         WORK, INFO )
332*
333*                 Check error code from ZLATMS and handle error.
334*
335                  IF( INFO.NE.0 ) THEN
336                     CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N,
337     $                            -1, -1, -1, IMAT, NFAIL, NERRS, NOUT )
338*
339*                    Skip all tests for this generated matrix
340*
341                     GO TO 250
342                  END IF
343*
344*                 For matrix types 3-6, zero one or more rows and
345*                 columns of the matrix to test that INFO is returned
346*                 correctly.
347*
348                  IF( ZEROT ) THEN
349                     IF( IMAT.EQ.3 ) THEN
350                        IZERO = 1
351                     ELSE IF( IMAT.EQ.4 ) THEN
352                        IZERO = N
353                     ELSE
354                        IZERO = N / 2 + 1
355                     END IF
356*
357                     IF( IMAT.LT.6 ) THEN
358*
359*                       Set row and column IZERO to zero.
360*
361                        IF( IUPLO.EQ.1 ) THEN
362                           IOFF = ( IZERO-1 )*LDA
363                           DO 20 I = 1, IZERO - 1
364                              A( IOFF+I ) = CZERO
365   20                      CONTINUE
366                           IOFF = IOFF + IZERO
367                           DO 30 I = IZERO, N
368                              A( IOFF ) = CZERO
369                              IOFF = IOFF + LDA
370   30                      CONTINUE
371                        ELSE
372                           IOFF = IZERO
373                           DO 40 I = 1, IZERO - 1
374                              A( IOFF ) = CZERO
375                              IOFF = IOFF + LDA
376   40                      CONTINUE
377                           IOFF = IOFF - IZERO
378                           DO 50 I = IZERO, N
379                              A( IOFF+I ) = CZERO
380   50                      CONTINUE
381                        END IF
382                     ELSE
383                        IF( IUPLO.EQ.1 ) THEN
384*
385*                          Set the first IZERO rows and columns to zero.
386*
387                           IOFF = 0
388                           DO 70 J = 1, N
389                              I2 = MIN( J, IZERO )
390                              DO 60 I = 1, I2
391                                 A( IOFF+I ) = CZERO
392   60                         CONTINUE
393                              IOFF = IOFF + LDA
394   70                      CONTINUE
395                        ELSE
396*
397*                          Set the last IZERO rows and columns to zero.
398*
399                           IOFF = 0
400                           DO 90 J = 1, N
401                              I1 = MAX( J, IZERO )
402                              DO 80 I = I1, N
403                                 A( IOFF+I ) = CZERO
404   80                         CONTINUE
405                              IOFF = IOFF + LDA
406   90                      CONTINUE
407                        END IF
408                     END IF
409                  ELSE
410                     IZERO = 0
411                  END IF
412*
413*                 End generate the test matrix A.
414*
415*
416*              Do for each value of NB in NBVAL
417*
418               DO 240 INB = 1, NNB
419*
420*                 Set the optimal blocksize, which will be later
421*                 returned by ILAENV.
422*
423                  NB = NBVAL( INB )
424                  CALL XLAENV( 1, NB )
425*
426*                 Copy the test matrix A into matrix AFAC which
427*                 will be factorized in place. This is needed to
428*                 preserve the test matrix A for subsequent tests.
429*
430                  CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
431*
432*                 Compute the L*D*L**T or U*D*U**T factorization of the
433*                 matrix. IWORK stores details of the interchanges and
434*                 the block structure of D. AINV is a work array for
435*                 block factorization, LWORK is the length of AINV.
436*
437                  LWORK = MAX( 2, NB )*LDA
438                  SRNAMT = 'ZHETRF_ROOK'
439                  CALL ZHETRF_ROOK( UPLO, N, AFAC, LDA, IWORK, AINV,
440     $                              LWORK, INFO )
441*
442*                 Adjust the expected value of INFO to account for
443*                 pivoting.
444*
445                  K = IZERO
446                  IF( K.GT.0 ) THEN
447  100                CONTINUE
448                     IF( IWORK( K ).LT.0 ) THEN
449                        IF( IWORK( K ).NE.-K ) THEN
450                           K = -IWORK( K )
451                           GO TO 100
452                        END IF
453                     ELSE IF( IWORK( K ).NE.K ) THEN
454                        K = IWORK( K )
455                        GO TO 100
456                     END IF
457                  END IF
458*
459*                 Check error code from ZHETRF_ROOK and handle error.
460*
461                  IF( INFO.NE.K)
462     $               CALL ALAERH( PATH, 'ZHETRF_ROOK', INFO, K,
463     $                            UPLO, N, N, -1, -1, NB, IMAT,
464     $                            NFAIL, NERRS, NOUT )
465*
466*                 Set the condition estimate flag if the INFO is not 0.
467*
468                  IF( INFO.NE.0 ) THEN
469                     TRFCON = .TRUE.
470                  ELSE
471                     TRFCON = .FALSE.
472                  END IF
473*
474*+    TEST 1
475*                 Reconstruct matrix from factors and compute residual.
476*
477                  CALL ZHET01_ROOK( UPLO, N, A, LDA, AFAC, LDA, IWORK,
478     $                              AINV, LDA, RWORK, RESULT( 1 ) )
479                  NT = 1
480*
481*+    TEST 2
482*                 Form the inverse and compute the residual,
483*                 if the factorization was competed without INFO > 0
484*                 (i.e. there is no zero rows and columns).
485*                 Do it only for the first block size.
486*
487                  IF( INB.EQ.1 .AND. .NOT.TRFCON ) THEN
488                     CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
489                     SRNAMT = 'ZHETRI_ROOK'
490                     CALL ZHETRI_ROOK( UPLO, N, AINV, LDA, IWORK, WORK,
491     $                                 INFO )
492*
493*                    Check error code from ZHETRI_ROOK and handle error.
494*
495                     IF( INFO.NE.0 )
496     $                  CALL ALAERH( PATH, 'ZHETRI_ROOK', INFO, -1,
497     $                               UPLO, N, N, -1, -1, -1, IMAT,
498     $                               NFAIL, NERRS, NOUT )
499*
500*                    Compute the residual for a Hermitian matrix times
501*                    its inverse.
502*
503                     CALL ZPOT03( UPLO, N, A, LDA, AINV, LDA, WORK, LDA,
504     $                            RWORK, RCONDC, RESULT( 2 ) )
505                     NT = 2
506                  END IF
507*
508*                 Print information about the tests that did not pass
509*                 the threshold.
510*
511                  DO 110 K = 1, NT
512                     IF( RESULT( K ).GE.THRESH ) THEN
513                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
514     $                     CALL ALAHD( NOUT, PATH )
515                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
516     $                     RESULT( K )
517                        NFAIL = NFAIL + 1
518                     END IF
519  110             CONTINUE
520                  NRUN = NRUN + NT
521*
522*+    TEST 3
523*                 Compute largest element in U or L
524*
525                  RESULT( 3 ) = ZERO
526                  DTEMP = ZERO
527*
528                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) ) /
529     $                    ( ONE-ALPHA )
530*
531                  IF( IUPLO.EQ.1 ) THEN
532*
533*                 Compute largest element in U
534*
535                     K = N
536  120                CONTINUE
537                     IF( K.LE.1 )
538     $                  GO TO 130
539*
540                     IF( IWORK( K ).GT.ZERO ) THEN
541*
542*                       Get max absolute value from elements
543*                       in column k in U
544*
545                        DTEMP = ZLANGE( 'M', K-1, 1,
546     $                          AFAC( ( K-1 )*LDA+1 ), LDA, RWORK )
547                     ELSE
548*
549*                       Get max absolute value from elements
550*                       in columns k and k-1 in U
551*
552                        DTEMP = ZLANGE( 'M', K-2, 2,
553     $                          AFAC( ( K-2 )*LDA+1 ), LDA, RWORK )
554                        K = K - 1
555*
556                     END IF
557*
558*                    DTEMP should be bounded by CONST
559*
560                     DTEMP = DTEMP - CONST + THRESH
561                     IF( DTEMP.GT.RESULT( 3 ) )
562     $                  RESULT( 3 ) = DTEMP
563*
564                     K = K - 1
565*
566                     GO TO 120
567  130                CONTINUE
568*
569                  ELSE
570*
571*                 Compute largest element in L
572*
573                     K = 1
574  140                CONTINUE
575                     IF( K.GE.N )
576     $                  GO TO 150
577*
578                     IF( IWORK( K ).GT.ZERO ) THEN
579*
580*                       Get max absolute value from elements
581*                       in column k in L
582*
583                        DTEMP = ZLANGE( 'M', N-K, 1,
584     $                          AFAC( ( K-1 )*LDA+K+1 ), LDA, RWORK )
585                     ELSE
586*
587*                       Get max absolute value from elements
588*                       in columns k and k+1 in L
589*
590                        DTEMP = ZLANGE( 'M', N-K-1, 2,
591     $                          AFAC( ( K-1 )*LDA+K+2 ), LDA, RWORK )
592                        K = K + 1
593*
594                     END IF
595*
596*                    DTEMP should be bounded by CONST
597*
598                     DTEMP = DTEMP - CONST + THRESH
599                     IF( DTEMP.GT.RESULT( 3 ) )
600     $                  RESULT( 3 ) = DTEMP
601*
602                     K = K + 1
603*
604                     GO TO 140
605  150                CONTINUE
606                  END IF
607*
608*
609*+    TEST 4
610*                 Compute largest 2-Norm of 2-by-2 diag blocks
611*
612                  RESULT( 4 ) = ZERO
613                  DTEMP = ZERO
614*
615                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) )*
616     $                    ( ( ONE + ALPHA ) / ( ONE - ALPHA ) )
617                  CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
618*
619                  IF( IUPLO.EQ.1 ) THEN
620*
621*                    Loop backward for UPLO = 'U'
622*
623                     K = N
624  160                CONTINUE
625                     IF( K.LE.1 )
626     $                  GO TO 170
627*
628                     IF( IWORK( K ).LT.ZERO ) THEN
629*
630*                       Get the two eigenvalues of a 2-by-2 block,
631*                       store them in WORK array
632*
633                        CALL ZHEEVX( 'N', 'A', UPLO, 2,
634     $                               AINV( ( K-2 )*LDA+K-1 ), LDA,DTEMP,
635     $                               DTEMP, ITEMP, ITEMP, ZERO, ITEMP,
636     $                               RWORK, CDUMMY, 1, WORK, 16,
637     $                               RWORK( 3 ), IWORK( N+1 ), IDUMMY,
638     $                               INFO )
639*
640                        LAM_MAX = MAX( ABS( RWORK( 1 ) ),
641     $                            ABS( RWORK( 2 ) ) )
642                        LAM_MIN = MIN( ABS( RWORK( 1 ) ),
643     $                            ABS( RWORK( 2 ) ) )
644*
645                        DTEMP = LAM_MAX / LAM_MIN
646*
647*                       DTEMP should be bounded by CONST
648*
649                        DTEMP = ABS( DTEMP ) - CONST + THRESH
650                        IF( DTEMP.GT.RESULT( 4 ) )
651     $                     RESULT( 4 ) = DTEMP
652                        K = K - 1
653*
654                     END IF
655*
656                     K = K - 1
657*
658                     GO TO 160
659  170                CONTINUE
660*
661                  ELSE
662*
663*                    Loop forward for UPLO = 'L'
664*
665                     K = 1
666  180                CONTINUE
667                     IF( K.GE.N )
668     $                  GO TO 190
669*
670                     IF( IWORK( K ).LT.ZERO ) THEN
671*
672*                       Get the two eigenvalues of a 2-by-2 block,
673*                       store them in WORK array
674*
675                        CALL ZHEEVX( 'N', 'A', UPLO, 2,
676     $                               AINV( ( K-1 )*LDA+K ), LDA, DTEMP,
677     $                               DTEMP, ITEMP, ITEMP, ZERO, ITEMP,
678     $                               RWORK, CDUMMY, 1, WORK, 16,
679     $                               RWORK( 3 ), IWORK( N+1 ), IDUMMY,
680     $                               INFO )
681*
682                        LAM_MAX = MAX( ABS( RWORK( 1 ) ),
683     $                            ABS( RWORK( 2 ) ) )
684                        LAM_MIN = MIN( ABS( RWORK( 1 ) ),
685     $                            ABS( RWORK( 2 ) ) )
686*
687                        DTEMP = LAM_MAX / LAM_MIN
688*
689*                       DTEMP should be bounded by CONST
690*
691                        DTEMP = ABS( DTEMP ) - CONST + THRESH
692                        IF( DTEMP.GT.RESULT( 4 ) )
693     $                     RESULT( 4 ) = DTEMP
694                        K = K + 1
695*
696                     END IF
697*
698                     K = K + 1
699*
700                     GO TO 180
701  190                CONTINUE
702                  END IF
703*
704*                 Print information about the tests that did not pass
705*                 the threshold.
706*
707                  DO 200 K = 3, 4
708                     IF( RESULT( K ).GE.THRESH ) THEN
709                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
710     $                     CALL ALAHD( NOUT, PATH )
711                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
712     $                     RESULT( K )
713                        NFAIL = NFAIL + 1
714                     END IF
715  200             CONTINUE
716                  NRUN = NRUN + 2
717*
718*                 Skip the other tests if this is not the first block
719*                 size.
720*
721                  IF( INB.GT.1 )
722     $               GO TO 240
723*
724*                 Do only the condition estimate if INFO is not 0.
725*
726                  IF( TRFCON ) THEN
727                     RCONDC = ZERO
728                     GO TO 230
729                  END IF
730*
731*                 Do for each value of NRHS in NSVAL.
732*
733                  DO 220 IRHS = 1, NNS
734                     NRHS = NSVAL( IRHS )
735*
736*                    Begin loop over NRHS values
737*
738*
739*+    TEST 5 ( Using TRS_ROOK)
740*                 Solve and compute residual for  A * X = B.
741*
742*                    Choose a set of NRHS random solution vectors
743*                    stored in XACT and set up the right hand side B
744*
745                     SRNAMT = 'ZLARHS'
746                     CALL ZLARHS( MATPATH, XTYPE, UPLO, ' ', N, N,
747     $                            KL, KU, NRHS, A, LDA, XACT, LDA,
748     $                            B, LDA, ISEED, INFO )
749                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
750*
751                     SRNAMT = 'ZHETRS_ROOK'
752                     CALL ZHETRS_ROOK( UPLO, N, NRHS, AFAC, LDA, IWORK,
753     $                                 X, LDA, INFO )
754*
755*                    Check error code from ZHETRS_ROOK and handle error.
756*
757                     IF( INFO.NE.0 )
758     $                  CALL ALAERH( PATH, 'ZHETRS_ROOK', INFO, 0,
759     $                               UPLO, N, N, -1, -1, NRHS, IMAT,
760     $                               NFAIL, NERRS, NOUT )
761*
762                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK, LDA )
763*
764*                    Compute the residual for the solution
765*
766                     CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA, WORK,
767     $                            LDA, RWORK, RESULT( 5 ) )
768*
769*+    TEST 6
770*                 Check solution from generated exact solution.
771*
772                     CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
773     $                            RESULT( 6 ) )
774*
775*                    Print information about the tests that did not pass
776*                    the threshold.
777*
778                     DO 210 K = 5, 6
779                        IF( RESULT( K ).GE.THRESH ) THEN
780                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
781     $                        CALL ALAHD( NOUT, PATH )
782                           WRITE( NOUT, FMT = 9998 )UPLO, N, NRHS,
783     $                        IMAT, K, RESULT( K )
784                           NFAIL = NFAIL + 1
785                        END IF
786  210                CONTINUE
787                     NRUN = NRUN + 2
788*
789*                 End do for each value of NRHS in NSVAL.
790*
791  220             CONTINUE
792*
793*+    TEST 7
794*                 Get an estimate of RCOND = 1/CNDNUM.
795*
796  230             CONTINUE
797                  ANORM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
798                  SRNAMT = 'ZHECON_ROOK'
799                  CALL ZHECON_ROOK( UPLO, N, AFAC, LDA, IWORK, ANORM,
800     $                              RCOND, WORK, INFO )
801*
802*                 Check error code from ZHECON_ROOK and handle error.
803*
804                  IF( INFO.NE.0 )
805     $               CALL ALAERH( PATH, 'ZHECON_ROOK', INFO, 0,
806     $                             UPLO, N, N, -1, -1, -1, IMAT,
807     $                             NFAIL, NERRS, NOUT )
808*
809*                 Compute the test ratio to compare values of RCOND
810*
811                  RESULT( 7 ) = DGET06( RCOND, RCONDC )
812*
813*                 Print information about the tests that did not pass
814*                 the threshold.
815*
816                  IF( RESULT( 7 ).GE.THRESH ) THEN
817                     IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
818     $                  CALL ALAHD( NOUT, PATH )
819                     WRITE( NOUT, FMT = 9997 )UPLO, N, IMAT, 7,
820     $                  RESULT( 7 )
821                     NFAIL = NFAIL + 1
822                  END IF
823                  NRUN = NRUN + 1
824  240          CONTINUE
825*
826  250       CONTINUE
827  260    CONTINUE
828  270 CONTINUE
829*
830*     Print a summary of the results.
831*
832      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
833*
834 9999 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NB =', I4, ', type ',
835     $      I2, ', test ', I2, ', ratio =', G12.5 )
836 9998 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NRHS=', I3, ', type ',
837     $      I2, ', test ', I2, ', ratio =', G12.5 )
838 9997 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ',', 10X, ' type ', I2,
839     $      ', test ', I2, ', ratio =', G12.5 )
840      RETURN
841*
842*     End of ZCHKHE_ROOK
843*
844      END
845