1*> \brief \b SCHKGB
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 SCHKGB( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
12*                          NSVAL, THRESH, TSTERR, A, LA, AFAC, LAFAC, B,
13*                          X, XACT, WORK, RWORK, IWORK, NOUT )
14*
15*       .. Scalar Arguments ..
16*       LOGICAL            TSTERR
17*       INTEGER            LA, LAFAC, NM, NN, NNB, NNS, NOUT
18*       REAL               THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
23*      $                   NVAL( * )
24*       REAL               A( * ), AFAC( * ), B( * ), RWORK( * ),
25*      $                   WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> SCHKGB tests SGBTRF, -TRS, -RFS, and -CON
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] NM
49*> \verbatim
50*>          NM is INTEGER
51*>          The number of values of M contained in the vector MVAL.
52*> \endverbatim
53*>
54*> \param[in] MVAL
55*> \verbatim
56*>          MVAL is INTEGER array, dimension (NM)
57*>          The values of the matrix row dimension M.
58*> \endverbatim
59*>
60*> \param[in] NN
61*> \verbatim
62*>          NN is INTEGER
63*>          The number of values of N contained in the vector NVAL.
64*> \endverbatim
65*>
66*> \param[in] NVAL
67*> \verbatim
68*>          NVAL is INTEGER array, dimension (NN)
69*>          The values of the matrix column dimension N.
70*> \endverbatim
71*>
72*> \param[in] NNB
73*> \verbatim
74*>          NNB is INTEGER
75*>          The number of values of NB contained in the vector NBVAL.
76*> \endverbatim
77*>
78*> \param[in] NBVAL
79*> \verbatim
80*>          NBVAL is INTEGER array, dimension (NNB)
81*>          The values of the blocksize NB.
82*> \endverbatim
83*>
84*> \param[in] NNS
85*> \verbatim
86*>          NNS is INTEGER
87*>          The number of values of NRHS contained in the vector NSVAL.
88*> \endverbatim
89*>
90*> \param[in] NSVAL
91*> \verbatim
92*>          NSVAL is INTEGER array, dimension (NNS)
93*>          The values of the number of right hand sides NRHS.
94*> \endverbatim
95*>
96*> \param[in] THRESH
97*> \verbatim
98*>          THRESH is REAL
99*>          The threshold value for the test ratios.  A result is
100*>          included in the output file if RESULT >= THRESH.  To have
101*>          every test ratio printed, use THRESH = 0.
102*> \endverbatim
103*>
104*> \param[in] TSTERR
105*> \verbatim
106*>          TSTERR is LOGICAL
107*>          Flag that indicates whether error exits are to be tested.
108*> \endverbatim
109*>
110*> \param[out] A
111*> \verbatim
112*>          A is REAL array, dimension (LA)
113*> \endverbatim
114*>
115*> \param[in] LA
116*> \verbatim
117*>          LA is INTEGER
118*>          The length of the array A.  LA >= (KLMAX+KUMAX+1)*NMAX
119*>          where KLMAX is the largest entry in the local array KLVAL,
120*>                KUMAX is the largest entry in the local array KUVAL and
121*>                NMAX is the largest entry in the input array NVAL.
122*> \endverbatim
123*>
124*> \param[out] AFAC
125*> \verbatim
126*>          AFAC is REAL array, dimension (LAFAC)
127*> \endverbatim
128*>
129*> \param[in] LAFAC
130*> \verbatim
131*>          LAFAC is INTEGER
132*>          The length of the array AFAC. LAFAC >= (2*KLMAX+KUMAX+1)*NMAX
133*>          where KLMAX is the largest entry in the local array KLVAL,
134*>                KUMAX is the largest entry in the local array KUVAL and
135*>                NMAX is the largest entry in the input array NVAL.
136*> \endverbatim
137*>
138*> \param[out] B
139*> \verbatim
140*>          B is REAL array, dimension (NMAX*NSMAX)
141*>          where NSMAX is the largest entry in NSVAL.
142*> \endverbatim
143*>
144*> \param[out] X
145*> \verbatim
146*>          X is REAL array, dimension (NMAX*NSMAX)
147*> \endverbatim
148*>
149*> \param[out] XACT
150*> \verbatim
151*>          XACT is REAL array, dimension (NMAX*NSMAX)
152*> \endverbatim
153*>
154*> \param[out] WORK
155*> \verbatim
156*>          WORK is REAL array, dimension
157*>                      (NMAX*max(3,NSMAX,NMAX))
158*> \endverbatim
159*>
160*> \param[out] RWORK
161*> \verbatim
162*>          RWORK is REAL array, dimension
163*>                      (NMAX+2*NSMAX)
164*> \endverbatim
165*>
166*> \param[out] IWORK
167*> \verbatim
168*>          IWORK is INTEGER array, dimension (2*NMAX)
169*> \endverbatim
170*>
171*> \param[in] NOUT
172*> \verbatim
173*>          NOUT is INTEGER
174*>          The unit number for output.
175*> \endverbatim
176*
177*  Authors:
178*  ========
179*
180*> \author Univ. of Tennessee
181*> \author Univ. of California Berkeley
182*> \author Univ. of Colorado Denver
183*> \author NAG Ltd.
184*
185*> \ingroup single_lin
186*
187*  =====================================================================
188      SUBROUTINE SCHKGB( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
189     $                   NSVAL, THRESH, TSTERR, A, LA, AFAC, LAFAC, B,
190     $                   X, XACT, WORK, RWORK, IWORK, NOUT )
191*
192*  -- LAPACK test routine --
193*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
194*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195*
196*     .. Scalar Arguments ..
197      LOGICAL            TSTERR
198      INTEGER            LA, LAFAC, NM, NN, NNB, NNS, NOUT
199      REAL               THRESH
200*     ..
201*     .. Array Arguments ..
202      LOGICAL            DOTYPE( * )
203      INTEGER            IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
204     $                   NVAL( * )
205      REAL               A( * ), AFAC( * ), B( * ), RWORK( * ),
206     $                   WORK( * ), X( * ), XACT( * )
207*     ..
208*
209*  =====================================================================
210*
211*     .. Parameters ..
212      REAL               ONE, ZERO
213      PARAMETER          ( ONE = 1.0E+0, ZERO = 0.0E+0 )
214      INTEGER            NTYPES, NTESTS
215      PARAMETER          ( NTYPES = 8, NTESTS = 7 )
216      INTEGER            NBW, NTRAN
217      PARAMETER          ( NBW = 4, NTRAN = 3 )
218*     ..
219*     .. Local Scalars ..
220      LOGICAL            TRFCON, ZEROT
221      CHARACTER          DIST, NORM, TRANS, TYPE, XTYPE
222      CHARACTER*3        PATH
223      INTEGER            I, I1, I2, IKL, IKU, IM, IMAT, IN, INB, INFO,
224     $                   IOFF, IRHS, ITRAN, IZERO, J, K, KL, KOFF, KU,
225     $                   LDA, LDAFAC, LDB, M, MODE, N, NB, NERRS, NFAIL,
226     $                   NIMAT, NKL, NKU, NRHS, NRUN
227      REAL               AINVNM, ANORM, ANORMI, ANORMO, CNDNUM, RCOND,
228     $                   RCONDC, RCONDI, RCONDO
229*     ..
230*     .. Local Arrays ..
231      CHARACTER          TRANSS( NTRAN )
232      INTEGER            ISEED( 4 ), ISEEDY( 4 ), KLVAL( NBW ),
233     $                   KUVAL( NBW )
234      REAL               RESULT( NTESTS )
235*     ..
236*     .. External Functions ..
237      REAL               SGET06, SLANGB, SLANGE
238      EXTERNAL           SGET06, SLANGB, SLANGE
239*     ..
240*     .. External Subroutines ..
241      EXTERNAL           ALAERH, ALAHD, ALASUM, SCOPY, SERRGE, SGBCON,
242     $                   SGBRFS, SGBT01, SGBT02, SGBT05, SGBTRF, SGBTRS,
243     $                   SGET04, SLACPY, SLARHS, SLASET, SLATB4, SLATMS,
244     $                   XLAENV
245*     ..
246*     .. Intrinsic Functions ..
247      INTRINSIC          MAX, MIN
248*     ..
249*     .. Scalars in Common ..
250      LOGICAL            LERR, OK
251      CHARACTER*32       SRNAMT
252      INTEGER            INFOT, NUNIT
253*     ..
254*     .. Common blocks ..
255      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
256      COMMON             / SRNAMC / SRNAMT
257*     ..
258*     .. Data statements ..
259      DATA               ISEEDY / 1988, 1989, 1990, 1991 / ,
260     $                   TRANSS / 'N', 'T', 'C' /
261*     ..
262*     .. Executable Statements ..
263*
264*     Initialize constants and the random number seed.
265*
266      PATH( 1: 1 ) = 'Single precision'
267      PATH( 2: 3 ) = 'GB'
268      NRUN = 0
269      NFAIL = 0
270      NERRS = 0
271      DO 10 I = 1, 4
272         ISEED( I ) = ISEEDY( I )
273   10 CONTINUE
274*
275*     Test the error exits
276*
277      IF( TSTERR )
278     $   CALL SERRGE( PATH, NOUT )
279      INFOT = 0
280      CALL XLAENV( 2, 2 )
281*
282*     Initialize the first value for the lower and upper bandwidths.
283*
284      KLVAL( 1 ) = 0
285      KUVAL( 1 ) = 0
286*
287*     Do for each value of M in MVAL
288*
289      DO 160 IM = 1, NM
290         M = MVAL( IM )
291*
292*        Set values to use for the lower bandwidth.
293*
294         KLVAL( 2 ) = M + ( M+1 ) / 4
295*
296*        KLVAL( 2 ) = MAX( M-1, 0 )
297*
298         KLVAL( 3 ) = ( 3*M-1 ) / 4
299         KLVAL( 4 ) = ( M+1 ) / 4
300*
301*        Do for each value of N in NVAL
302*
303         DO 150 IN = 1, NN
304            N = NVAL( IN )
305            XTYPE = 'N'
306*
307*           Set values to use for the upper bandwidth.
308*
309            KUVAL( 2 ) = N + ( N+1 ) / 4
310*
311*           KUVAL( 2 ) = MAX( N-1, 0 )
312*
313            KUVAL( 3 ) = ( 3*N-1 ) / 4
314            KUVAL( 4 ) = ( N+1 ) / 4
315*
316*           Set limits on the number of loop iterations.
317*
318            NKL = MIN( M+1, 4 )
319            IF( N.EQ.0 )
320     $         NKL = 2
321            NKU = MIN( N+1, 4 )
322            IF( M.EQ.0 )
323     $         NKU = 2
324            NIMAT = NTYPES
325            IF( M.LE.0 .OR. N.LE.0 )
326     $         NIMAT = 1
327*
328            DO 140 IKL = 1, NKL
329*
330*              Do for KL = 0, (5*M+1)/4, (3M-1)/4, and (M+1)/4. This
331*              order makes it easier to skip redundant values for small
332*              values of M.
333*
334               KL = KLVAL( IKL )
335               DO 130 IKU = 1, NKU
336*
337*                 Do for KU = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This
338*                 order makes it easier to skip redundant values for
339*                 small values of N.
340*
341                  KU = KUVAL( IKU )
342*
343*                 Check that A and AFAC are big enough to generate this
344*                 matrix.
345*
346                  LDA = KL + KU + 1
347                  LDAFAC = 2*KL + KU + 1
348                  IF( ( LDA*N ).GT.LA .OR. ( LDAFAC*N ).GT.LAFAC ) THEN
349                     IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
350     $                  CALL ALAHD( NOUT, PATH )
351                     IF( N*( KL+KU+1 ).GT.LA ) THEN
352                        WRITE( NOUT, FMT = 9999 )LA, M, N, KL, KU,
353     $                     N*( KL+KU+1 )
354                        NERRS = NERRS + 1
355                     END IF
356                     IF( N*( 2*KL+KU+1 ).GT.LAFAC ) THEN
357                        WRITE( NOUT, FMT = 9998 )LAFAC, M, N, KL, KU,
358     $                     N*( 2*KL+KU+1 )
359                        NERRS = NERRS + 1
360                     END IF
361                     GO TO 130
362                  END IF
363*
364                  DO 120 IMAT = 1, NIMAT
365*
366*                    Do the tests only if DOTYPE( IMAT ) is true.
367*
368                     IF( .NOT.DOTYPE( IMAT ) )
369     $                  GO TO 120
370*
371*                    Skip types 2, 3, or 4 if the matrix size is too
372*                    small.
373*
374                     ZEROT = IMAT.GE.2 .AND. IMAT.LE.4
375                     IF( ZEROT .AND. N.LT.IMAT-1 )
376     $                  GO TO 120
377*
378                     IF( .NOT.ZEROT .OR. .NOT.DOTYPE( 1 ) ) THEN
379*
380*                       Set up parameters with SLATB4 and generate a
381*                       test matrix with SLATMS.
382*
383                        CALL SLATB4( PATH, IMAT, M, N, TYPE, KL, KU,
384     $                               ANORM, MODE, CNDNUM, DIST )
385*
386                        KOFF = MAX( 1, KU+2-N )
387                        DO 20 I = 1, KOFF - 1
388                           A( I ) = ZERO
389   20                   CONTINUE
390                        SRNAMT = 'SLATMS'
391                        CALL SLATMS( M, N, DIST, ISEED, TYPE, RWORK,
392     $                               MODE, CNDNUM, ANORM, KL, KU, 'Z',
393     $                               A( KOFF ), LDA, WORK, INFO )
394*
395*                       Check the error code from SLATMS.
396*
397                        IF( INFO.NE.0 ) THEN
398                           CALL ALAERH( PATH, 'SLATMS', INFO, 0, ' ', M,
399     $                                  N, KL, KU, -1, IMAT, NFAIL,
400     $                                  NERRS, NOUT )
401                           GO TO 120
402                        END IF
403                     ELSE IF( IZERO.GT.0 ) THEN
404*
405*                       Use the same matrix for types 3 and 4 as for
406*                       type 2 by copying back the zeroed out column.
407*
408                        CALL SCOPY( I2-I1+1, B, 1, A( IOFF+I1 ), 1 )
409                     END IF
410*
411*                    For types 2, 3, and 4, zero one or more columns of
412*                    the matrix to test that INFO is returned correctly.
413*
414                     IZERO = 0
415                     IF( ZEROT ) THEN
416                        IF( IMAT.EQ.2 ) THEN
417                           IZERO = 1
418                        ELSE IF( IMAT.EQ.3 ) THEN
419                           IZERO = MIN( M, N )
420                        ELSE
421                           IZERO = MIN( M, N ) / 2 + 1
422                        END IF
423                        IOFF = ( IZERO-1 )*LDA
424                        IF( IMAT.LT.4 ) THEN
425*
426*                          Store the column to be zeroed out in B.
427*
428                           I1 = MAX( 1, KU+2-IZERO )
429                           I2 = MIN( KL+KU+1, KU+1+( M-IZERO ) )
430                           CALL SCOPY( I2-I1+1, A( IOFF+I1 ), 1, B, 1 )
431*
432                           DO 30 I = I1, I2
433                              A( IOFF+I ) = ZERO
434   30                      CONTINUE
435                        ELSE
436                           DO 50 J = IZERO, N
437                              DO 40 I = MAX( 1, KU+2-J ),
438     $                                MIN( KL+KU+1, KU+1+( M-J ) )
439                                 A( IOFF+I ) = ZERO
440   40                         CONTINUE
441                              IOFF = IOFF + LDA
442   50                      CONTINUE
443                        END IF
444                     END IF
445*
446*                    These lines, if used in place of the calls in the
447*                    loop over INB, cause the code to bomb on a Sun
448*                    SPARCstation.
449*
450*                     ANORMO = SLANGB( 'O', N, KL, KU, A, LDA, RWORK )
451*                     ANORMI = SLANGB( 'I', N, KL, KU, A, LDA, RWORK )
452*
453*                    Do for each blocksize in NBVAL
454*
455                     DO 110 INB = 1, NNB
456                        NB = NBVAL( INB )
457                        CALL XLAENV( 1, NB )
458*
459*                       Compute the LU factorization of the band matrix.
460*
461                        IF( M.GT.0 .AND. N.GT.0 )
462     $                     CALL SLACPY( 'Full', KL+KU+1, N, A, LDA,
463     $                                  AFAC( KL+1 ), LDAFAC )
464                        SRNAMT = 'SGBTRF'
465                        CALL SGBTRF( M, N, KL, KU, AFAC, LDAFAC, IWORK,
466     $                               INFO )
467*
468*                       Check error code from SGBTRF.
469*
470                        IF( INFO.NE.IZERO )
471     $                     CALL ALAERH( PATH, 'SGBTRF', INFO, IZERO,
472     $                                  ' ', M, N, KL, KU, NB, IMAT,
473     $                                  NFAIL, NERRS, NOUT )
474                        TRFCON = .FALSE.
475*
476*+    TEST 1
477*                       Reconstruct matrix from factors and compute
478*                       residual.
479*
480                        CALL SGBT01( M, N, KL, KU, A, LDA, AFAC, LDAFAC,
481     $                               IWORK, WORK, RESULT( 1 ) )
482*
483*                       Print information about the tests so far that
484*                       did not pass the threshold.
485*
486                        IF( RESULT( 1 ).GE.THRESH ) THEN
487                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
488     $                        CALL ALAHD( NOUT, PATH )
489                           WRITE( NOUT, FMT = 9997 )M, N, KL, KU, NB,
490     $                        IMAT, 1, RESULT( 1 )
491                           NFAIL = NFAIL + 1
492                        END IF
493                        NRUN = NRUN + 1
494*
495*                       Skip the remaining tests if this is not the
496*                       first block size or if M .ne. N.
497*
498                        IF( INB.GT.1 .OR. M.NE.N )
499     $                     GO TO 110
500*
501                        ANORMO = SLANGB( 'O', N, KL, KU, A, LDA, RWORK )
502                        ANORMI = SLANGB( 'I', N, KL, KU, A, LDA, RWORK )
503*
504                        IF( INFO.EQ.0 ) THEN
505*
506*                          Form the inverse of A so we can get a good
507*                          estimate of CNDNUM = norm(A) * norm(inv(A)).
508*
509                           LDB = MAX( 1, N )
510                           CALL SLASET( 'Full', N, N, ZERO, ONE, WORK,
511     $                                  LDB )
512                           SRNAMT = 'SGBTRS'
513                           CALL SGBTRS( 'No transpose', N, KL, KU, N,
514     $                                  AFAC, LDAFAC, IWORK, WORK, LDB,
515     $                                  INFO )
516*
517*                          Compute the 1-norm condition number of A.
518*
519                           AINVNM = SLANGE( 'O', N, N, WORK, LDB,
520     $                              RWORK )
521                           IF( ANORMO.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
522                              RCONDO = ONE
523                           ELSE
524                              RCONDO = ( ONE / ANORMO ) / AINVNM
525                           END IF
526*
527*                          Compute the infinity-norm condition number of
528*                          A.
529*
530                           AINVNM = SLANGE( 'I', N, N, WORK, LDB,
531     $                              RWORK )
532                           IF( ANORMI.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
533                              RCONDI = ONE
534                           ELSE
535                              RCONDI = ( ONE / ANORMI ) / AINVNM
536                           END IF
537                        ELSE
538*
539*                          Do only the condition estimate if INFO.NE.0.
540*
541                           TRFCON = .TRUE.
542                           RCONDO = ZERO
543                           RCONDI = ZERO
544                        END IF
545*
546*                       Skip the solve tests if the matrix is singular.
547*
548                        IF( TRFCON )
549     $                     GO TO 90
550*
551                        DO 80 IRHS = 1, NNS
552                           NRHS = NSVAL( IRHS )
553                           XTYPE = 'N'
554*
555                           DO 70 ITRAN = 1, NTRAN
556                              TRANS = TRANSS( ITRAN )
557                              IF( ITRAN.EQ.1 ) THEN
558                                 RCONDC = RCONDO
559                                 NORM = 'O'
560                              ELSE
561                                 RCONDC = RCONDI
562                                 NORM = 'I'
563                              END IF
564*
565*+    TEST 2:
566*                             Solve and compute residual for op(A) * X = B.
567*
568                              SRNAMT = 'SLARHS'
569                              CALL SLARHS( PATH, XTYPE, ' ', TRANS, N,
570     $                                     N, KL, KU, NRHS, A, LDA,
571     $                                     XACT, LDB, B, LDB, ISEED,
572     $                                     INFO )
573                              XTYPE = 'C'
574                              CALL SLACPY( 'Full', N, NRHS, B, LDB, X,
575     $                                     LDB )
576*
577                              SRNAMT = 'SGBTRS'
578                              CALL SGBTRS( TRANS, N, KL, KU, NRHS, AFAC,
579     $                                     LDAFAC, IWORK, X, LDB, INFO )
580*
581*                             Check error code from SGBTRS.
582*
583                              IF( INFO.NE.0 )
584     $                           CALL ALAERH( PATH, 'SGBTRS', INFO, 0,
585     $                                        TRANS, N, N, KL, KU, -1,
586     $                                        IMAT, NFAIL, NERRS, NOUT )
587*
588                              CALL SLACPY( 'Full', N, NRHS, B, LDB,
589     $                                     WORK, LDB )
590                              CALL SGBT02( TRANS, M, N, KL, KU, NRHS, A,
591     $                                     LDA, X, LDB, WORK, LDB,
592     $                                     RWORK, RESULT( 2 ) )
593*
594*+    TEST 3:
595*                             Check solution from generated exact
596*                             solution.
597*
598                              CALL SGET04( N, NRHS, X, LDB, XACT, LDB,
599     $                                     RCONDC, RESULT( 3 ) )
600*
601*+    TESTS 4, 5, 6:
602*                             Use iterative refinement to improve the
603*                             solution.
604*
605                              SRNAMT = 'SGBRFS'
606                              CALL SGBRFS( TRANS, N, KL, KU, NRHS, A,
607     $                                     LDA, AFAC, LDAFAC, IWORK, B,
608     $                                     LDB, X, LDB, RWORK,
609     $                                     RWORK( NRHS+1 ), WORK,
610     $                                     IWORK( N+1 ), INFO )
611*
612*                             Check error code from SGBRFS.
613*
614                              IF( INFO.NE.0 )
615     $                           CALL ALAERH( PATH, 'SGBRFS', INFO, 0,
616     $                                        TRANS, N, N, KL, KU, NRHS,
617     $                                        IMAT, NFAIL, NERRS, NOUT )
618*
619                              CALL SGET04( N, NRHS, X, LDB, XACT, LDB,
620     $                                     RCONDC, RESULT( 4 ) )
621                              CALL SGBT05( TRANS, N, KL, KU, NRHS, A,
622     $                                     LDA, B, LDB, X, LDB, XACT,
623     $                                     LDB, RWORK, RWORK( NRHS+1 ),
624     $                                     RESULT( 5 ) )
625                              DO 60 K = 2, 6
626                                 IF( RESULT( K ).GE.THRESH ) THEN
627                                    IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
628     $                                 CALL ALAHD( NOUT, PATH )
629                                    WRITE( NOUT, FMT = 9996 )TRANS, N,
630     $                                 KL, KU, NRHS, IMAT, K,
631     $                                 RESULT( K )
632                                    NFAIL = NFAIL + 1
633                                 END IF
634   60                         CONTINUE
635                              NRUN = NRUN + 5
636   70                      CONTINUE
637   80                   CONTINUE
638*
639*+    TEST 7:
640*                          Get an estimate of RCOND = 1/CNDNUM.
641*
642   90                   CONTINUE
643                        DO 100 ITRAN = 1, 2
644                           IF( ITRAN.EQ.1 ) THEN
645                              ANORM = ANORMO
646                              RCONDC = RCONDO
647                              NORM = 'O'
648                           ELSE
649                              ANORM = ANORMI
650                              RCONDC = RCONDI
651                              NORM = 'I'
652                           END IF
653                           SRNAMT = 'SGBCON'
654                           CALL SGBCON( NORM, N, KL, KU, AFAC, LDAFAC,
655     $                                  IWORK, ANORM, RCOND, WORK,
656     $                                  IWORK( N+1 ), INFO )
657*
658*                             Check error code from SGBCON.
659*
660                           IF( INFO.NE.0 )
661     $                        CALL ALAERH( PATH, 'SGBCON', INFO, 0,
662     $                                     NORM, N, N, KL, KU, -1, IMAT,
663     $                                     NFAIL, NERRS, NOUT )
664*
665                           RESULT( 7 ) = SGET06( RCOND, RCONDC )
666*
667*                          Print information about the tests that did
668*                          not pass the threshold.
669*
670                           IF( RESULT( 7 ).GE.THRESH ) THEN
671                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
672     $                           CALL ALAHD( NOUT, PATH )
673                              WRITE( NOUT, FMT = 9995 )NORM, N, KL, KU,
674     $                           IMAT, 7, RESULT( 7 )
675                              NFAIL = NFAIL + 1
676                           END IF
677                           NRUN = NRUN + 1
678  100                   CONTINUE
679*
680  110                CONTINUE
681  120             CONTINUE
682  130          CONTINUE
683  140       CONTINUE
684  150    CONTINUE
685  160 CONTINUE
686*
687*     Print a summary of the results.
688*
689      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
690*
691 9999 FORMAT( ' *** In SCHKGB, LA=', I5, ' is too small for M=', I5,
692     $      ', N=', I5, ', KL=', I4, ', KU=', I4,
693     $      / ' ==> Increase LA to at least ', I5 )
694 9998 FORMAT( ' *** In SCHKGB, LAFAC=', I5, ' is too small for M=', I5,
695     $      ', N=', I5, ', KL=', I4, ', KU=', I4,
696     $      / ' ==> Increase LAFAC to at least ', I5 )
697 9997 FORMAT( ' M =', I5, ', N =', I5, ', KL=', I5, ', KU=', I5,
698     $      ', NB =', I4, ', type ', I1, ', test(', I1, ')=', G12.5 )
699 9996 FORMAT( ' TRANS=''', A1, ''', N=', I5, ', KL=', I5, ', KU=', I5,
700     $      ', NRHS=', I3, ', type ', I1, ', test(', I1, ')=', G12.5 )
701 9995 FORMAT( ' NORM =''', A1, ''', N=', I5, ', KL=', I5, ', KU=', I5,
702     $      ',', 10X, ' type ', I1, ', test(', I1, ')=', G12.5 )
703*
704      RETURN
705*
706*     End of SCHKGB
707*
708      END
709