LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dchkbb.f
Go to the documentation of this file.
1*> \brief \b DCHKBB
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 DCHKBB( 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* DOUBLE PRECISION THRESH
20* ..
21* .. Array Arguments ..
22* LOGICAL DOTYPE( * )
23* INTEGER ISEED( 4 ), KK( * ), MVAL( * ), NVAL( * )
24* DOUBLE PRECISION 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*> DCHKBB tests the reduction of a general real rectangular band
36*> matrix to bidiagonal form.
37*>
38*> DGBBRD 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*> DGBBRD 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, DCHKBB 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*> DCHKBB 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, DCHKBB
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 DCHKBB to continue the same random number
181*> sequence.
182*> \endverbatim
183*>
184*> \param[in] THRESH
185*> \verbatim
186*> THRESH is DOUBLE PRECISION
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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (max(NN))
232*> Used to hold the diagonal of the bidiagonal matrix computed
233*> by DGBBRD.
234*> \endverbatim
235*>
236*> \param[out] BE
237*> \verbatim
238*> BE is DOUBLE PRECISION array, dimension (max(NN))
239*> Used to hold the off-diagonal of the bidiagonal matrix
240*> computed by DGBBRD.
241*> \endverbatim
242*>
243*> \param[out] Q
244*> \verbatim
245*> Q is DOUBLE PRECISION array, dimension (LDQ, max(NN))
246*> Used to hold the orthogonal matrix Q computed by DGBBRD.
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 DOUBLE PRECISION array, dimension (LDP, max(NN))
259*> Used to hold the orthogonal matrix P computed by DGBBRD.
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 DOUBLE PRECISION array, dimension (LDC, max(NN))
272*> Used to hold the matrix C updated by DGBBRD.
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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 double_eig
349*
350* =====================================================================
351 SUBROUTINE dchkbb( 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 DOUBLE PRECISION THRESH
364* ..
365* .. Array Arguments ..
366 LOGICAL DOTYPE( * )
367 INTEGER ISEED( 4 ), KK( * ), MVAL( * ), NVAL( * )
368 DOUBLE PRECISION 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 DOUBLE PRECISION ZERO, ONE
377 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
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 DOUBLE PRECISION 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 DOUBLE PRECISION DLAMCH
396 EXTERNAL DLAMCH
397* ..
398* .. External Subroutines ..
399 EXTERNAL dbdt01, dbdt02, dgbbrd, dlacpy, dlahd2, dlaset,
401* ..
402* .. Intrinsic Functions ..
403 INTRINSIC abs, dble, max, min, 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( 'DCHKBB', -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 = dlamch( 'Safe minimum' )
486 ovfl = one / unfl
487 ulp = dlamch( 'Epsilon' )*dlamch( '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 / dble( 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 dlaset( 'Full', lda, n, zero, zero, a, lda )
566 CALL dlaset( '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 dlatms( 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 dlatms( 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 dlatmr( 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 dlatmr( 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 dlacpy( 'Full', m, nrhs, c, ldc, cc, ldc )
643*
644* Call DGBBRD to compute B, Q and P, and to update C.
645*
646 CALL dgbbrd( '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 )'DGBBRD', 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 dbdt01( m, n, -1, a, lda, q, ldq, bd, be, p, ldp,
667 $ work, result( 1 ) )
668 CALL dort01( 'Columns', m, m, q, ldq, work, lwork,
669 $ result( 2 ) )
670 CALL dort01( 'Rows', n, n, p, ldp, work, lwork,
671 $ result( 3 ) )
672 CALL dbdt02( 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 dlahd2( nounit, 'DBB' )
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 dlasum( 'DBB', nounit, nerrs, ntestt )
700 RETURN
701*
702 9999 FORMAT( ' DCHKBB: ', 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 DCHKBB
709*
710 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dbdt01(m, n, kd, a, lda, q, ldq, d, e, pt, ldpt, work, resid)
DBDT01
Definition dbdt01.f:141
subroutine dbdt02(m, n, b, ldb, c, ldc, u, ldu, work, resid)
DBDT02
Definition dbdt02.f:112
subroutine dchkbb(nsizes, mval, nval, nwdths, kk, ntypes, dotype, nrhs, iseed, thresh, nounit, a, lda, ab, ldab, bd, be, q, ldq, p, ldp, c, ldc, cc, work, lwork, result, info)
DCHKBB
Definition dchkbb.f:355
subroutine dlahd2(iounit, path)
DLAHD2
Definition dlahd2.f:65
subroutine dlasum(type, iounit, ie, nrun)
DLASUM
Definition dlasum.f:43
subroutine dlatmr(m, n, dist, iseed, sym, d, mode, cond, dmax, rsign, grade, dl, model, condl, dr, moder, condr, pivtng, ipivot, kl, ku, sparse, anorm, pack, a, lda, iwork, info)
DLATMR
Definition dlatmr.f:471
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dort01(rowcol, m, n, u, ldu, work, lwork, resid)
DORT01
Definition dort01.f:116
subroutine dgbbrd(vect, m, n, ncc, kl, ku, ab, ldab, d, e, q, ldq, pt, ldpt, c, ldc, work, info)
DGBBRD
Definition dgbbrd.f:187
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:110