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 (NNB)
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*> \ingroup complex16_lin
167*
168*  =====================================================================
169      SUBROUTINE ZCHKHE_ROOK( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
170     $                        THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
171     $                        XACT, WORK, RWORK, IWORK, NOUT )
172*
173*  -- LAPACK test routine --
174*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
175*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
176*
177*     .. Scalar Arguments ..
178      LOGICAL            TSTERR
179      INTEGER            NMAX, NN, NNB, NNS, NOUT
180      DOUBLE PRECISION   THRESH
181*     ..
182*     .. Array Arguments ..
183      LOGICAL            DOTYPE( * )
184      INTEGER            IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
185      DOUBLE PRECISION   RWORK( * )
186      COMPLEX*16         A( * ), AFAC( * ), AINV( * ), B( * ),
187     $                   WORK( * ), X( * ), XACT( * )
188*     ..
189*
190*  =====================================================================
191*
192*     .. Parameters ..
193      DOUBLE PRECISION   ZERO, ONE
194      PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
195      DOUBLE PRECISION   ONEHALF
196      PARAMETER          ( ONEHALF = 0.5D+0 )
197      DOUBLE PRECISION   EIGHT, SEVTEN
198      PARAMETER          ( EIGHT = 8.0D+0, SEVTEN = 17.0D+0 )
199      COMPLEX*16         CZERO
200      PARAMETER          ( CZERO = ( 0.0D+0, 0.0D+0 ) )
201      INTEGER            NTYPES
202      PARAMETER          ( NTYPES = 10 )
203      INTEGER            NTESTS
204      PARAMETER          ( NTESTS = 7 )
205*     ..
206*     .. Local Scalars ..
207      LOGICAL            TRFCON, ZEROT
208      CHARACTER          DIST, TYPE, UPLO, XTYPE
209      CHARACTER*3        PATH, MATPATH
210      INTEGER            I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
211     $                   IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
212     $                   N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
213      DOUBLE PRECISION   ALPHA, ANORM, CNDNUM, CONST, SING_MAX,
214     $                   SING_MIN, RCOND, RCONDC, DTEMP
215*     ..
216*     .. Local Arrays ..
217      CHARACTER          UPLOS( 2 )
218      INTEGER            ISEED( 4 ), ISEEDY( 4 )
219      DOUBLE PRECISION   RESULT( NTESTS )
220      COMPLEX*16         BLOCK( 2, 2 ), ZDUMMY( 1 )
221*     ..
222*     .. External Functions ..
223      DOUBLE PRECISION   ZLANGE, ZLANHE, DGET06
224      EXTERNAL           ZLANGE, ZLANHE, DGET06
225*     ..
226*     .. External Subroutines ..
227      EXTERNAL           ALAERH, ALAHD, ALASUM, ZERRHE, ZGESVD, ZGET04,
228     $                   ZLACPY, ZLARHS, ZLATB4, ZLATMS, ZPOT02,
229     $                   ZPOT03, ZHECON_ROOK, ZHET01_ROOK, ZHETRF_ROOK,
230     $                   ZHETRI_ROOK, ZHETRS_ROOK, XLAENV
231*     ..
232*     .. Intrinsic Functions ..
233      INTRINSIC          CONJG, MAX, MIN, SQRT
234*     ..
235*     .. Scalars in Common ..
236      LOGICAL            LERR, OK
237      CHARACTER*32       SRNAMT
238      INTEGER            INFOT, NUNIT
239*     ..
240*     .. Common blocks ..
241      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
242      COMMON             / SRNAMC / SRNAMT
243*     ..
244*     .. Data statements ..
245      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
246      DATA               UPLOS / 'U', 'L' /
247*     ..
248*     .. Executable Statements ..
249*
250*     Initialize constants and the random number seed.
251*
252      ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT
253*
254*     Test path
255*
256      PATH( 1: 1 ) = 'Zomplex precision'
257      PATH( 2: 3 ) = 'HR'
258*
259*     Path to generate matrices
260*
261      MATPATH( 1: 1 ) = 'Zomplex precision'
262      MATPATH( 2: 3 ) = 'HE'
263*
264      NRUN = 0
265      NFAIL = 0
266      NERRS = 0
267      DO 10 I = 1, 4
268         ISEED( I ) = ISEEDY( I )
269   10 CONTINUE
270*
271*     Test the error exits
272*
273      IF( TSTERR )
274     $   CALL ZERRHE( PATH, NOUT )
275      INFOT = 0
276*
277*     Set the minimum block size for which the block routine should
278*     be used, which will be later returned by ILAENV
279*
280      CALL XLAENV( 2, 2 )
281*
282*     Do for each value of N in NVAL
283*
284      DO 270 IN = 1, NN
285         N = NVAL( IN )
286         LDA = MAX( N, 1 )
287         XTYPE = 'N'
288         NIMAT = NTYPES
289         IF( N.LE.0 )
290     $      NIMAT = 1
291*
292         IZERO = 0
293*
294*        Do for each value of matrix type IMAT
295*
296         DO 260 IMAT = 1, NIMAT
297*
298*           Do the tests only if DOTYPE( IMAT ) is true.
299*
300            IF( .NOT.DOTYPE( IMAT ) )
301     $         GO TO 260
302*
303*           Skip types 3, 4, 5, or 6 if the matrix size is too small.
304*
305            ZEROT = IMAT.GE.3 .AND. IMAT.LE.6
306            IF( ZEROT .AND. N.LT.IMAT-2 )
307     $         GO TO 260
308*
309*           Do first for UPLO = 'U', then for UPLO = 'L'
310*
311            DO 250 IUPLO = 1, 2
312               UPLO = UPLOS( IUPLO )
313*
314*                 Begin generate the test matrix A.
315*
316*                 Set up parameters with ZLATB4 for the matrix generator
317*                 based on the type of matrix to be generated.
318*
319                  CALL ZLATB4( MATPATH, IMAT, N, N, TYPE, KL, KU, ANORM,
320     $                         MODE, CNDNUM, DIST )
321*
322*                 Generate a matrix with ZLATMS.
323*
324                  SRNAMT = 'ZLATMS'
325                  CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
326     $                         CNDNUM, ANORM, KL, KU, UPLO, A, LDA,
327     $                         WORK, INFO )
328*
329*                 Check error code from ZLATMS and handle error.
330*
331                  IF( INFO.NE.0 ) THEN
332                     CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N,
333     $                            -1, -1, -1, IMAT, NFAIL, NERRS, NOUT )
334*
335*                    Skip all tests for this generated matrix
336*
337                     GO TO 250
338                  END IF
339*
340*                 For matrix types 3-6, zero one or more rows and
341*                 columns of the matrix to test that INFO is returned
342*                 correctly.
343*
344                  IF( ZEROT ) THEN
345                     IF( IMAT.EQ.3 ) THEN
346                        IZERO = 1
347                     ELSE IF( IMAT.EQ.4 ) THEN
348                        IZERO = N
349                     ELSE
350                        IZERO = N / 2 + 1
351                     END IF
352*
353                     IF( IMAT.LT.6 ) THEN
354*
355*                       Set row and column IZERO to zero.
356*
357                        IF( IUPLO.EQ.1 ) THEN
358                           IOFF = ( IZERO-1 )*LDA
359                           DO 20 I = 1, IZERO - 1
360                              A( IOFF+I ) = CZERO
361   20                      CONTINUE
362                           IOFF = IOFF + IZERO
363                           DO 30 I = IZERO, N
364                              A( IOFF ) = CZERO
365                              IOFF = IOFF + LDA
366   30                      CONTINUE
367                        ELSE
368                           IOFF = IZERO
369                           DO 40 I = 1, IZERO - 1
370                              A( IOFF ) = CZERO
371                              IOFF = IOFF + LDA
372   40                      CONTINUE
373                           IOFF = IOFF - IZERO
374                           DO 50 I = IZERO, N
375                              A( IOFF+I ) = CZERO
376   50                      CONTINUE
377                        END IF
378                     ELSE
379                        IF( IUPLO.EQ.1 ) THEN
380*
381*                          Set the first IZERO rows and columns to zero.
382*
383                           IOFF = 0
384                           DO 70 J = 1, N
385                              I2 = MIN( J, IZERO )
386                              DO 60 I = 1, I2
387                                 A( IOFF+I ) = CZERO
388   60                         CONTINUE
389                              IOFF = IOFF + LDA
390   70                      CONTINUE
391                        ELSE
392*
393*                          Set the last IZERO rows and columns to zero.
394*
395                           IOFF = 0
396                           DO 90 J = 1, N
397                              I1 = MAX( J, IZERO )
398                              DO 80 I = I1, N
399                                 A( IOFF+I ) = CZERO
400   80                         CONTINUE
401                              IOFF = IOFF + LDA
402   90                      CONTINUE
403                        END IF
404                     END IF
405                  ELSE
406                     IZERO = 0
407                  END IF
408*
409*                 End generate the test matrix A.
410*
411*
412*              Do for each value of NB in NBVAL
413*
414               DO 240 INB = 1, NNB
415*
416*                 Set the optimal blocksize, which will be later
417*                 returned by ILAENV.
418*
419                  NB = NBVAL( INB )
420                  CALL XLAENV( 1, NB )
421*
422*                 Copy the test matrix A into matrix AFAC which
423*                 will be factorized in place. This is needed to
424*                 preserve the test matrix A for subsequent tests.
425*
426                  CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
427*
428*                 Compute the L*D*L**T or U*D*U**T factorization of the
429*                 matrix. IWORK stores details of the interchanges and
430*                 the block structure of D. AINV is a work array for
431*                 block factorization, LWORK is the length of AINV.
432*
433                  LWORK = MAX( 2, NB )*LDA
434                  SRNAMT = 'ZHETRF_ROOK'
435                  CALL ZHETRF_ROOK( UPLO, N, AFAC, LDA, IWORK, AINV,
436     $                              LWORK, INFO )
437*
438*                 Adjust the expected value of INFO to account for
439*                 pivoting.
440*
441                  K = IZERO
442                  IF( K.GT.0 ) THEN
443  100                CONTINUE
444                     IF( IWORK( K ).LT.0 ) THEN
445                        IF( IWORK( K ).NE.-K ) THEN
446                           K = -IWORK( K )
447                           GO TO 100
448                        END IF
449                     ELSE IF( IWORK( K ).NE.K ) THEN
450                        K = IWORK( K )
451                        GO TO 100
452                     END IF
453                  END IF
454*
455*                 Check error code from ZHETRF_ROOK and handle error.
456*
457                  IF( INFO.NE.K)
458     $               CALL ALAERH( PATH, 'ZHETRF_ROOK', INFO, K,
459     $                            UPLO, N, N, -1, -1, NB, IMAT,
460     $                            NFAIL, NERRS, NOUT )
461*
462*                 Set the condition estimate flag if the INFO is not 0.
463*
464                  IF( INFO.NE.0 ) THEN
465                     TRFCON = .TRUE.
466                  ELSE
467                     TRFCON = .FALSE.
468                  END IF
469*
470*+    TEST 1
471*                 Reconstruct matrix from factors and compute residual.
472*
473                  CALL ZHET01_ROOK( UPLO, N, A, LDA, AFAC, LDA, IWORK,
474     $                              AINV, LDA, RWORK, RESULT( 1 ) )
475                  NT = 1
476*
477*+    TEST 2
478*                 Form the inverse and compute the residual,
479*                 if the factorization was competed without INFO > 0
480*                 (i.e. there is no zero rows and columns).
481*                 Do it only for the first block size.
482*
483                  IF( INB.EQ.1 .AND. .NOT.TRFCON ) THEN
484                     CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
485                     SRNAMT = 'ZHETRI_ROOK'
486                     CALL ZHETRI_ROOK( UPLO, N, AINV, LDA, IWORK, WORK,
487     $                                 INFO )
488*
489*                    Check error code from ZHETRI_ROOK and handle error.
490*
491                     IF( INFO.NE.0 )
492     $                  CALL ALAERH( PATH, 'ZHETRI_ROOK', INFO, -1,
493     $                               UPLO, N, N, -1, -1, -1, IMAT,
494     $                               NFAIL, NERRS, NOUT )
495*
496*                    Compute the residual for a Hermitian matrix times
497*                    its inverse.
498*
499                     CALL ZPOT03( UPLO, N, A, LDA, AINV, LDA, WORK, LDA,
500     $                            RWORK, RCONDC, RESULT( 2 ) )
501                     NT = 2
502                  END IF
503*
504*                 Print information about the tests that did not pass
505*                 the threshold.
506*
507                  DO 110 K = 1, NT
508                     IF( RESULT( K ).GE.THRESH ) THEN
509                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
510     $                     CALL ALAHD( NOUT, PATH )
511                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
512     $                     RESULT( K )
513                        NFAIL = NFAIL + 1
514                     END IF
515  110             CONTINUE
516                  NRUN = NRUN + NT
517*
518*+    TEST 3
519*                 Compute largest element in U or L
520*
521                  RESULT( 3 ) = ZERO
522                  DTEMP = ZERO
523*
524                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) ) /
525     $                    ( ONE-ALPHA )
526*
527                  IF( IUPLO.EQ.1 ) THEN
528*
529*                 Compute largest element in U
530*
531                     K = N
532  120                CONTINUE
533                     IF( K.LE.1 )
534     $                  GO TO 130
535*
536                     IF( IWORK( K ).GT.ZERO ) THEN
537*
538*                       Get max absolute value from elements
539*                       in column k in U
540*
541                        DTEMP = ZLANGE( 'M', K-1, 1,
542     $                          AFAC( ( K-1 )*LDA+1 ), LDA, RWORK )
543                     ELSE
544*
545*                       Get max absolute value from elements
546*                       in columns k and k-1 in U
547*
548                        DTEMP = ZLANGE( 'M', K-2, 2,
549     $                          AFAC( ( K-2 )*LDA+1 ), LDA, RWORK )
550                        K = K - 1
551*
552                     END IF
553*
554*                    DTEMP should be bounded by CONST
555*
556                     DTEMP = DTEMP - CONST + THRESH
557                     IF( DTEMP.GT.RESULT( 3 ) )
558     $                  RESULT( 3 ) = DTEMP
559*
560                     K = K - 1
561*
562                     GO TO 120
563  130                CONTINUE
564*
565                  ELSE
566*
567*                 Compute largest element in L
568*
569                     K = 1
570  140                CONTINUE
571                     IF( K.GE.N )
572     $                  GO TO 150
573*
574                     IF( IWORK( K ).GT.ZERO ) THEN
575*
576*                       Get max absolute value from elements
577*                       in column k in L
578*
579                        DTEMP = ZLANGE( 'M', N-K, 1,
580     $                          AFAC( ( K-1 )*LDA+K+1 ), LDA, RWORK )
581                     ELSE
582*
583*                       Get max absolute value from elements
584*                       in columns k and k+1 in L
585*
586                        DTEMP = ZLANGE( 'M', N-K-1, 2,
587     $                          AFAC( ( K-1 )*LDA+K+2 ), LDA, RWORK )
588                        K = K + 1
589*
590                     END IF
591*
592*                    DTEMP should be bounded by CONST
593*
594                     DTEMP = DTEMP - CONST + THRESH
595                     IF( DTEMP.GT.RESULT( 3 ) )
596     $                  RESULT( 3 ) = DTEMP
597*
598                     K = K + 1
599*
600                     GO TO 140
601  150                CONTINUE
602                  END IF
603*
604*
605*+    TEST 4
606*                 Compute largest 2-Norm (condition number)
607*                 of 2-by-2 diag blocks
608*
609                  RESULT( 4 ) = ZERO
610                  DTEMP = ZERO
611*
612                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) )*
613     $                    ( ( ONE + ALPHA ) / ( ONE - ALPHA ) )
614                  CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
615*
616                  IF( IUPLO.EQ.1 ) THEN
617*
618*                    Loop backward for UPLO = 'U'
619*
620                     K = N
621  160                CONTINUE
622                     IF( K.LE.1 )
623     $                  GO TO 170
624*
625                     IF( IWORK( K ).LT.ZERO ) THEN
626*
627*                       Get the two singular values
628*                       (real and non-negative) of a 2-by-2 block,
629*                       store them in RWORK array
630*
631                        BLOCK( 1, 1 ) = AFAC( ( K-2 )*LDA+K-1 )
632                        BLOCK( 1, 2 ) = AFAC( (K-1)*LDA+K-1 )
633                        BLOCK( 2, 1 ) = CONJG( BLOCK( 1, 2 ) )
634                        BLOCK( 2, 2 ) = AFAC( (K-1)*LDA+K )
635*
636                        CALL ZGESVD( 'N', 'N', 2, 2, BLOCK, 2, RWORK,
637     $                               ZDUMMY, 1, ZDUMMY, 1,
638     $                               WORK, 6, RWORK( 3 ), INFO )
639*
640*
641                        SING_MAX = RWORK( 1 )
642                        SING_MIN = RWORK( 2 )
643*
644                        DTEMP = SING_MAX / SING_MIN
645*
646*                       DTEMP should be bounded by CONST
647*
648                        DTEMP = DTEMP - CONST + THRESH
649                        IF( DTEMP.GT.RESULT( 4 ) )
650     $                     RESULT( 4 ) = DTEMP
651                        K = K - 1
652*
653                     END IF
654*
655                     K = K - 1
656*
657                     GO TO 160
658  170                CONTINUE
659*
660                  ELSE
661*
662*                    Loop forward for UPLO = 'L'
663*
664                     K = 1
665  180                CONTINUE
666                     IF( K.GE.N )
667     $                  GO TO 190
668*
669                     IF( IWORK( K ).LT.ZERO ) THEN
670*
671*                       Get the two singular values
672*                       (real and non-negative) of a 2-by-2 block,
673*                       store them in RWORK array
674*
675                        BLOCK( 1, 1 ) = AFAC( ( K-1 )*LDA+K )
676                        BLOCK( 2, 1 ) = AFAC( ( K-1 )*LDA+K+1 )
677                        BLOCK( 1, 2 ) = CONJG( BLOCK( 2, 1 ) )
678                        BLOCK( 2, 2 ) = AFAC( K*LDA+K+1 )
679*
680                        CALL ZGESVD( 'N', 'N', 2, 2, BLOCK, 2, RWORK,
681     $                               ZDUMMY, 1, ZDUMMY, 1,
682     $                               WORK, 6, RWORK(3), INFO )
683*
684                        SING_MAX = RWORK( 1 )
685                        SING_MIN = RWORK( 2 )
686*
687                        DTEMP = SING_MAX / SING_MIN
688*
689*                       DTEMP should be bounded by CONST
690*
691                        DTEMP = 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