LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ cchktb()

subroutine cchktb ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nns,
integer, dimension( * )  nsval,
real  thresh,
logical  tsterr,
integer  nmax,
complex, dimension( * )  ab,
complex, dimension( * )  ainv,
complex, dimension( * )  b,
complex, dimension( * )  x,
complex, dimension( * )  xact,
complex, dimension( * )  work,
real, dimension( * )  rwork,
integer  nout 
)

CCHKTB

Purpose:
 CCHKTB tests CTBTRS, -RFS, and -CON, and CLATBS.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix column dimension N.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[in]THRESH
          THRESH is REAL
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The leading dimension of the work arrays.
          NMAX >= the maximum value of N in NVAL.
[out]AB
          AB is COMPLEX array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL array, dimension
                      (max(NMAX,2*NSMAX))
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 147 of file cchktb.f.

149*
150* -- LAPACK test routine --
151* -- LAPACK is a software package provided by Univ. of Tennessee, --
152* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
153*
154* .. Scalar Arguments ..
155 LOGICAL TSTERR
156 INTEGER NMAX, NN, NNS, NOUT
157 REAL THRESH
158* ..
159* .. Array Arguments ..
160 LOGICAL DOTYPE( * )
161 INTEGER NSVAL( * ), NVAL( * )
162 REAL RWORK( * )
163 COMPLEX AB( * ), AINV( * ), B( * ), WORK( * ), X( * ),
164 $ XACT( * )
165* ..
166*
167* =====================================================================
168*
169* .. Parameters ..
170 INTEGER NTYPE1, NTYPES
171 parameter( ntype1 = 9, ntypes = 17 )
172 INTEGER NTESTS
173 parameter( ntests = 8 )
174 INTEGER NTRAN
175 parameter( ntran = 3 )
176 REAL ONE, ZERO
177 parameter( one = 1.0e+0, zero = 0.0e+0 )
178* ..
179* .. Local Scalars ..
180 CHARACTER DIAG, NORM, TRANS, UPLO, XTYPE
181 CHARACTER*3 PATH
182 INTEGER I, IDIAG, IK, IMAT, IN, INFO, IRHS, ITRAN,
183 $ IUPLO, J, K, KD, LDA, LDAB, N, NERRS, NFAIL,
184 $ NIMAT, NIMAT2, NK, NRHS, NRUN
185 REAL AINVNM, ANORM, RCOND, RCONDC, RCONDI, RCONDO,
186 $ SCALE
187* ..
188* .. Local Arrays ..
189 CHARACTER TRANSS( NTRAN ), UPLOS( 2 )
190 INTEGER ISEED( 4 ), ISEEDY( 4 )
191 REAL RESULT( NTESTS )
192* ..
193* .. External Functions ..
194 LOGICAL LSAME
195 REAL CLANTB, CLANTR
196 EXTERNAL lsame, clantb, clantr
197* ..
198* .. External Subroutines ..
199 EXTERNAL alaerh, alahd, alasum, ccopy, cerrtr, cget04,
202 $ ctbtrs
203* ..
204* .. Scalars in Common ..
205 LOGICAL LERR, OK
206 CHARACTER*32 SRNAMT
207 INTEGER INFOT, IOUNIT
208* ..
209* .. Common blocks ..
210 COMMON / infoc / infot, iounit, ok, lerr
211 COMMON / srnamc / srnamt
212* ..
213* .. Intrinsic Functions ..
214 INTRINSIC cmplx, max, min
215* ..
216* .. Data statements ..
217 DATA iseedy / 1988, 1989, 1990, 1991 /
218 DATA uplos / 'U', 'L' / , transs / 'N', 'T', 'C' /
219* ..
220* .. Executable Statements ..
221*
222* Initialize constants and the random number seed.
223*
224 path( 1: 1 ) = 'Complex precision'
225 path( 2: 3 ) = 'TB'
226 nrun = 0
227 nfail = 0
228 nerrs = 0
229 DO 10 i = 1, 4
230 iseed( i ) = iseedy( i )
231 10 CONTINUE
232*
233* Test the error exits
234*
235 IF( tsterr )
236 $ CALL cerrtr( path, nout )
237 infot = 0
238*
239 DO 140 in = 1, nn
240*
241* Do for each value of N in NVAL
242*
243 n = nval( in )
244 lda = max( 1, n )
245 xtype = 'N'
246 nimat = ntype1
247 nimat2 = ntypes
248 IF( n.LE.0 ) THEN
249 nimat = 1
250 nimat2 = ntype1 + 1
251 END IF
252*
253 nk = min( n+1, 4 )
254 DO 130 ik = 1, nk
255*
256* Do for KD = 0, N, (3N-1)/4, and (N+1)/4. This order makes
257* it easier to skip redundant values for small values of N.
258*
259 IF( ik.EQ.1 ) THEN
260 kd = 0
261 ELSE IF( ik.EQ.2 ) THEN
262 kd = max( n, 0 )
263 ELSE IF( ik.EQ.3 ) THEN
264 kd = ( 3*n-1 ) / 4
265 ELSE IF( ik.EQ.4 ) THEN
266 kd = ( n+1 ) / 4
267 END IF
268 ldab = kd + 1
269*
270 DO 90 imat = 1, nimat
271*
272* Do the tests only if DOTYPE( IMAT ) is true.
273*
274 IF( .NOT.dotype( imat ) )
275 $ GO TO 90
276*
277 DO 80 iuplo = 1, 2
278*
279* Do first for UPLO = 'U', then for UPLO = 'L'
280*
281 uplo = uplos( iuplo )
282*
283* Call CLATTB to generate a triangular test matrix.
284*
285 srnamt = 'CLATTB'
286 CALL clattb( imat, uplo, 'No transpose', diag, iseed,
287 $ n, kd, ab, ldab, x, work, rwork, info )
288*
289* Set IDIAG = 1 for non-unit matrices, 2 for unit.
290*
291 IF( lsame( diag, 'N' ) ) THEN
292 idiag = 1
293 ELSE
294 idiag = 2
295 END IF
296*
297* Form the inverse of A so we can get a good estimate
298* of RCONDC = 1/(norm(A) * norm(inv(A))).
299*
300 CALL claset( 'Full', n, n, cmplx( zero ),
301 $ cmplx( one ), ainv, lda )
302 IF( lsame( uplo, 'U' ) ) THEN
303 DO 20 j = 1, n
304 CALL ctbsv( uplo, 'No transpose', diag, j, kd,
305 $ ab, ldab, ainv( ( j-1 )*lda+1 ), 1 )
306 20 CONTINUE
307 ELSE
308 DO 30 j = 1, n
309 CALL ctbsv( uplo, 'No transpose', diag, n-j+1,
310 $ kd, ab( ( j-1 )*ldab+1 ), ldab,
311 $ ainv( ( j-1 )*lda+j ), 1 )
312 30 CONTINUE
313 END IF
314*
315* Compute the 1-norm condition number of A.
316*
317 anorm = clantb( '1', uplo, diag, n, kd, ab, ldab,
318 $ rwork )
319 ainvnm = clantr( '1', uplo, diag, n, n, ainv, lda,
320 $ rwork )
321 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
322 rcondo = one
323 ELSE
324 rcondo = ( one / anorm ) / ainvnm
325 END IF
326*
327* Compute the infinity-norm condition number of A.
328*
329 anorm = clantb( 'I', uplo, diag, n, kd, ab, ldab,
330 $ rwork )
331 ainvnm = clantr( 'I', uplo, diag, n, n, ainv, lda,
332 $ rwork )
333 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
334 rcondi = one
335 ELSE
336 rcondi = ( one / anorm ) / ainvnm
337 END IF
338*
339 DO 60 irhs = 1, nns
340 nrhs = nsval( irhs )
341 xtype = 'N'
342*
343 DO 50 itran = 1, ntran
344*
345* Do for op(A) = A, A**T, or A**H.
346*
347 trans = transs( itran )
348 IF( itran.EQ.1 ) THEN
349 norm = 'O'
350 rcondc = rcondo
351 ELSE
352 norm = 'I'
353 rcondc = rcondi
354 END IF
355*
356*+ TEST 1
357* Solve and compute residual for op(A)*x = b.
358*
359 srnamt = 'CLARHS'
360 CALL clarhs( path, xtype, uplo, trans, n, n, kd,
361 $ idiag, nrhs, ab, ldab, xact, lda,
362 $ b, lda, iseed, info )
363 xtype = 'C'
364 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
365*
366 srnamt = 'CTBTRS'
367 CALL ctbtrs( uplo, trans, diag, n, kd, nrhs, ab,
368 $ ldab, x, lda, info )
369*
370* Check error code from CTBTRS.
371*
372 IF( info.NE.0 )
373 $ CALL alaerh( path, 'CTBTRS', info, 0,
374 $ uplo // trans // diag, n, n, kd,
375 $ kd, nrhs, imat, nfail, nerrs,
376 $ nout )
377*
378 CALL ctbt02( uplo, trans, diag, n, kd, nrhs, ab,
379 $ ldab, x, lda, b, lda, work, rwork,
380 $ result( 1 ) )
381*
382*+ TEST 2
383* Check solution from generated exact solution.
384*
385 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
386 $ result( 2 ) )
387*
388*+ TESTS 3, 4, and 5
389* Use iterative refinement to improve the solution
390* and compute error bounds.
391*
392 srnamt = 'CTBRFS'
393 CALL ctbrfs( uplo, trans, diag, n, kd, nrhs, ab,
394 $ ldab, b, lda, x, lda, rwork,
395 $ rwork( nrhs+1 ), work,
396 $ rwork( 2*nrhs+1 ), info )
397*
398* Check error code from CTBRFS.
399*
400 IF( info.NE.0 )
401 $ CALL alaerh( path, 'CTBRFS', info, 0,
402 $ uplo // trans // diag, n, n, kd,
403 $ kd, nrhs, imat, nfail, nerrs,
404 $ nout )
405*
406 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
407 $ result( 3 ) )
408 CALL ctbt05( uplo, trans, diag, n, kd, nrhs, ab,
409 $ ldab, b, lda, x, lda, xact, lda,
410 $ rwork, rwork( nrhs+1 ),
411 $ result( 4 ) )
412*
413* Print information about the tests that did not
414* pass the threshold.
415*
416 DO 40 k = 1, 5
417 IF( result( k ).GE.thresh ) THEN
418 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
419 $ CALL alahd( nout, path )
420 WRITE( nout, fmt = 9999 )uplo, trans,
421 $ diag, n, kd, nrhs, imat, k, result( k )
422 nfail = nfail + 1
423 END IF
424 40 CONTINUE
425 nrun = nrun + 5
426 50 CONTINUE
427 60 CONTINUE
428*
429*+ TEST 6
430* Get an estimate of RCOND = 1/CNDNUM.
431*
432 DO 70 itran = 1, 2
433 IF( itran.EQ.1 ) THEN
434 norm = 'O'
435 rcondc = rcondo
436 ELSE
437 norm = 'I'
438 rcondc = rcondi
439 END IF
440 srnamt = 'CTBCON'
441 CALL ctbcon( norm, uplo, diag, n, kd, ab, ldab,
442 $ rcond, work, rwork, info )
443*
444* Check error code from CTBCON.
445*
446 IF( info.NE.0 )
447 $ CALL alaerh( path, 'CTBCON', info, 0,
448 $ norm // uplo // diag, n, n, kd, kd,
449 $ -1, imat, nfail, nerrs, nout )
450*
451 CALL ctbt06( rcond, rcondc, uplo, diag, n, kd, ab,
452 $ ldab, rwork, result( 6 ) )
453*
454* Print the test ratio if it is .GE. THRESH.
455*
456 IF( result( 6 ).GE.thresh ) THEN
457 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
458 $ CALL alahd( nout, path )
459 WRITE( nout, fmt = 9998 ) 'CTBCON', norm, uplo,
460 $ diag, n, kd, imat, 6, result( 6 )
461 nfail = nfail + 1
462 END IF
463 nrun = nrun + 1
464 70 CONTINUE
465 80 CONTINUE
466 90 CONTINUE
467*
468* Use pathological test matrices to test CLATBS.
469*
470 DO 120 imat = ntype1 + 1, nimat2
471*
472* Do the tests only if DOTYPE( IMAT ) is true.
473*
474 IF( .NOT.dotype( imat ) )
475 $ GO TO 120
476*
477 DO 110 iuplo = 1, 2
478*
479* Do first for UPLO = 'U', then for UPLO = 'L'
480*
481 uplo = uplos( iuplo )
482 DO 100 itran = 1, ntran
483*
484* Do for op(A) = A, A**T, and A**H.
485*
486 trans = transs( itran )
487*
488* Call CLATTB to generate a triangular test matrix.
489*
490 srnamt = 'CLATTB'
491 CALL clattb( imat, uplo, trans, diag, iseed, n, kd,
492 $ ab, ldab, x, work, rwork, info )
493*
494*+ TEST 7
495* Solve the system op(A)*x = b
496*
497 srnamt = 'CLATBS'
498 CALL ccopy( n, x, 1, b, 1 )
499 CALL clatbs( uplo, trans, diag, 'N', n, kd, ab,
500 $ ldab, b, scale, rwork, info )
501*
502* Check error code from CLATBS.
503*
504 IF( info.NE.0 )
505 $ CALL alaerh( path, 'CLATBS', info, 0,
506 $ uplo // trans // diag // 'N', n, n,
507 $ kd, kd, -1, imat, nfail, nerrs,
508 $ nout )
509*
510 CALL ctbt03( uplo, trans, diag, n, kd, 1, ab, ldab,
511 $ scale, rwork, one, b, lda, x, lda,
512 $ work, result( 7 ) )
513*
514*+ TEST 8
515* Solve op(A)*x = b again with NORMIN = 'Y'.
516*
517 CALL ccopy( n, x, 1, b, 1 )
518 CALL clatbs( uplo, trans, diag, 'Y', n, kd, ab,
519 $ ldab, b, scale, rwork, info )
520*
521* Check error code from CLATBS.
522*
523 IF( info.NE.0 )
524 $ CALL alaerh( path, 'CLATBS', info, 0,
525 $ uplo // trans // diag // 'Y', n, n,
526 $ kd, kd, -1, imat, nfail, nerrs,
527 $ nout )
528*
529 CALL ctbt03( uplo, trans, diag, n, kd, 1, ab, ldab,
530 $ scale, rwork, one, b, lda, x, lda,
531 $ work, result( 8 ) )
532*
533* Print information about the tests that did not pass
534* the threshold.
535*
536 IF( result( 7 ).GE.thresh ) THEN
537 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
538 $ CALL alahd( nout, path )
539 WRITE( nout, fmt = 9997 )'CLATBS', uplo, trans,
540 $ diag, 'N', n, kd, imat, 7, result( 7 )
541 nfail = nfail + 1
542 END IF
543 IF( result( 8 ).GE.thresh ) THEN
544 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
545 $ CALL alahd( nout, path )
546 WRITE( nout, fmt = 9997 )'CLATBS', uplo, trans,
547 $ diag, 'Y', n, kd, imat, 8, result( 8 )
548 nfail = nfail + 1
549 END IF
550 nrun = nrun + 2
551 100 CONTINUE
552 110 CONTINUE
553 120 CONTINUE
554 130 CONTINUE
555 140 CONTINUE
556*
557* Print a summary of the results.
558*
559 CALL alasum( path, nout, nfail, nrun, nerrs )
560*
561 9999 FORMAT( ' UPLO=''', a1, ''', TRANS=''', a1, ''',
562 $ DIAG=''', a1, ''', N=', i5, ', KD=', i5, ', NRHS=', i5,
563 $ ', type ', i2, ', test(', i2, ')=', g12.5 )
564 9998 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ''', a1, ''',',
565 $ i5, ',', i5, ', ... ), type ', i2, ', test(', i2, ')=',
566 $ g12.5 )
567 9997 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ''', a1, ''', ''',
568 $ a1, ''',', i5, ',', i5, ', ... ), type ', i2, ', test(',
569 $ i1, ')=', g12.5 )
570 RETURN
571*
572* End of CCHKTB
573*
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 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 cerrtr(path, nunit)
CERRTR
Definition cerrtr.f:54
subroutine cget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
CGET04
Definition cget04.f:102
subroutine clattb(imat, uplo, trans, diag, iseed, n, kd, ab, ldab, b, work, rwork, info)
CLATTB
Definition clattb.f:141
subroutine ctbt02(uplo, trans, diag, n, kd, nrhs, ab, ldab, x, ldx, b, ldb, work, rwork, resid)
CTBT02
Definition ctbt02.f:159
subroutine ctbt03(uplo, trans, diag, n, kd, nrhs, ab, ldab, scale, cnorm, tscal, x, ldx, b, ldb, work, resid)
CTBT03
Definition ctbt03.f:177
subroutine ctbt05(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
CTBT05
Definition ctbt05.f:189
subroutine ctbt06(rcond, rcondc, uplo, diag, n, kd, ab, ldab, rwork, rat)
CTBT06
Definition ctbt06.f:126
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
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
real function clantb(norm, uplo, diag, n, k, ab, ldab, work)
CLANTB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clantb.f:141
real function clantr(norm, uplo, diag, m, n, a, lda, work)
CLANTR returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clantr.f:142
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 clatbs(uplo, trans, diag, normin, n, kd, ab, ldab, x, scale, cnorm, info)
CLATBS solves a triangular banded system of equations.
Definition clatbs.f:243
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine ctbcon(norm, uplo, diag, n, kd, ab, ldab, rcond, work, rwork, info)
CTBCON
Definition ctbcon.f:143
subroutine ctbrfs(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CTBRFS
Definition ctbrfs.f:188
subroutine ctbsv(uplo, trans, diag, n, k, a, lda, x, incx)
CTBSV
Definition ctbsv.f:189
subroutine ctbtrs(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, info)
CTBTRS
Definition ctbtrs.f:146
Here is the call graph for this function:
Here is the caller graph for this function: