LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cchkgb.f
Go to the documentation of this file.
1*> \brief \b CCHKGB
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE CCHKGB( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
12* NSVAL, THRESH, TSTERR, A, LA, AFAC, LAFAC, B,
13* X, XACT, WORK, RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER LA, LAFAC, NM, NN, NNB, NNS, NOUT
18* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
23* $ NVAL( * )
24* REAL RWORK( * )
25* COMPLEX A( * ), AFAC( * ), B( * ), WORK( * ), X( * ),
26* $ XACT( * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> CCHKGB tests CGBTRF, -TRS, -RFS, and -CON
36*> \endverbatim
37*
38* Arguments:
39* ==========
40*
41*> \param[in] DOTYPE
42*> \verbatim
43*> DOTYPE is LOGICAL array, dimension (NTYPES)
44*> The matrix types to be used for testing. Matrices of type j
45*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
46*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
47*> \endverbatim
48*>
49*> \param[in] NM
50*> \verbatim
51*> NM is INTEGER
52*> The number of values of M contained in the vector MVAL.
53*> \endverbatim
54*>
55*> \param[in] MVAL
56*> \verbatim
57*> MVAL is INTEGER array, dimension (NM)
58*> The values of the matrix row dimension M.
59*> \endverbatim
60*>
61*> \param[in] NN
62*> \verbatim
63*> NN is INTEGER
64*> The number of values of N contained in the vector NVAL.
65*> \endverbatim
66*>
67*> \param[in] NVAL
68*> \verbatim
69*> NVAL is INTEGER array, dimension (NN)
70*> The values of the matrix column dimension N.
71*> \endverbatim
72*>
73*> \param[in] NNB
74*> \verbatim
75*> NNB is INTEGER
76*> The number of values of NB contained in the vector NBVAL.
77*> \endverbatim
78*>
79*> \param[in] NBVAL
80*> \verbatim
81*> NBVAL is INTEGER array, dimension (NNB)
82*> The values of the blocksize NB.
83*> \endverbatim
84*>
85*> \param[in] NNS
86*> \verbatim
87*> NNS is INTEGER
88*> The number of values of NRHS contained in the vector NSVAL.
89*> \endverbatim
90*>
91*> \param[in] NSVAL
92*> \verbatim
93*> NSVAL is INTEGER array, dimension (NNS)
94*> The values of the number of right hand sides NRHS.
95*> \endverbatim
96*>
97*> \param[in] THRESH
98*> \verbatim
99*> THRESH is REAL
100*> The threshold value for the test ratios. A result is
101*> included in the output file if RESULT >= THRESH. To have
102*> every test ratio printed, use THRESH = 0.
103*> \endverbatim
104*>
105*> \param[in] TSTERR
106*> \verbatim
107*> TSTERR is LOGICAL
108*> Flag that indicates whether error exits are to be tested.
109*> \endverbatim
110*>
111*> \param[out] A
112*> \verbatim
113*> A is COMPLEX array, dimension (LA)
114*> \endverbatim
115*>
116*> \param[in] LA
117*> \verbatim
118*> LA is INTEGER
119*> The length of the array A. LA >= (KLMAX+KUMAX+1)*NMAX
120*> where KLMAX is the largest entry in the local array KLVAL,
121*> KUMAX is the largest entry in the local array KUVAL and
122*> NMAX is the largest entry in the input array NVAL.
123*> \endverbatim
124*>
125*> \param[out] AFAC
126*> \verbatim
127*> AFAC is COMPLEX array, dimension (LAFAC)
128*> \endverbatim
129*>
130*> \param[in] LAFAC
131*> \verbatim
132*> LAFAC is INTEGER
133*> The length of the array AFAC. LAFAC >= (2*KLMAX+KUMAX+1)*NMAX
134*> where KLMAX is the largest entry in the local array KLVAL,
135*> KUMAX is the largest entry in the local array KUVAL and
136*> NMAX is the largest entry in the input array NVAL.
137*> \endverbatim
138*>
139*> \param[out] B
140*> \verbatim
141*> B is COMPLEX array, dimension (NMAX*NSMAX)
142*> \endverbatim
143*>
144*> \param[out] X
145*> \verbatim
146*> X is COMPLEX array, dimension (NMAX*NSMAX)
147*> \endverbatim
148*>
149*> \param[out] XACT
150*> \verbatim
151*> XACT is COMPLEX array, dimension (NMAX*NSMAX)
152*> \endverbatim
153*>
154*> \param[out] WORK
155*> \verbatim
156*> WORK is COMPLEX array, dimension
157*> (NMAX*max(3,NSMAX,NMAX))
158*> \endverbatim
159*>
160*> \param[out] RWORK
161*> \verbatim
162*> RWORK is REAL array, dimension
163*> (NMAX+2*NSMAX)
164*> \endverbatim
165*>
166*> \param[out] IWORK
167*> \verbatim
168*> IWORK is INTEGER array, dimension (NMAX)
169*> \endverbatim
170*>
171*> \param[in] NOUT
172*> \verbatim
173*> NOUT is INTEGER
174*> The unit number for output.
175*> \endverbatim
176*
177* Authors:
178* ========
179*
180*> \author Univ. of Tennessee
181*> \author Univ. of California Berkeley
182*> \author Univ. of Colorado Denver
183*> \author NAG Ltd.
184*
185*> \ingroup complex_lin
186*
187* =====================================================================
188 SUBROUTINE cchkgb( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NNS,
189 $ NSVAL, THRESH, TSTERR, A, LA, AFAC, LAFAC, B,
190 $ X, XACT, WORK, RWORK, IWORK, NOUT )
191*
192* -- LAPACK test routine --
193* -- LAPACK is a software package provided by Univ. of Tennessee, --
194* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195*
196* .. Scalar Arguments ..
197 LOGICAL TSTERR
198 INTEGER LA, LAFAC, NM, NN, NNB, NNS, NOUT
199 REAL THRESH
200* ..
201* .. Array Arguments ..
202 LOGICAL DOTYPE( * )
203 INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NSVAL( * ),
204 $ nval( * )
205 REAL RWORK( * )
206 COMPLEX A( * ), AFAC( * ), B( * ), WORK( * ), X( * ),
207 $ xact( * )
208* ..
209*
210* =====================================================================
211*
212* .. Parameters ..
213 REAL ONE, ZERO
214 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
215 INTEGER NTYPES, NTESTS
216 parameter( ntypes = 8, ntests = 7 )
217 INTEGER NBW, NTRAN
218 parameter( nbw = 4, ntran = 3 )
219* ..
220* .. Local Scalars ..
221 LOGICAL TRFCON, ZEROT
222 CHARACTER DIST, NORM, TRANS, TYPE, XTYPE
223 CHARACTER*3 PATH
224 INTEGER I, I1, I2, IKL, IKU, IM, IMAT, IN, INB, INFO,
225 $ ioff, irhs, itran, izero, j, k, kl, koff, ku,
226 $ lda, ldafac, ldb, m, mode, n, nb, nerrs, nfail,
227 $ nimat, nkl, nku, nrhs, nrun
228 REAL AINVNM, ANORM, ANORMI, ANORMO, CNDNUM, RCOND,
229 $ RCONDC, RCONDI, RCONDO
230* ..
231* .. Local Arrays ..
232 CHARACTER TRANSS( NTRAN )
233 INTEGER ISEED( 4 ), ISEEDY( 4 ), KLVAL( NBW ),
234 $ kuval( nbw )
235 REAL RESULT( NTESTS )
236* ..
237* .. External Functions ..
238 REAL CLANGB, CLANGE, SGET06
239 EXTERNAL CLANGB, CLANGE, SGET06
240* ..
241* .. External Subroutines ..
242 EXTERNAL alaerh, alahd, alasum, ccopy, cerrge, cgbcon,
245 $ xlaenv
246* ..
247* .. Intrinsic Functions ..
248 INTRINSIC cmplx, max, min
249* ..
250* .. Scalars in Common ..
251 LOGICAL LERR, OK
252 CHARACTER*32 SRNAMT
253 INTEGER INFOT, NUNIT
254* ..
255* .. Common blocks ..
256 COMMON / infoc / infot, nunit, ok, lerr
257 COMMON / srnamc / srnamt
258* ..
259* .. Data statements ..
260 DATA iseedy / 1988, 1989, 1990, 1991 / ,
261 $ transs / 'N', 'T', 'C' /
262* ..
263* .. Executable Statements ..
264*
265* Initialize constants and the random number seed.
266*
267 path( 1: 1 ) = 'Complex precision'
268 path( 2: 3 ) = 'GB'
269 nrun = 0
270 nfail = 0
271 nerrs = 0
272 DO 10 i = 1, 4
273 iseed( i ) = iseedy( i )
274 10 CONTINUE
275*
276* Test the error exits
277*
278 IF( tsterr )
279 $ CALL cerrge( path, nout )
280 infot = 0
281*
282* Initialize the first value for the lower and upper bandwidths.
283*
284 klval( 1 ) = 0
285 kuval( 1 ) = 0
286*
287* Do for each value of M in MVAL
288*
289 DO 160 im = 1, nm
290 m = mval( im )
291*
292* Set values to use for the lower bandwidth.
293*
294 klval( 2 ) = m + ( m+1 ) / 4
295*
296* KLVAL( 2 ) = MAX( M-1, 0 )
297*
298 klval( 3 ) = ( 3*m-1 ) / 4
299 klval( 4 ) = ( m+1 ) / 4
300*
301* Do for each value of N in NVAL
302*
303 DO 150 in = 1, nn
304 n = nval( in )
305 xtype = 'N'
306*
307* Set values to use for the upper bandwidth.
308*
309 kuval( 2 ) = n + ( n+1 ) / 4
310*
311* KUVAL( 2 ) = MAX( N-1, 0 )
312*
313 kuval( 3 ) = ( 3*n-1 ) / 4
314 kuval( 4 ) = ( n+1 ) / 4
315*
316* Set limits on the number of loop iterations.
317*
318 nkl = min( m+1, 4 )
319 IF( n.EQ.0 )
320 $ nkl = 2
321 nku = min( n+1, 4 )
322 IF( m.EQ.0 )
323 $ nku = 2
324 nimat = ntypes
325 IF( m.LE.0 .OR. n.LE.0 )
326 $ nimat = 1
327*
328 DO 140 ikl = 1, nkl
329*
330* Do for KL = 0, (5*M+1)/4, (3M-1)/4, and (M+1)/4. This
331* order makes it easier to skip redundant values for small
332* values of M.
333*
334 kl = klval( ikl )
335 DO 130 iku = 1, nku
336*
337* Do for KU = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This
338* order makes it easier to skip redundant values for
339* small values of N.
340*
341 ku = kuval( iku )
342*
343* Check that A and AFAC are big enough to generate this
344* matrix.
345*
346 lda = kl + ku + 1
347 ldafac = 2*kl + ku + 1
348 IF( ( lda*n ).GT.la .OR. ( ldafac*n ).GT.lafac ) THEN
349 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
350 $ CALL alahd( nout, path )
351 IF( n*( kl+ku+1 ).GT.la ) THEN
352 WRITE( nout, fmt = 9999 )la, m, n, kl, ku,
353 $ n*( kl+ku+1 )
354 nerrs = nerrs + 1
355 END IF
356 IF( n*( 2*kl+ku+1 ).GT.lafac ) THEN
357 WRITE( nout, fmt = 9998 )lafac, m, n, kl, ku,
358 $ n*( 2*kl+ku+1 )
359 nerrs = nerrs + 1
360 END IF
361 GO TO 130
362 END IF
363*
364 DO 120 imat = 1, nimat
365*
366* Do the tests only if DOTYPE( IMAT ) is true.
367*
368 IF( .NOT.dotype( imat ) )
369 $ GO TO 120
370*
371* Skip types 2, 3, or 4 if the matrix size is too
372* small.
373*
374 zerot = imat.GE.2 .AND. imat.LE.4
375 IF( zerot .AND. n.LT.imat-1 )
376 $ GO TO 120
377*
378 IF( .NOT.zerot .OR. .NOT.dotype( 1 ) ) THEN
379*
380* Set up parameters with CLATB4 and generate a
381* test matrix with CLATMS.
382*
383 CALL clatb4( path, imat, m, n, TYPE, kl, ku,
384 $ anorm, mode, cndnum, dist )
385*
386 koff = max( 1, ku+2-n )
387 DO 20 i = 1, koff - 1
388 a( i ) = zero
389 20 CONTINUE
390 srnamt = 'CLATMS'
391 CALL clatms( m, n, dist, iseed, TYPE, rwork,
392 $ mode, cndnum, anorm, kl, ku, 'Z',
393 $ a( koff ), lda, work, info )
394*
395* Check the error code from CLATMS.
396*
397 IF( info.NE.0 ) THEN
398 CALL alaerh( path, 'CLATMS', info, 0, ' ', m,
399 $ n, kl, ku, -1, imat, nfail,
400 $ nerrs, nout )
401 GO TO 120
402 END IF
403 ELSE IF( izero.GT.0 ) THEN
404*
405* Use the same matrix for types 3 and 4 as for
406* type 2 by copying back the zeroed out column.
407*
408 CALL ccopy( i2-i1+1, b, 1, a( ioff+i1 ), 1 )
409 END IF
410*
411* For types 2, 3, and 4, zero one or more columns of
412* the matrix to test that INFO is returned correctly.
413*
414 izero = 0
415 IF( zerot ) THEN
416 IF( imat.EQ.2 ) THEN
417 izero = 1
418 ELSE IF( imat.EQ.3 ) THEN
419 izero = min( m, n )
420 ELSE
421 izero = min( m, n ) / 2 + 1
422 END IF
423 ioff = ( izero-1 )*lda
424 IF( imat.LT.4 ) THEN
425*
426* Store the column to be zeroed out in B.
427*
428 i1 = max( 1, ku+2-izero )
429 i2 = min( kl+ku+1, ku+1+( m-izero ) )
430 CALL ccopy( i2-i1+1, a( ioff+i1 ), 1, b, 1 )
431*
432 DO 30 i = i1, i2
433 a( ioff+i ) = zero
434 30 CONTINUE
435 ELSE
436 DO 50 j = izero, n
437 DO 40 i = max( 1, ku+2-j ),
438 $ min( kl+ku+1, ku+1+( m-j ) )
439 a( ioff+i ) = zero
440 40 CONTINUE
441 ioff = ioff + lda
442 50 CONTINUE
443 END IF
444 END IF
445*
446* These lines, if used in place of the calls in the
447* loop over INB, cause the code to bomb on a Sun
448* SPARCstation.
449*
450* ANORMO = CLANGB( 'O', N, KL, KU, A, LDA, RWORK )
451* ANORMI = CLANGB( 'I', N, KL, KU, A, LDA, RWORK )
452*
453* Do for each blocksize in NBVAL
454*
455 DO 110 inb = 1, nnb
456 nb = nbval( inb )
457 CALL xlaenv( 1, nb )
458*
459* Compute the LU factorization of the band matrix.
460*
461 IF( m.GT.0 .AND. n.GT.0 )
462 $ CALL clacpy( 'Full', kl+ku+1, n, a, lda,
463 $ afac( kl+1 ), ldafac )
464 srnamt = 'CGBTRF'
465 CALL cgbtrf( m, n, kl, ku, afac, ldafac, iwork,
466 $ info )
467*
468* Check error code from CGBTRF.
469*
470 IF( info.NE.izero )
471 $ CALL alaerh( path, 'CGBTRF', info, izero,
472 $ ' ', m, n, kl, ku, nb, imat,
473 $ nfail, nerrs, nout )
474 trfcon = .false.
475*
476*+ TEST 1
477* Reconstruct matrix from factors and compute
478* residual.
479*
480 CALL cgbt01( m, n, kl, ku, a, lda, afac, ldafac,
481 $ iwork, work, result( 1 ) )
482*
483* Print information about the tests so far that
484* did not pass the threshold.
485*
486 IF( result( 1 ).GE.thresh ) THEN
487 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
488 $ CALL alahd( nout, path )
489 WRITE( nout, fmt = 9997 )m, n, kl, ku, nb,
490 $ imat, 1, result( 1 )
491 nfail = nfail + 1
492 END IF
493 nrun = nrun + 1
494*
495* Skip the remaining tests if this is not the
496* first block size or if M .ne. N.
497*
498 IF( inb.GT.1 .OR. m.NE.n )
499 $ GO TO 110
500*
501 anormo = clangb( 'O', n, kl, ku, a, lda, rwork )
502 anormi = clangb( 'I', n, kl, ku, a, lda, rwork )
503*
504 IF( info.EQ.0 ) THEN
505*
506* Form the inverse of A so we can get a good
507* estimate of CNDNUM = norm(A) * norm(inv(A)).
508*
509 ldb = max( 1, n )
510 CALL claset( 'Full', n, n, cmplx( zero ),
511 $ cmplx( one ), work, ldb )
512 srnamt = 'CGBTRS'
513 CALL cgbtrs( 'No transpose', n, kl, ku, n,
514 $ afac, ldafac, iwork, work, ldb,
515 $ info )
516*
517* Compute the 1-norm condition number of A.
518*
519 ainvnm = clange( 'O', n, n, work, ldb,
520 $ rwork )
521 IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
522 rcondo = one
523 ELSE
524 rcondo = ( one / anormo ) / ainvnm
525 END IF
526*
527* Compute the infinity-norm condition number of
528* A.
529*
530 ainvnm = clange( 'I', n, n, work, ldb,
531 $ rwork )
532 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
533 rcondi = one
534 ELSE
535 rcondi = ( one / anormi ) / ainvnm
536 END IF
537 ELSE
538*
539* Do only the condition estimate if INFO.NE.0.
540*
541 trfcon = .true.
542 rcondo = zero
543 rcondi = zero
544 END IF
545*
546* Skip the solve tests if the matrix is singular.
547*
548 IF( trfcon )
549 $ GO TO 90
550*
551 DO 80 irhs = 1, nns
552 nrhs = nsval( irhs )
553 xtype = 'N'
554*
555 DO 70 itran = 1, ntran
556 trans = transs( itran )
557 IF( itran.EQ.1 ) THEN
558 rcondc = rcondo
559 norm = 'O'
560 ELSE
561 rcondc = rcondi
562 norm = 'I'
563 END IF
564*
565*+ TEST 2:
566* Solve and compute residual for op(A) * X = B.
567*
568 srnamt = 'CLARHS'
569 CALL clarhs( path, xtype, ' ', trans, n,
570 $ n, kl, ku, nrhs, a, lda,
571 $ xact, ldb, b, ldb, iseed,
572 $ info )
573 xtype = 'C'
574 CALL clacpy( 'Full', n, nrhs, b, ldb, x,
575 $ ldb )
576*
577 srnamt = 'CGBTRS'
578 CALL cgbtrs( trans, n, kl, ku, nrhs, afac,
579 $ ldafac, iwork, x, ldb, info )
580*
581* Check error code from CGBTRS.
582*
583 IF( info.NE.0 )
584 $ CALL alaerh( path, 'CGBTRS', info, 0,
585 $ trans, n, n, kl, ku, -1,
586 $ imat, nfail, nerrs, nout )
587*
588 CALL clacpy( 'Full', n, nrhs, b, ldb,
589 $ work, ldb )
590 CALL cgbt02( trans, m, n, kl, ku, nrhs, a,
591 $ lda, x, ldb, work, ldb,
592 $ rwork, result( 2 ) )
593*
594*+ TEST 3:
595* Check solution from generated exact
596* solution.
597*
598 CALL cget04( n, nrhs, x, ldb, xact, ldb,
599 $ rcondc, result( 3 ) )
600*
601*+ TESTS 4, 5, 6:
602* Use iterative refinement to improve the
603* solution.
604*
605 srnamt = 'CGBRFS'
606 CALL cgbrfs( trans, n, kl, ku, nrhs, a,
607 $ lda, afac, ldafac, iwork, b,
608 $ ldb, x, ldb, rwork,
609 $ rwork( nrhs+1 ), work,
610 $ rwork( 2*nrhs+1 ), info )
611*
612* Check error code from CGBRFS.
613*
614 IF( info.NE.0 )
615 $ CALL alaerh( path, 'CGBRFS', info, 0,
616 $ trans, n, n, kl, ku, nrhs,
617 $ imat, nfail, nerrs, nout )
618*
619 CALL cget04( n, nrhs, x, ldb, xact, ldb,
620 $ rcondc, result( 4 ) )
621 CALL cgbt05( trans, n, kl, ku, nrhs, a,
622 $ lda, b, ldb, x, ldb, xact,
623 $ ldb, rwork, rwork( nrhs+1 ),
624 $ result( 5 ) )
625*
626* Print information about the tests that did
627* not pass the threshold.
628*
629 DO 60 k = 2, 6
630 IF( result( k ).GE.thresh ) THEN
631 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
632 $ CALL alahd( nout, path )
633 WRITE( nout, fmt = 9996 )trans, n,
634 $ kl, ku, nrhs, imat, k,
635 $ result( k )
636 nfail = nfail + 1
637 END IF
638 60 CONTINUE
639 nrun = nrun + 5
640 70 CONTINUE
641 80 CONTINUE
642*
643*+ TEST 7:
644* Get an estimate of RCOND = 1/CNDNUM.
645*
646 90 CONTINUE
647 DO 100 itran = 1, 2
648 IF( itran.EQ.1 ) THEN
649 anorm = anormo
650 rcondc = rcondo
651 norm = 'O'
652 ELSE
653 anorm = anormi
654 rcondc = rcondi
655 norm = 'I'
656 END IF
657 srnamt = 'CGBCON'
658 CALL cgbcon( norm, n, kl, ku, afac, ldafac,
659 $ iwork, anorm, rcond, work,
660 $ rwork, info )
661*
662* Check error code from CGBCON.
663*
664 IF( info.NE.0 )
665 $ CALL alaerh( path, 'CGBCON', info, 0,
666 $ norm, n, n, kl, ku, -1, imat,
667 $ nfail, nerrs, nout )
668*
669 result( 7 ) = sget06( rcond, rcondc )
670*
671* Print information about the tests that did
672* not pass the threshold.
673*
674 IF( result( 7 ).GE.thresh ) THEN
675 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
676 $ CALL alahd( nout, path )
677 WRITE( nout, fmt = 9995 )norm, n, kl, ku,
678 $ imat, 7, result( 7 )
679 nfail = nfail + 1
680 END IF
681 nrun = nrun + 1
682 100 CONTINUE
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 CCHKGB, 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 CCHKGB, 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 CCHKGB
710*
711 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine clarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
CLARHS
Definition clarhs.f:208
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine cchkgb(dotype, nm, mval, nn, nval, nnb, nbval, nns, nsval, thresh, tsterr, a, la, afac, lafac, b, x, xact, work, rwork, iwork, nout)
CCHKGB
Definition cchkgb.f:191
subroutine cerrge(path, nunit)
CERRGE
Definition cerrge.f:55
subroutine cgbt01(m, n, kl, ku, a, lda, afac, ldafac, ipiv, work, resid)
CGBT01
Definition cgbt01.f:126
subroutine cgbt02(trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
CGBT02
Definition cgbt02.f:148
subroutine cgbt05(trans, n, kl, ku, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
CGBT05
Definition cgbt05.f:176
subroutine cget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
CGET04
Definition cget04.f:102
subroutine clatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
CLATB4
Definition clatb4.f:121
subroutine clatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
CLATMS
Definition clatms.f:332
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cgbcon(norm, n, kl, ku, ab, ldab, ipiv, anorm, rcond, work, rwork, info)
CGBCON
Definition cgbcon.f:147
subroutine cgbrfs(trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CGBRFS
Definition cgbrfs.f:206
subroutine cgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
CGBTRF
Definition cgbtrf.f:144
subroutine cgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
CGBTRS
Definition cgbtrs.f:138
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