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