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

◆ dchktb()

subroutine dchktb ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nns,
integer, dimension( * )  nsval,
double precision  thresh,
logical  tsterr,
integer  nmax,
double precision, dimension( * )  ab,
double precision, dimension( * )  ainv,
double precision, dimension( * )  b,
double precision, dimension( * )  x,
double precision, dimension( * )  xact,
double precision, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

DCHKTB

Purpose:
 DCHKTB tests DTBTRS, -RFS, and -CON, and DLATBS.
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 DOUBLE PRECISION
          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 DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AINV
          AINV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(NMAX,2*NSMAX))
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[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 152 of file dchktb.f.

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