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 December 2016
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.7.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*     December 2016
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     $                   IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
215     $                   N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
216      DOUBLE PRECISION   ALPHA, ANORM, CNDNUM, CONST, SING_MAX,
217     $                   SING_MIN, RCOND, RCONDC, DTEMP
218*     ..
219*     .. Local Arrays ..
220      CHARACTER          UPLOS( 2 )
221      INTEGER            ISEED( 4 ), ISEEDY( 4 )
222      DOUBLE PRECISION   RESULT( NTESTS )
223      COMPLEX*16         BLOCK( 2, 2 ), ZDUMMY( 1 )
224*     ..
225*     .. External Functions ..
226      DOUBLE PRECISION   ZLANGE, ZLANHE, DGET06
227      EXTERNAL           ZLANGE, ZLANHE, DGET06
228*     ..
229*     .. External Subroutines ..
230      EXTERNAL           ALAERH, ALAHD, ALASUM, ZERRHE, ZGESVD, ZGET04,
231     $                   ZLACPY, ZLARHS, ZLATB4, ZLATMS, ZPOT02,
232     $                   ZPOT03, ZHECON_ROOK, ZHET01_ROOK, ZHETRF_ROOK,
233     $                   ZHETRI_ROOK, ZHETRS_ROOK, XLAENV
234*     ..
235*     .. Intrinsic Functions ..
236      INTRINSIC          CONJG, MAX, MIN, SQRT
237*     ..
238*     .. Scalars in Common ..
239      LOGICAL            LERR, OK
240      CHARACTER*32       SRNAMT
241      INTEGER            INFOT, NUNIT
242*     ..
243*     .. Common blocks ..
244      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
245      COMMON             / SRNAMC / SRNAMT
246*     ..
247*     .. Data statements ..
248      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
249      DATA               UPLOS / 'U', 'L' /
250*     ..
251*     .. Executable Statements ..
252*
253*     Initialize constants and the random number seed.
254*
255      ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT
256*
257*     Test path
258*
259      PATH( 1: 1 ) = 'Zomplex precision'
260      PATH( 2: 3 ) = 'HR'
261*
262*     Path to generate matrices
263*
264      MATPATH( 1: 1 ) = 'Zomplex precision'
265      MATPATH( 2: 3 ) = 'HE'
266*
267      NRUN = 0
268      NFAIL = 0
269      NERRS = 0
270      DO 10 I = 1, 4
271         ISEED( I ) = ISEEDY( I )
272   10 CONTINUE
273*
274*     Test the error exits
275*
276      IF( TSTERR )
277     $   CALL ZERRHE( PATH, NOUT )
278      INFOT = 0
279*
280*     Set the minimum block size for which the block routine should
281*     be used, which will be later returned by ILAENV
282*
283      CALL XLAENV( 2, 2 )
284*
285*     Do for each value of N in NVAL
286*
287      DO 270 IN = 1, NN
288         N = NVAL( IN )
289         LDA = MAX( N, 1 )
290         XTYPE = 'N'
291         NIMAT = NTYPES
292         IF( N.LE.0 )
293     $      NIMAT = 1
294*
295         IZERO = 0
296*
297*        Do for each value of matrix type IMAT
298*
299         DO 260 IMAT = 1, NIMAT
300*
301*           Do the tests only if DOTYPE( IMAT ) is true.
302*
303            IF( .NOT.DOTYPE( IMAT ) )
304     $         GO TO 260
305*
306*           Skip types 3, 4, 5, or 6 if the matrix size is too small.
307*
308            ZEROT = IMAT.GE.3 .AND. IMAT.LE.6
309            IF( ZEROT .AND. N.LT.IMAT-2 )
310     $         GO TO 260
311*
312*           Do first for UPLO = 'U', then for UPLO = 'L'
313*
314            DO 250 IUPLO = 1, 2
315               UPLO = UPLOS( IUPLO )
316*
317*                 Begin generate the test matrix A.
318*
319*                 Set up parameters with ZLATB4 for the matrix generator
320*                 based on the type of matrix to be generated.
321*
322                  CALL ZLATB4( MATPATH, IMAT, N, N, TYPE, KL, KU, ANORM,
323     $                         MODE, CNDNUM, DIST )
324*
325*                 Generate a matrix with ZLATMS.
326*
327                  SRNAMT = 'ZLATMS'
328                  CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
329     $                         CNDNUM, ANORM, KL, KU, UPLO, A, LDA,
330     $                         WORK, INFO )
331*
332*                 Check error code from ZLATMS and handle error.
333*
334                  IF( INFO.NE.0 ) THEN
335                     CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N,
336     $                            -1, -1, -1, IMAT, NFAIL, NERRS, NOUT )
337*
338*                    Skip all tests for this generated matrix
339*
340                     GO TO 250
341                  END IF
342*
343*                 For matrix types 3-6, zero one or more rows and
344*                 columns of the matrix to test that INFO is returned
345*                 correctly.
346*
347                  IF( ZEROT ) THEN
348                     IF( IMAT.EQ.3 ) THEN
349                        IZERO = 1
350                     ELSE IF( IMAT.EQ.4 ) THEN
351                        IZERO = N
352                     ELSE
353                        IZERO = N / 2 + 1
354                     END IF
355*
356                     IF( IMAT.LT.6 ) THEN
357*
358*                       Set row and column IZERO to zero.
359*
360                        IF( IUPLO.EQ.1 ) THEN
361                           IOFF = ( IZERO-1 )*LDA
362                           DO 20 I = 1, IZERO - 1
363                              A( IOFF+I ) = CZERO
364   20                      CONTINUE
365                           IOFF = IOFF + IZERO
366                           DO 30 I = IZERO, N
367                              A( IOFF ) = CZERO
368                              IOFF = IOFF + LDA
369   30                      CONTINUE
370                        ELSE
371                           IOFF = IZERO
372                           DO 40 I = 1, IZERO - 1
373                              A( IOFF ) = CZERO
374                              IOFF = IOFF + LDA
375   40                      CONTINUE
376                           IOFF = IOFF - IZERO
377                           DO 50 I = IZERO, N
378                              A( IOFF+I ) = CZERO
379   50                      CONTINUE
380                        END IF
381                     ELSE
382                        IF( IUPLO.EQ.1 ) THEN
383*
384*                          Set the first IZERO rows and columns to zero.
385*
386                           IOFF = 0
387                           DO 70 J = 1, N
388                              I2 = MIN( J, IZERO )
389                              DO 60 I = 1, I2
390                                 A( IOFF+I ) = CZERO
391   60                         CONTINUE
392                              IOFF = IOFF + LDA
393   70                      CONTINUE
394                        ELSE
395*
396*                          Set the last IZERO rows and columns to zero.
397*
398                           IOFF = 0
399                           DO 90 J = 1, N
400                              I1 = MAX( J, IZERO )
401                              DO 80 I = I1, N
402                                 A( IOFF+I ) = CZERO
403   80                         CONTINUE
404                              IOFF = IOFF + LDA
405   90                      CONTINUE
406                        END IF
407                     END IF
408                  ELSE
409                     IZERO = 0
410                  END IF
411*
412*                 End generate the test matrix A.
413*
414*
415*              Do for each value of NB in NBVAL
416*
417               DO 240 INB = 1, NNB
418*
419*                 Set the optimal blocksize, which will be later
420*                 returned by ILAENV.
421*
422                  NB = NBVAL( INB )
423                  CALL XLAENV( 1, NB )
424*
425*                 Copy the test matrix A into matrix AFAC which
426*                 will be factorized in place. This is needed to
427*                 preserve the test matrix A for subsequent tests.
428*
429                  CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
430*
431*                 Compute the L*D*L**T or U*D*U**T factorization of the
432*                 matrix. IWORK stores details of the interchanges and
433*                 the block structure of D. AINV is a work array for
434*                 block factorization, LWORK is the length of AINV.
435*
436                  LWORK = MAX( 2, NB )*LDA
437                  SRNAMT = 'ZHETRF_ROOK'
438                  CALL ZHETRF_ROOK( UPLO, N, AFAC, LDA, IWORK, AINV,
439     $                              LWORK, INFO )
440*
441*                 Adjust the expected value of INFO to account for
442*                 pivoting.
443*
444                  K = IZERO
445                  IF( K.GT.0 ) THEN
446  100                CONTINUE
447                     IF( IWORK( K ).LT.0 ) THEN
448                        IF( IWORK( K ).NE.-K ) THEN
449                           K = -IWORK( K )
450                           GO TO 100
451                        END IF
452                     ELSE IF( IWORK( K ).NE.K ) THEN
453                        K = IWORK( K )
454                        GO TO 100
455                     END IF
456                  END IF
457*
458*                 Check error code from ZHETRF_ROOK and handle error.
459*
460                  IF( INFO.NE.K)
461     $               CALL ALAERH( PATH, 'ZHETRF_ROOK', INFO, K,
462     $                            UPLO, N, N, -1, -1, NB, IMAT,
463     $                            NFAIL, NERRS, NOUT )
464*
465*                 Set the condition estimate flag if the INFO is not 0.
466*
467                  IF( INFO.NE.0 ) THEN
468                     TRFCON = .TRUE.
469                  ELSE
470                     TRFCON = .FALSE.
471                  END IF
472*
473*+    TEST 1
474*                 Reconstruct matrix from factors and compute residual.
475*
476                  CALL ZHET01_ROOK( UPLO, N, A, LDA, AFAC, LDA, IWORK,
477     $                              AINV, LDA, RWORK, RESULT( 1 ) )
478                  NT = 1
479*
480*+    TEST 2
481*                 Form the inverse and compute the residual,
482*                 if the factorization was competed without INFO > 0
483*                 (i.e. there is no zero rows and columns).
484*                 Do it only for the first block size.
485*
486                  IF( INB.EQ.1 .AND. .NOT.TRFCON ) THEN
487                     CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
488                     SRNAMT = 'ZHETRI_ROOK'
489                     CALL ZHETRI_ROOK( UPLO, N, AINV, LDA, IWORK, WORK,
490     $                                 INFO )
491*
492*                    Check error code from ZHETRI_ROOK and handle error.
493*
494                     IF( INFO.NE.0 )
495     $                  CALL ALAERH( PATH, 'ZHETRI_ROOK', INFO, -1,
496     $                               UPLO, N, N, -1, -1, -1, IMAT,
497     $                               NFAIL, NERRS, NOUT )
498*
499*                    Compute the residual for a Hermitian matrix times
500*                    its inverse.
501*
502                     CALL ZPOT03( UPLO, N, A, LDA, AINV, LDA, WORK, LDA,
503     $                            RWORK, RCONDC, RESULT( 2 ) )
504                     NT = 2
505                  END IF
506*
507*                 Print information about the tests that did not pass
508*                 the threshold.
509*
510                  DO 110 K = 1, NT
511                     IF( RESULT( K ).GE.THRESH ) THEN
512                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
513     $                     CALL ALAHD( NOUT, PATH )
514                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
515     $                     RESULT( K )
516                        NFAIL = NFAIL + 1
517                     END IF
518  110             CONTINUE
519                  NRUN = NRUN + NT
520*
521*+    TEST 3
522*                 Compute largest element in U or L
523*
524                  RESULT( 3 ) = ZERO
525                  DTEMP = ZERO
526*
527                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) ) /
528     $                    ( ONE-ALPHA )
529*
530                  IF( IUPLO.EQ.1 ) THEN
531*
532*                 Compute largest element in U
533*
534                     K = N
535  120                CONTINUE
536                     IF( K.LE.1 )
537     $                  GO TO 130
538*
539                     IF( IWORK( K ).GT.ZERO ) THEN
540*
541*                       Get max absolute value from elements
542*                       in column k in U
543*
544                        DTEMP = ZLANGE( 'M', K-1, 1,
545     $                          AFAC( ( K-1 )*LDA+1 ), LDA, RWORK )
546                     ELSE
547*
548*                       Get max absolute value from elements
549*                       in columns k and k-1 in U
550*
551                        DTEMP = ZLANGE( 'M', K-2, 2,
552     $                          AFAC( ( K-2 )*LDA+1 ), LDA, RWORK )
553                        K = K - 1
554*
555                     END IF
556*
557*                    DTEMP should be bounded by CONST
558*
559                     DTEMP = DTEMP - CONST + THRESH
560                     IF( DTEMP.GT.RESULT( 3 ) )
561     $                  RESULT( 3 ) = DTEMP
562*
563                     K = K - 1
564*
565                     GO TO 120
566  130                CONTINUE
567*
568                  ELSE
569*
570*                 Compute largest element in L
571*
572                     K = 1
573  140                CONTINUE
574                     IF( K.GE.N )
575     $                  GO TO 150
576*
577                     IF( IWORK( K ).GT.ZERO ) THEN
578*
579*                       Get max absolute value from elements
580*                       in column k in L
581*
582                        DTEMP = ZLANGE( 'M', N-K, 1,
583     $                          AFAC( ( K-1 )*LDA+K+1 ), LDA, RWORK )
584                     ELSE
585*
586*                       Get max absolute value from elements
587*                       in columns k and k+1 in L
588*
589                        DTEMP = ZLANGE( 'M', N-K-1, 2,
590     $                          AFAC( ( K-1 )*LDA+K+2 ), LDA, RWORK )
591                        K = K + 1
592*
593                     END IF
594*
595*                    DTEMP should be bounded by CONST
596*
597                     DTEMP = DTEMP - CONST + THRESH
598                     IF( DTEMP.GT.RESULT( 3 ) )
599     $                  RESULT( 3 ) = DTEMP
600*
601                     K = K + 1
602*
603                     GO TO 140
604  150                CONTINUE
605                  END IF
606*
607*
608*+    TEST 4
609*                 Compute largest 2-Norm (condition number)
610*                 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 singular values
631*                       (real and non-negative) of a 2-by-2 block,
632*                       store them in RWORK array
633*
634                        BLOCK( 1, 1 ) = AFAC( ( K-2 )*LDA+K-1 )
635                        BLOCK( 1, 2 ) = AFAC( (K-1)*LDA+K-1 )
636                        BLOCK( 2, 1 ) = CONJG( BLOCK( 1, 2 ) )
637                        BLOCK( 2, 2 ) = AFAC( (K-1)*LDA+K )
638*
639                        CALL ZGESVD( 'N', 'N', 2, 2, BLOCK, 2, RWORK,
640     $                               ZDUMMY, 1, ZDUMMY, 1,
641     $                               WORK, 6, RWORK( 3 ), INFO )
642*
643*
644                        SING_MAX = RWORK( 1 )
645                        SING_MIN = RWORK( 2 )
646*
647                        DTEMP = SING_MAX / SING_MIN
648*
649*                       DTEMP should be bounded by CONST
650*
651                        DTEMP = DTEMP - CONST + THRESH
652                        IF( DTEMP.GT.RESULT( 4 ) )
653     $                     RESULT( 4 ) = DTEMP
654                        K = K - 1
655*
656                     END IF
657*
658                     K = K - 1
659*
660                     GO TO 160
661  170                CONTINUE
662*
663                  ELSE
664*
665*                    Loop forward for UPLO = 'L'
666*
667                     K = 1
668  180                CONTINUE
669                     IF( K.GE.N )
670     $                  GO TO 190
671*
672                     IF( IWORK( K ).LT.ZERO ) THEN
673*
674*                       Get the two singular values
675*                       (real and non-negative) of a 2-by-2 block,
676*                       store them in RWORK array
677*
678                        BLOCK( 1, 1 ) = AFAC( ( K-1 )*LDA+K )
679                        BLOCK( 2, 1 ) = AFAC( ( K-1 )*LDA+K+1 )
680                        BLOCK( 1, 2 ) = CONJG( BLOCK( 2, 1 ) )
681                        BLOCK( 2, 2 ) = AFAC( K*LDA+K+1 )
682*
683                        CALL ZGESVD( 'N', 'N', 2, 2, BLOCK, 2, RWORK,
684     $                               ZDUMMY, 1, ZDUMMY, 1,
685     $                               WORK, 6, RWORK(3), INFO )
686*
687                        SING_MAX = RWORK( 1 )
688                        SING_MIN = RWORK( 2 )
689*
690                        DTEMP = SING_MAX / SING_MIN
691*
692*                       DTEMP should be bounded by CONST
693*
694                        DTEMP = DTEMP - CONST + THRESH
695                        IF( DTEMP.GT.RESULT( 4 ) )
696     $                     RESULT( 4 ) = DTEMP
697                        K = K + 1
698*
699                     END IF
700*
701                     K = K + 1
702*
703                     GO TO 180
704  190                CONTINUE
705                  END IF
706*
707*                 Print information about the tests that did not pass
708*                 the threshold.
709*
710                  DO 200 K = 3, 4
711                     IF( RESULT( K ).GE.THRESH ) THEN
712                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
713     $                     CALL ALAHD( NOUT, PATH )
714                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
715     $                     RESULT( K )
716                        NFAIL = NFAIL + 1
717                     END IF
718  200             CONTINUE
719                  NRUN = NRUN + 2
720*
721*                 Skip the other tests if this is not the first block
722*                 size.
723*
724                  IF( INB.GT.1 )
725     $               GO TO 240
726*
727*                 Do only the condition estimate if INFO is not 0.
728*
729                  IF( TRFCON ) THEN
730                     RCONDC = ZERO
731                     GO TO 230
732                  END IF
733*
734*                 Do for each value of NRHS in NSVAL.
735*
736                  DO 220 IRHS = 1, NNS
737                     NRHS = NSVAL( IRHS )
738*
739*                    Begin loop over NRHS values
740*
741*
742*+    TEST 5 ( Using TRS_ROOK)
743*                 Solve and compute residual for  A * X = B.
744*
745*                    Choose a set of NRHS random solution vectors
746*                    stored in XACT and set up the right hand side B
747*
748                     SRNAMT = 'ZLARHS'
749                     CALL ZLARHS( MATPATH, XTYPE, UPLO, ' ', N, N,
750     $                            KL, KU, NRHS, A, LDA, XACT, LDA,
751     $                            B, LDA, ISEED, INFO )
752                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
753*
754                     SRNAMT = 'ZHETRS_ROOK'
755                     CALL ZHETRS_ROOK( UPLO, N, NRHS, AFAC, LDA, IWORK,
756     $                                 X, LDA, INFO )
757*
758*                    Check error code from ZHETRS_ROOK and handle error.
759*
760                     IF( INFO.NE.0 )
761     $                  CALL ALAERH( PATH, 'ZHETRS_ROOK', INFO, 0,
762     $                               UPLO, N, N, -1, -1, NRHS, IMAT,
763     $                               NFAIL, NERRS, NOUT )
764*
765                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK, LDA )
766*
767*                    Compute the residual for the solution
768*
769                     CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA, WORK,
770     $                            LDA, RWORK, RESULT( 5 ) )
771*
772*+    TEST 6
773*                 Check solution from generated exact solution.
774*
775                     CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
776     $                            RESULT( 6 ) )
777*
778*                    Print information about the tests that did not pass
779*                    the threshold.
780*
781                     DO 210 K = 5, 6
782                        IF( RESULT( K ).GE.THRESH ) THEN
783                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
784     $                        CALL ALAHD( NOUT, PATH )
785                           WRITE( NOUT, FMT = 9998 )UPLO, N, NRHS,
786     $                        IMAT, K, RESULT( K )
787                           NFAIL = NFAIL + 1
788                        END IF
789  210                CONTINUE
790                     NRUN = NRUN + 2
791*
792*                 End do for each value of NRHS in NSVAL.
793*
794  220             CONTINUE
795*
796*+    TEST 7
797*                 Get an estimate of RCOND = 1/CNDNUM.
798*
799  230             CONTINUE
800                  ANORM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
801                  SRNAMT = 'ZHECON_ROOK'
802                  CALL ZHECON_ROOK( UPLO, N, AFAC, LDA, IWORK, ANORM,
803     $                              RCOND, WORK, INFO )
804*
805*                 Check error code from ZHECON_ROOK and handle error.
806*
807                  IF( INFO.NE.0 )
808     $               CALL ALAERH( PATH, 'ZHECON_ROOK', INFO, 0,
809     $                             UPLO, N, N, -1, -1, -1, IMAT,
810     $                             NFAIL, NERRS, NOUT )
811*
812*                 Compute the test ratio to compare values of RCOND
813*
814                  RESULT( 7 ) = DGET06( RCOND, RCONDC )
815*
816*                 Print information about the tests that did not pass
817*                 the threshold.
818*
819                  IF( RESULT( 7 ).GE.THRESH ) THEN
820                     IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
821     $                  CALL ALAHD( NOUT, PATH )
822                     WRITE( NOUT, FMT = 9997 )UPLO, N, IMAT, 7,
823     $                  RESULT( 7 )
824                     NFAIL = NFAIL + 1
825                  END IF
826                  NRUN = NRUN + 1
827  240          CONTINUE
828*
829  250       CONTINUE
830  260    CONTINUE
831  270 CONTINUE
832*
833*     Print a summary of the results.
834*
835      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
836*
837 9999 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NB =', I4, ', type ',
838     $      I2, ', test ', I2, ', ratio =', G12.5 )
839 9998 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NRHS=', I3, ', type ',
840     $      I2, ', test ', I2, ', ratio =', G12.5 )
841 9997 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ',', 10X, ' type ', I2,
842     $      ', test ', I2, ', ratio =', G12.5 )
843      RETURN
844*
845*     End of ZCHKHE_ROOK
846*
847      END
848