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