1*> \brief \b DCHKGB
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 DCHKGB( 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*       DOUBLE PRECISION   THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
23*      $                   NVAL( * )
24*       DOUBLE PRECISION   A( * ), AFAC( * ), B( * ), RWORK( * ),
25*      $                   WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> DCHKGB tests DGBTRF, -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 DOUBLE PRECISION
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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (NMAX*NSMAX)
147*> \endverbatim
148*>
149*> \param[out] XACT
150*> \verbatim
151*>          XACT is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
152*> \endverbatim
153*>
154*> \param[out] WORK
155*> \verbatim
156*>          WORK is DOUBLE PRECISION array, dimension
157*>                      (NMAX*max(3,NSMAX,NMAX))
158*> \endverbatim
159*>
160*> \param[out] RWORK
161*> \verbatim
162*>          RWORK is DOUBLE PRECISION array, dimension
163*>                      (max(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*> \date November 2011
186*
187*> \ingroup double_lin
188*
189*  =====================================================================
190      SUBROUTINE DCHKGB( 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.4.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*     November 2011
198*
199*     .. Scalar Arguments ..
200      LOGICAL            TSTERR
201      INTEGER            LA, LAFAC, NM, NN, NNB, NNS, NOUT
202      DOUBLE PRECISION   THRESH
203*     ..
204*     .. Array Arguments ..
205      LOGICAL            DOTYPE( * )
206      INTEGER            IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
207     $                   NVAL( * )
208      DOUBLE PRECISION   A( * ), AFAC( * ), B( * ), RWORK( * ),
209     $                   WORK( * ), X( * ), XACT( * )
210*     ..
211*
212*  =====================================================================
213*
214*     .. Parameters ..
215      DOUBLE PRECISION   ONE, ZERO
216      PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
217      INTEGER            NTYPES, NTESTS
218      PARAMETER          ( NTYPES = 8, NTESTS = 7 )
219      INTEGER            NBW, NTRAN
220      PARAMETER          ( NBW = 4, NTRAN = 3 )
221*     ..
222*     .. Local Scalars ..
223      LOGICAL            TRFCON, ZEROT
224      CHARACTER          DIST, NORM, TRANS, TYPE, XTYPE
225      CHARACTER*3        PATH
226      INTEGER            I, I1, I2, IKL, IKU, IM, IMAT, IN, INB, INFO,
227     $                   IOFF, IRHS, ITRAN, IZERO, J, K, KL, KOFF, KU,
228     $                   LDA, LDAFAC, LDB, M, MODE, N, NB, NERRS, NFAIL,
229     $                   NIMAT, NKL, NKU, NRHS, NRUN
230      DOUBLE PRECISION   AINVNM, ANORM, ANORMI, ANORMO, CNDNUM, RCOND,
231     $                   RCONDC, RCONDI, RCONDO
232*     ..
233*     .. Local Arrays ..
234      CHARACTER          TRANSS( NTRAN )
235      INTEGER            ISEED( 4 ), ISEEDY( 4 ), KLVAL( NBW ),
236     $                   KUVAL( NBW )
237      DOUBLE PRECISION   RESULT( NTESTS )
238*     ..
239*     .. External Functions ..
240      DOUBLE PRECISION   DGET06, DLANGB, DLANGE
241      EXTERNAL           DGET06, DLANGB, DLANGE
242*     ..
243*     .. External Subroutines ..
244      EXTERNAL           ALAERH, ALAHD, ALASUM, DCOPY, DERRGE, DGBCON,
245     $                   DGBRFS, DGBT01, DGBT02, DGBT05, DGBTRF, DGBTRS,
246     $                   DGET04, DLACPY, DLARHS, DLASET, DLATB4, DLATMS,
247     $                   XLAENV
248*     ..
249*     .. Intrinsic Functions ..
250      INTRINSIC          MAX, MIN
251*     ..
252*     .. Scalars in Common ..
253      LOGICAL            LERR, OK
254      CHARACTER*32       SRNAMT
255      INTEGER            INFOT, NUNIT
256*     ..
257*     .. Common blocks ..
258      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
259      COMMON             / SRNAMC / SRNAMT
260*     ..
261*     .. Data statements ..
262      DATA               ISEEDY / 1988, 1989, 1990, 1991 / ,
263     $                   TRANSS / 'N', 'T', 'C' /
264*     ..
265*     .. Executable Statements ..
266*
267*     Initialize constants and the random number seed.
268*
269      PATH( 1: 1 ) = 'Double precision'
270      PATH( 2: 3 ) = 'GB'
271      NRUN = 0
272      NFAIL = 0
273      NERRS = 0
274      DO 10 I = 1, 4
275         ISEED( I ) = ISEEDY( I )
276   10 CONTINUE
277*
278*     Test the error exits
279*
280      IF( TSTERR )
281     $   CALL DERRGE( PATH, NOUT )
282      INFOT = 0
283      CALL XLAENV( 2, 2 )
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 DLATB4 and generate a
384*                       test matrix with DLATMS.
385*
386                        CALL DLATB4( 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 = 'DLATMS'
394                        CALL DLATMS( 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 DLATMS.
399*
400                        IF( INFO.NE.0 ) THEN
401                           CALL ALAERH( PATH, 'DLATMS', 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 DCOPY( 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 DCOPY( 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 = DLANGB( 'O', N, KL, KU, A, LDA, RWORK )
454*                     ANORMI = DLANGB( '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 DLACPY( 'Full', KL+KU+1, N, A, LDA,
466     $                                  AFAC( KL+1 ), LDAFAC )
467                        SRNAMT = 'DGBTRF'
468                        CALL DGBTRF( M, N, KL, KU, AFAC, LDAFAC, IWORK,
469     $                               INFO )
470*
471*                       Check error code from DGBTRF.
472*
473                        IF( INFO.NE.IZERO )
474     $                     CALL ALAERH( PATH, 'DGBTRF', 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 DGBT01( 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 = DLANGB( 'O', N, KL, KU, A, LDA, RWORK )
505                        ANORMI = DLANGB( '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 DLASET( 'Full', N, N, ZERO, ONE, WORK,
514     $                                  LDB )
515                           SRNAMT = 'DGBTRS'
516                           CALL DGBTRS( '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 = DLANGE( '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 = DLANGE( '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 = 'DLARHS'
572                              CALL DLARHS( PATH, XTYPE, ' ', TRANS, N,
573     $                                     N, KL, KU, NRHS, A, LDA,
574     $                                     XACT, LDB, B, LDB, ISEED,
575     $                                     INFO )
576                              XTYPE = 'C'
577                              CALL DLACPY( 'Full', N, NRHS, B, LDB, X,
578     $                                     LDB )
579*
580                              SRNAMT = 'DGBTRS'
581                              CALL DGBTRS( TRANS, N, KL, KU, NRHS, AFAC,
582     $                                     LDAFAC, IWORK, X, LDB, INFO )
583*
584*                             Check error code from DGBTRS.
585*
586                              IF( INFO.NE.0 )
587     $                           CALL ALAERH( PATH, 'DGBTRS', INFO, 0,
588     $                                        TRANS, N, N, KL, KU, -1,
589     $                                        IMAT, NFAIL, NERRS, NOUT )
590*
591                              CALL DLACPY( 'Full', N, NRHS, B, LDB,
592     $                                     WORK, LDB )
593                              CALL DGBT02( 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 DGET04( 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 = 'DGBRFS'
609                              CALL DGBRFS( TRANS, N, KL, KU, NRHS, A,
610     $                                     LDA, AFAC, LDAFAC, IWORK, B,
611     $                                     LDB, X, LDB, RWORK,
612     $                                     RWORK( NRHS+1 ), WORK,
613     $                                     IWORK( N+1 ), INFO )
614*
615*                             Check error code from DGBRFS.
616*
617                              IF( INFO.NE.0 )
618     $                           CALL ALAERH( PATH, 'DGBRFS', INFO, 0,
619     $                                        TRANS, N, N, KL, KU, NRHS,
620     $                                        IMAT, NFAIL, NERRS, NOUT )
621*
622                              CALL DGET04( N, NRHS, X, LDB, XACT, LDB,
623     $                                     RCONDC, RESULT( 4 ) )
624                              CALL DGBT05( TRANS, N, KL, KU, NRHS, A,
625     $                                     LDA, B, LDB, X, LDB, XACT,
626     $                                     LDB, RWORK, RWORK( NRHS+1 ),
627     $                                     RESULT( 5 ) )
628                              DO 60 K = 2, 6
629                                 IF( RESULT( K ).GE.THRESH ) THEN
630                                    IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
631     $                                 CALL ALAHD( NOUT, PATH )
632                                    WRITE( NOUT, FMT = 9996 )TRANS, N,
633     $                                 KL, KU, NRHS, IMAT, K,
634     $                                 RESULT( K )
635                                    NFAIL = NFAIL + 1
636                                 END IF
637   60                         CONTINUE
638                              NRUN = NRUN + 5
639   70                      CONTINUE
640   80                   CONTINUE
641*
642*+    TEST 7:
643*                          Get an estimate of RCOND = 1/CNDNUM.
644*
645   90                   CONTINUE
646                        DO 100 ITRAN = 1, 2
647                           IF( ITRAN.EQ.1 ) THEN
648                              ANORM = ANORMO
649                              RCONDC = RCONDO
650                              NORM = 'O'
651                           ELSE
652                              ANORM = ANORMI
653                              RCONDC = RCONDI
654                              NORM = 'I'
655                           END IF
656                           SRNAMT = 'DGBCON'
657                           CALL DGBCON( NORM, N, KL, KU, AFAC, LDAFAC,
658     $                                  IWORK, ANORM, RCOND, WORK,
659     $                                  IWORK( N+1 ), INFO )
660*
661*                             Check error code from DGBCON.
662*
663                           IF( INFO.NE.0 )
664     $                        CALL ALAERH( PATH, 'DGBCON', INFO, 0,
665     $                                     NORM, N, N, KL, KU, -1, IMAT,
666     $                                     NFAIL, NERRS, NOUT )
667*
668                           RESULT( 7 ) = DGET06( RCOND, RCONDC )
669*
670*                          Print information about the tests that did
671*                          not pass the threshold.
672*
673                           IF( RESULT( 7 ).GE.THRESH ) THEN
674                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
675     $                           CALL ALAHD( NOUT, PATH )
676                              WRITE( NOUT, FMT = 9995 )NORM, N, KL, KU,
677     $                           IMAT, 7, RESULT( 7 )
678                              NFAIL = NFAIL + 1
679                           END IF
680                           NRUN = NRUN + 1
681  100                   CONTINUE
682*
683  110                CONTINUE
684  120             CONTINUE
685  130          CONTINUE
686  140       CONTINUE
687  150    CONTINUE
688  160 CONTINUE
689*
690*     Print a summary of the results.
691*
692      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
693*
694 9999 FORMAT( ' *** In DCHKGB, LA=', I5, ' is too small for M=', I5,
695     $      ', N=', I5, ', KL=', I4, ', KU=', I4,
696     $      / ' ==> Increase LA to at least ', I5 )
697 9998 FORMAT( ' *** In DCHKGB, LAFAC=', I5, ' is too small for M=', I5,
698     $      ', N=', I5, ', KL=', I4, ', KU=', I4,
699     $      / ' ==> Increase LAFAC to at least ', I5 )
700 9997 FORMAT( ' M =', I5, ', N =', I5, ', KL=', I5, ', KU=', I5,
701     $      ', NB =', I4, ', type ', I1, ', test(', I1, ')=', G12.5 )
702 9996 FORMAT( ' TRANS=''', A1, ''', N=', I5, ', KL=', I5, ', KU=', I5,
703     $      ', NRHS=', I3, ', type ', I1, ', test(', I1, ')=', G12.5 )
704 9995 FORMAT( ' NORM =''', A1, ''', N=', I5, ', KL=', I5, ', KU=', I5,
705     $      ',', 10X, ' type ', I1, ', test(', I1, ')=', G12.5 )
706*
707      RETURN
708*
709*     End of DCHKGB
710*
711      END
712