1*> \brief \b SCHKBB
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 SCHKBB( NSIZES, MVAL, NVAL, NWDTHS, KK, NTYPES, DOTYPE,
12*                          NRHS, ISEED, THRESH, NOUNIT, A, LDA, AB, LDAB,
13*                          BD, BE, Q, LDQ, P, LDP, C, LDC, CC, WORK,
14*                          LWORK, RESULT, INFO )
15*
16*       .. Scalar Arguments ..
17*       INTEGER            INFO, LDA, LDAB, LDC, LDP, LDQ, LWORK, NOUNIT,
18*      $                   NRHS, NSIZES, NTYPES, NWDTHS
19*       REAL               THRESH
20*       ..
21*       .. Array Arguments ..
22*       LOGICAL            DOTYPE( * )
23*       INTEGER            ISEED( 4 ), KK( * ), MVAL( * ), NVAL( * )
24*       REAL               A( LDA, * ), AB( LDAB, * ), BD( * ), BE( * ),
25*      $                   C( LDC, * ), CC( LDC, * ), P( LDP, * ),
26*      $                   Q( LDQ, * ), RESULT( * ), WORK( * )
27*       ..
28*
29*
30*> \par Purpose:
31*  =============
32*>
33*> \verbatim
34*>
35*> SCHKBB tests the reduction of a general real rectangular band
36*> matrix to bidiagonal form.
37*>
38*> SGBBRD factors a general band matrix A as  Q B P* , where * means
39*> transpose, B is upper bidiagonal, and Q and P are orthogonal;
40*> SGBBRD can also overwrite a given matrix C with Q* C .
41*>
42*> For each pair of matrix dimensions (M,N) and each selected matrix
43*> type, an M by N matrix A and an M by NRHS matrix C are generated.
44*> The problem dimensions are as follows
45*>    A:          M x N
46*>    Q:          M x M
47*>    P:          N x N
48*>    B:          min(M,N) x min(M,N)
49*>    C:          M x NRHS
50*>
51*> For each generated matrix, 4 tests are performed:
52*>
53*> (1)   | A - Q B PT | / ( |A| max(M,N) ulp ), PT = P'
54*>
55*> (2)   | I - Q' Q | / ( M ulp )
56*>
57*> (3)   | I - PT PT' | / ( N ulp )
58*>
59*> (4)   | Y - Q' C | / ( |Y| max(M,NRHS) ulp ), where Y = Q' C.
60*>
61*> The "types" are specified by a logical array DOTYPE( 1:NTYPES );
62*> if DOTYPE(j) is .TRUE., then matrix type "j" will be generated.
63*> Currently, the list of possible types is:
64*>
65*> The possible matrix types are
66*>
67*> (1)  The zero matrix.
68*> (2)  The identity matrix.
69*>
70*> (3)  A diagonal matrix with evenly spaced entries
71*>      1, ..., ULP  and random signs.
72*>      (ULP = (first number larger than 1) - 1 )
73*> (4)  A diagonal matrix with geometrically spaced entries
74*>      1, ..., ULP  and random signs.
75*> (5)  A diagonal matrix with "clustered" entries 1, ULP, ..., ULP
76*>      and random signs.
77*>
78*> (6)  Same as (3), but multiplied by SQRT( overflow threshold )
79*> (7)  Same as (3), but multiplied by SQRT( underflow threshold )
80*>
81*> (8)  A matrix of the form  U D V, where U and V are orthogonal and
82*>      D has evenly spaced entries 1, ..., ULP with random signs
83*>      on the diagonal.
84*>
85*> (9)  A matrix of the form  U D V, where U and V are orthogonal and
86*>      D has geometrically spaced entries 1, ..., ULP with random
87*>      signs on the diagonal.
88*>
89*> (10) A matrix of the form  U D V, where U and V are orthogonal and
90*>      D has "clustered" entries 1, ULP,..., ULP with random
91*>      signs on the diagonal.
92*>
93*> (11) Same as (8), but multiplied by SQRT( overflow threshold )
94*> (12) Same as (8), but multiplied by SQRT( underflow threshold )
95*>
96*> (13) Rectangular matrix with random entries chosen from (-1,1).
97*> (14) Same as (13), but multiplied by SQRT( overflow threshold )
98*> (15) Same as (13), but multiplied by SQRT( underflow threshold )
99*> \endverbatim
100*
101*  Arguments:
102*  ==========
103*
104*> \param[in] NSIZES
105*> \verbatim
106*>          NSIZES is INTEGER
107*>          The number of values of M and N contained in the vectors
108*>          MVAL and NVAL.  The matrix sizes are used in pairs (M,N).
109*>          If NSIZES is zero, SCHKBB does nothing.  NSIZES must be at
110*>          least zero.
111*> \endverbatim
112*>
113*> \param[in] MVAL
114*> \verbatim
115*>          MVAL is INTEGER array, dimension (NSIZES)
116*>          The values of the matrix row dimension M.
117*> \endverbatim
118*>
119*> \param[in] NVAL
120*> \verbatim
121*>          NVAL is INTEGER array, dimension (NSIZES)
122*>          The values of the matrix column dimension N.
123*> \endverbatim
124*>
125*> \param[in] NWDTHS
126*> \verbatim
127*>          NWDTHS is INTEGER
128*>          The number of bandwidths to use.  If it is zero,
129*>          SCHKBB does nothing.  It must be at least zero.
130*> \endverbatim
131*>
132*> \param[in] KK
133*> \verbatim
134*>          KK is INTEGER array, dimension (NWDTHS)
135*>          An array containing the bandwidths to be used for the band
136*>          matrices.  The values must be at least zero.
137*> \endverbatim
138*>
139*> \param[in] NTYPES
140*> \verbatim
141*>          NTYPES is INTEGER
142*>          The number of elements in DOTYPE.   If it is zero, SCHKBB
143*>          does nothing.  It must be at least zero.  If it is MAXTYP+1
144*>          and NSIZES is 1, then an additional type, MAXTYP+1 is
145*>          defined, which is to use whatever matrix is in A.  This
146*>          is only useful if DOTYPE(1:MAXTYP) is .FALSE. and
147*>          DOTYPE(MAXTYP+1) is .TRUE. .
148*> \endverbatim
149*>
150*> \param[in] DOTYPE
151*> \verbatim
152*>          DOTYPE is LOGICAL array, dimension (NTYPES)
153*>          If DOTYPE(j) is .TRUE., then for each size in NN a
154*>          matrix of that size and of type j will be generated.
155*>          If NTYPES is smaller than the maximum number of types
156*>          defined (PARAMETER MAXTYP), then types NTYPES+1 through
157*>          MAXTYP will not be generated.  If NTYPES is larger
158*>          than MAXTYP, DOTYPE(MAXTYP+1) through DOTYPE(NTYPES)
159*>          will be ignored.
160*> \endverbatim
161*>
162*> \param[in] NRHS
163*> \verbatim
164*>          NRHS is INTEGER
165*>          The number of columns in the "right-hand side" matrix C.
166*>          If NRHS = 0, then the operations on the right-hand side will
167*>          not be tested. NRHS must be at least 0.
168*> \endverbatim
169*>
170*> \param[in,out] ISEED
171*> \verbatim
172*>          ISEED is INTEGER array, dimension (4)
173*>          On entry ISEED specifies the seed of the random number
174*>          generator. The array elements should be between 0 and 4095;
175*>          if not they will be reduced mod 4096.  Also, ISEED(4) must
176*>          be odd.  The random number generator uses a linear
177*>          congruential sequence limited to small integers, and so
178*>          should produce machine independent random numbers. The
179*>          values of ISEED are changed on exit, and can be used in the
180*>          next call to SCHKBB to continue the same random number
181*>          sequence.
182*> \endverbatim
183*>
184*> \param[in] THRESH
185*> \verbatim
186*>          THRESH is REAL
187*>          A test will count as "failed" if the "error", computed as
188*>          described above, exceeds THRESH.  Note that the error
189*>          is scaled to be O(1), so THRESH should be a reasonably
190*>          small multiple of 1, e.g., 10 or 100.  In particular,
191*>          it should not depend on the precision (single vs. double)
192*>          or the size of the matrix.  It must be at least zero.
193*> \endverbatim
194*>
195*> \param[in] NOUNIT
196*> \verbatim
197*>          NOUNIT is INTEGER
198*>          The FORTRAN unit number for printing out error messages
199*>          (e.g., if a routine returns IINFO not equal to 0.)
200*> \endverbatim
201*>
202*> \param[in,out] A
203*> \verbatim
204*>          A is REAL array, dimension
205*>                            (LDA, max(NN))
206*>          Used to hold the matrix A.
207*> \endverbatim
208*>
209*> \param[in] LDA
210*> \verbatim
211*>          LDA is INTEGER
212*>          The leading dimension of A.  It must be at least 1
213*>          and at least max( NN ).
214*> \endverbatim
215*>
216*> \param[out] AB
217*> \verbatim
218*>          AB is REAL array, dimension (LDAB, max(NN))
219*>          Used to hold A in band storage format.
220*> \endverbatim
221*>
222*> \param[in] LDAB
223*> \verbatim
224*>          LDAB is INTEGER
225*>          The leading dimension of AB.  It must be at least 2 (not 1!)
226*>          and at least max( KK )+1.
227*> \endverbatim
228*>
229*> \param[out] BD
230*> \verbatim
231*>          BD is REAL array, dimension (max(NN))
232*>          Used to hold the diagonal of the bidiagonal matrix computed
233*>          by SGBBRD.
234*> \endverbatim
235*>
236*> \param[out] BE
237*> \verbatim
238*>          BE is REAL array, dimension (max(NN))
239*>          Used to hold the off-diagonal of the bidiagonal matrix
240*>          computed by SGBBRD.
241*> \endverbatim
242*>
243*> \param[out] Q
244*> \verbatim
245*>          Q is REAL array, dimension (LDQ, max(NN))
246*>          Used to hold the orthogonal matrix Q computed by SGBBRD.
247*> \endverbatim
248*>
249*> \param[in] LDQ
250*> \verbatim
251*>          LDQ is INTEGER
252*>          The leading dimension of Q.  It must be at least 1
253*>          and at least max( NN ).
254*> \endverbatim
255*>
256*> \param[out] P
257*> \verbatim
258*>          P is REAL array, dimension (LDP, max(NN))
259*>          Used to hold the orthogonal matrix P computed by SGBBRD.
260*> \endverbatim
261*>
262*> \param[in] LDP
263*> \verbatim
264*>          LDP is INTEGER
265*>          The leading dimension of P.  It must be at least 1
266*>          and at least max( NN ).
267*> \endverbatim
268*>
269*> \param[out] C
270*> \verbatim
271*>          C is REAL array, dimension (LDC, max(NN))
272*>          Used to hold the matrix C updated by SGBBRD.
273*> \endverbatim
274*>
275*> \param[in] LDC
276*> \verbatim
277*>          LDC is INTEGER
278*>          The leading dimension of U.  It must be at least 1
279*>          and at least max( NN ).
280*> \endverbatim
281*>
282*> \param[out] CC
283*> \verbatim
284*>          CC is REAL array, dimension (LDC, max(NN))
285*>          Used to hold a copy of the matrix C.
286*> \endverbatim
287*>
288*> \param[out] WORK
289*> \verbatim
290*>          WORK is REAL array, dimension (LWORK)
291*> \endverbatim
292*>
293*> \param[in] LWORK
294*> \verbatim
295*>          LWORK is INTEGER
296*>          The number of entries in WORK.  This must be at least
297*>          max( LDA+1, max(NN)+1 )*max(NN).
298*> \endverbatim
299*>
300*> \param[out] RESULT
301*> \verbatim
302*>          RESULT is REAL array, dimension (4)
303*>          The values computed by the tests described above.
304*>          The values are currently limited to 1/ulp, to avoid
305*>          overflow.
306*> \endverbatim
307*>
308*> \param[out] INFO
309*> \verbatim
310*>          INFO is INTEGER
311*>          If 0, then everything ran OK.
312*>
313*>-----------------------------------------------------------------------
314*>
315*>       Some Local Variables and Parameters:
316*>       ---- ----- --------- --- ----------
317*>       ZERO, ONE       Real 0 and 1.
318*>       MAXTYP          The number of types defined.
319*>       NTEST           The number of tests performed, or which can
320*>                       be performed so far, for the current matrix.
321*>       NTESTT          The total number of tests performed so far.
322*>       NMAX            Largest value in NN.
323*>       NMATS           The number of matrices generated so far.
324*>       NERRS           The number of tests which have exceeded THRESH
325*>                       so far.
326*>       COND, IMODE     Values to be passed to the matrix generators.
327*>       ANORM           Norm of A; passed to matrix generators.
328*>
329*>       OVFL, UNFL      Overflow and underflow thresholds.
330*>       ULP, ULPINV     Finest relative precision and its inverse.
331*>       RTOVFL, RTUNFL  Square roots of the previous 2 values.
332*>               The following four arrays decode JTYPE:
333*>       KTYPE(j)        The general type (1-10) for type "j".
334*>       KMODE(j)        The MODE value to be passed to the matrix
335*>                       generator for type "j".
336*>       KMAGN(j)        The order of magnitude ( O(1),
337*>                       O(overflow^(1/2) ), O(underflow^(1/2) )
338*> \endverbatim
339*
340*  Authors:
341*  ========
342*
343*> \author Univ. of Tennessee
344*> \author Univ. of California Berkeley
345*> \author Univ. of Colorado Denver
346*> \author NAG Ltd.
347*
348*> \ingroup single_eig
349*
350*  =====================================================================
351      SUBROUTINE SCHKBB( NSIZES, MVAL, NVAL, NWDTHS, KK, NTYPES, DOTYPE,
352     $                   NRHS, ISEED, THRESH, NOUNIT, A, LDA, AB, LDAB,
353     $                   BD, BE, Q, LDQ, P, LDP, C, LDC, CC, WORK,
354     $                   LWORK, RESULT, INFO )
355*
356*  -- LAPACK test routine (input) --
357*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
358*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
359*
360*     .. Scalar Arguments ..
361      INTEGER            INFO, LDA, LDAB, LDC, LDP, LDQ, LWORK, NOUNIT,
362     $                   NRHS, NSIZES, NTYPES, NWDTHS
363      REAL               THRESH
364*     ..
365*     .. Array Arguments ..
366      LOGICAL            DOTYPE( * )
367      INTEGER            ISEED( 4 ), KK( * ), MVAL( * ), NVAL( * )
368      REAL               A( LDA, * ), AB( LDAB, * ), BD( * ), BE( * ),
369     $                   C( LDC, * ), CC( LDC, * ), P( LDP, * ),
370     $                   Q( LDQ, * ), RESULT( * ), WORK( * )
371*     ..
372*
373*  =====================================================================
374*
375*     .. Parameters ..
376      REAL               ZERO, ONE
377      PARAMETER          ( ZERO = 0.0E0, ONE = 1.0E0 )
378      INTEGER            MAXTYP
379      PARAMETER          ( MAXTYP = 15 )
380*     ..
381*     .. Local Scalars ..
382      LOGICAL            BADMM, BADNN, BADNNB
383      INTEGER            I, IINFO, IMODE, ITYPE, J, JCOL, JR, JSIZE,
384     $                   JTYPE, JWIDTH, K, KL, KMAX, KU, M, MMAX, MNMAX,
385     $                   MNMIN, MTYPES, N, NERRS, NMATS, NMAX, NTEST,
386     $                   NTESTT
387      REAL               AMNINV, ANORM, COND, OVFL, RTOVFL, RTUNFL, ULP,
388     $                   ULPINV, UNFL
389*     ..
390*     .. Local Arrays ..
391      INTEGER            IDUMMA( 1 ), IOLDSD( 4 ), KMAGN( MAXTYP ),
392     $                   KMODE( MAXTYP ), KTYPE( MAXTYP )
393*     ..
394*     .. External Functions ..
395      REAL               SLAMCH
396      EXTERNAL           SLAMCH
397*     ..
398*     .. External Subroutines ..
399      EXTERNAL           SBDT01, SBDT02, SGBBRD, SLACPY, SLAHD2, SLASET,
400     $                   SLASUM, SLATMR, SLATMS, SORT01, XERBLA
401*     ..
402*     .. Intrinsic Functions ..
403      INTRINSIC          ABS, MAX, MIN, REAL, SQRT
404*     ..
405*     .. Data statements ..
406      DATA               KTYPE / 1, 2, 5*4, 5*6, 3*9 /
407      DATA               KMAGN / 2*1, 3*1, 2, 3, 3*1, 2, 3, 1, 2, 3 /
408      DATA               KMODE / 2*0, 4, 3, 1, 4, 4, 4, 3, 1, 4, 4, 0,
409     $                   0, 0 /
410*     ..
411*     .. Executable Statements ..
412*
413*     Check for errors
414*
415      NTESTT = 0
416      INFO = 0
417*
418*     Important constants
419*
420      BADMM = .FALSE.
421      BADNN = .FALSE.
422      MMAX = 1
423      NMAX = 1
424      MNMAX = 1
425      DO 10 J = 1, NSIZES
426         MMAX = MAX( MMAX, MVAL( J ) )
427         IF( MVAL( J ).LT.0 )
428     $      BADMM = .TRUE.
429         NMAX = MAX( NMAX, NVAL( J ) )
430         IF( NVAL( J ).LT.0 )
431     $      BADNN = .TRUE.
432         MNMAX = MAX( MNMAX, MIN( MVAL( J ), NVAL( J ) ) )
433   10 CONTINUE
434*
435      BADNNB = .FALSE.
436      KMAX = 0
437      DO 20 J = 1, NWDTHS
438         KMAX = MAX( KMAX, KK( J ) )
439         IF( KK( J ).LT.0 )
440     $      BADNNB = .TRUE.
441   20 CONTINUE
442*
443*     Check for errors
444*
445      IF( NSIZES.LT.0 ) THEN
446         INFO = -1
447      ELSE IF( BADMM ) THEN
448         INFO = -2
449      ELSE IF( BADNN ) THEN
450         INFO = -3
451      ELSE IF( NWDTHS.LT.0 ) THEN
452         INFO = -4
453      ELSE IF( BADNNB ) THEN
454         INFO = -5
455      ELSE IF( NTYPES.LT.0 ) THEN
456         INFO = -6
457      ELSE IF( NRHS.LT.0 ) THEN
458         INFO = -8
459      ELSE IF( LDA.LT.NMAX ) THEN
460         INFO = -13
461      ELSE IF( LDAB.LT.2*KMAX+1 ) THEN
462         INFO = -15
463      ELSE IF( LDQ.LT.NMAX ) THEN
464         INFO = -19
465      ELSE IF( LDP.LT.NMAX ) THEN
466         INFO = -21
467      ELSE IF( LDC.LT.NMAX ) THEN
468         INFO = -23
469      ELSE IF( ( MAX( LDA, NMAX )+1 )*NMAX.GT.LWORK ) THEN
470         INFO = -26
471      END IF
472*
473      IF( INFO.NE.0 ) THEN
474         CALL XERBLA( 'SCHKBB', -INFO )
475         RETURN
476      END IF
477*
478*     Quick return if possible
479*
480      IF( NSIZES.EQ.0 .OR. NTYPES.EQ.0 .OR. NWDTHS.EQ.0 )
481     $   RETURN
482*
483*     More Important constants
484*
485      UNFL = SLAMCH( 'Safe minimum' )
486      OVFL = ONE / UNFL
487      ULP = SLAMCH( 'Epsilon' )*SLAMCH( 'Base' )
488      ULPINV = ONE / ULP
489      RTUNFL = SQRT( UNFL )
490      RTOVFL = SQRT( OVFL )
491*
492*     Loop over sizes, widths, types
493*
494      NERRS = 0
495      NMATS = 0
496*
497      DO 160 JSIZE = 1, NSIZES
498         M = MVAL( JSIZE )
499         N = NVAL( JSIZE )
500         MNMIN = MIN( M, N )
501         AMNINV = ONE / REAL( MAX( 1, M, N ) )
502*
503         DO 150 JWIDTH = 1, NWDTHS
504            K = KK( JWIDTH )
505            IF( K.GE.M .AND. K.GE.N )
506     $         GO TO 150
507            KL = MAX( 0, MIN( M-1, K ) )
508            KU = MAX( 0, MIN( N-1, K ) )
509*
510            IF( NSIZES.NE.1 ) THEN
511               MTYPES = MIN( MAXTYP, NTYPES )
512            ELSE
513               MTYPES = MIN( MAXTYP+1, NTYPES )
514            END IF
515*
516            DO 140 JTYPE = 1, MTYPES
517               IF( .NOT.DOTYPE( JTYPE ) )
518     $            GO TO 140
519               NMATS = NMATS + 1
520               NTEST = 0
521*
522               DO 30 J = 1, 4
523                  IOLDSD( J ) = ISEED( J )
524   30          CONTINUE
525*
526*              Compute "A".
527*
528*              Control parameters:
529*
530*                  KMAGN  KMODE        KTYPE
531*              =1  O(1)   clustered 1  zero
532*              =2  large  clustered 2  identity
533*              =3  small  exponential  (none)
534*              =4         arithmetic   diagonal, (w/ singular values)
535*              =5         random log   (none)
536*              =6         random       nonhermitian, w/ singular values
537*              =7                      (none)
538*              =8                      (none)
539*              =9                      random nonhermitian
540*
541               IF( MTYPES.GT.MAXTYP )
542     $            GO TO 90
543*
544               ITYPE = KTYPE( JTYPE )
545               IMODE = KMODE( JTYPE )
546*
547*              Compute norm
548*
549               GO TO ( 40, 50, 60 )KMAGN( JTYPE )
550*
551   40          CONTINUE
552               ANORM = ONE
553               GO TO 70
554*
555   50          CONTINUE
556               ANORM = ( RTOVFL*ULP )*AMNINV
557               GO TO 70
558*
559   60          CONTINUE
560               ANORM = RTUNFL*MAX( M, N )*ULPINV
561               GO TO 70
562*
563   70          CONTINUE
564*
565               CALL SLASET( 'Full', LDA, N, ZERO, ZERO, A, LDA )
566               CALL SLASET( 'Full', LDAB, N, ZERO, ZERO, AB, LDAB )
567               IINFO = 0
568               COND = ULPINV
569*
570*              Special Matrices -- Identity & Jordan block
571*
572*                 Zero
573*
574               IF( ITYPE.EQ.1 ) THEN
575                  IINFO = 0
576*
577               ELSE IF( ITYPE.EQ.2 ) THEN
578*
579*                 Identity
580*
581                  DO 80 JCOL = 1, N
582                     A( JCOL, JCOL ) = ANORM
583   80             CONTINUE
584*
585               ELSE IF( ITYPE.EQ.4 ) THEN
586*
587*                 Diagonal Matrix, singular values specified
588*
589                  CALL SLATMS( M, N, 'S', ISEED, 'N', WORK, IMODE, COND,
590     $                         ANORM, 0, 0, 'N', A, LDA, WORK( M+1 ),
591     $                         IINFO )
592*
593               ELSE IF( ITYPE.EQ.6 ) THEN
594*
595*                 Nonhermitian, singular values specified
596*
597                  CALL SLATMS( M, N, 'S', ISEED, 'N', WORK, IMODE, COND,
598     $                         ANORM, KL, KU, 'N', A, LDA, WORK( M+1 ),
599     $                         IINFO )
600*
601               ELSE IF( ITYPE.EQ.9 ) THEN
602*
603*                 Nonhermitian, random entries
604*
605                  CALL SLATMR( M, N, 'S', ISEED, 'N', WORK, 6, ONE, ONE,
606     $                         'T', 'N', WORK( N+1 ), 1, ONE,
607     $                         WORK( 2*N+1 ), 1, ONE, 'N', IDUMMA, KL,
608     $                         KU, ZERO, ANORM, 'N', A, LDA, IDUMMA,
609     $                         IINFO )
610*
611               ELSE
612*
613                  IINFO = 1
614               END IF
615*
616*              Generate Right-Hand Side
617*
618               CALL SLATMR( M, NRHS, 'S', ISEED, 'N', WORK, 6, ONE, ONE,
619     $                      'T', 'N', WORK( M+1 ), 1, ONE,
620     $                      WORK( 2*M+1 ), 1, ONE, 'N', IDUMMA, M, NRHS,
621     $                      ZERO, ONE, 'NO', C, LDC, IDUMMA, IINFO )
622*
623               IF( IINFO.NE.0 ) THEN
624                  WRITE( NOUNIT, FMT = 9999 )'Generator', IINFO, N,
625     $               JTYPE, IOLDSD
626                  INFO = ABS( IINFO )
627                  RETURN
628               END IF
629*
630   90          CONTINUE
631*
632*              Copy A to band storage.
633*
634               DO 110 J = 1, N
635                  DO 100 I = MAX( 1, J-KU ), MIN( M, J+KL )
636                     AB( KU+1+I-J, J ) = A( I, J )
637  100             CONTINUE
638  110          CONTINUE
639*
640*              Copy C
641*
642               CALL SLACPY( 'Full', M, NRHS, C, LDC, CC, LDC )
643*
644*              Call SGBBRD to compute B, Q and P, and to update C.
645*
646               CALL SGBBRD( 'B', M, N, NRHS, KL, KU, AB, LDAB, BD, BE,
647     $                      Q, LDQ, P, LDP, CC, LDC, WORK, IINFO )
648*
649               IF( IINFO.NE.0 ) THEN
650                  WRITE( NOUNIT, FMT = 9999 )'SGBBRD', IINFO, N, JTYPE,
651     $               IOLDSD
652                  INFO = ABS( IINFO )
653                  IF( IINFO.LT.0 ) THEN
654                     RETURN
655                  ELSE
656                     RESULT( 1 ) = ULPINV
657                     GO TO 120
658                  END IF
659               END IF
660*
661*              Test 1:  Check the decomposition A := Q * B * P'
662*                   2:  Check the orthogonality of Q
663*                   3:  Check the orthogonality of P
664*                   4:  Check the computation of Q' * C
665*
666               CALL SBDT01( M, N, -1, A, LDA, Q, LDQ, BD, BE, P, LDP,
667     $                      WORK, RESULT( 1 ) )
668               CALL SORT01( 'Columns', M, M, Q, LDQ, WORK, LWORK,
669     $                      RESULT( 2 ) )
670               CALL SORT01( 'Rows', N, N, P, LDP, WORK, LWORK,
671     $                      RESULT( 3 ) )
672               CALL SBDT02( M, NRHS, C, LDC, CC, LDC, Q, LDQ, WORK,
673     $                      RESULT( 4 ) )
674*
675*              End of Loop -- Check for RESULT(j) > THRESH
676*
677               NTEST = 4
678  120          CONTINUE
679               NTESTT = NTESTT + NTEST
680*
681*              Print out tests which fail.
682*
683               DO 130 JR = 1, NTEST
684                  IF( RESULT( JR ).GE.THRESH ) THEN
685                     IF( NERRS.EQ.0 )
686     $                  CALL SLAHD2( NOUNIT, 'SBB' )
687                     NERRS = NERRS + 1
688                     WRITE( NOUNIT, FMT = 9998 )M, N, K, IOLDSD, JTYPE,
689     $                  JR, RESULT( JR )
690                  END IF
691  130          CONTINUE
692*
693  140       CONTINUE
694  150    CONTINUE
695  160 CONTINUE
696*
697*     Summary
698*
699      CALL SLASUM( 'SBB', NOUNIT, NERRS, NTESTT )
700      RETURN
701*
702 9999 FORMAT( ' SCHKBB: ', A, ' returned INFO=', I5, '.', / 9X, 'M=',
703     $      I5, ' N=', I5, ' K=', I5, ', JTYPE=', I5, ', ISEED=(',
704     $      3( I5, ',' ), I5, ')' )
705 9998 FORMAT( ' M =', I4, ' N=', I4, ', K=', I3, ', seed=',
706     $      4( I4, ',' ), ' type ', I2, ', test(', I2, ')=', G10.3 )
707*
708*     End of SCHKBB
709*
710      END
711