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