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