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