LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchktb.f
Go to the documentation of this file.
1*> \brief \b ZCHKTB
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 ZCHKTB( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
12* NMAX, AB, AINV, B, X, XACT, WORK, RWORK, NOUT )
13*
14* .. Scalar Arguments ..
15* LOGICAL TSTERR
16* INTEGER NMAX, NN, NNS, NOUT
17* DOUBLE PRECISION THRESH
18* ..
19* .. Array Arguments ..
20* LOGICAL DOTYPE( * )
21* INTEGER NSVAL( * ), NVAL( * )
22* DOUBLE PRECISION RWORK( * )
23* COMPLEX*16 AB( * ), AINV( * ), B( * ), WORK( * ), X( * ),
24* $ XACT( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> ZCHKTB tests ZTBTRS, -RFS, and -CON, and ZLATBS.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] DOTYPE
40*> \verbatim
41*> DOTYPE is LOGICAL array, dimension (NTYPES)
42*> The matrix types to be used for testing. Matrices of type j
43*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
44*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
45*> \endverbatim
46*>
47*> \param[in] NN
48*> \verbatim
49*> NN is INTEGER
50*> The number of values of N contained in the vector NVAL.
51*> \endverbatim
52*>
53*> \param[in] NVAL
54*> \verbatim
55*> NVAL is INTEGER array, dimension (NN)
56*> The values of the matrix column dimension N.
57*> \endverbatim
58*>
59*> \param[in] NNS
60*> \verbatim
61*> NNS is INTEGER
62*> The number of values of NRHS contained in the vector NSVAL.
63*> \endverbatim
64*>
65*> \param[in] NSVAL
66*> \verbatim
67*> NSVAL is INTEGER array, dimension (NNS)
68*> The values of the number of right hand sides NRHS.
69*> \endverbatim
70*>
71*> \param[in] THRESH
72*> \verbatim
73*> THRESH is DOUBLE PRECISION
74*> The threshold value for the test ratios. A result is
75*> included in the output file if RESULT >= THRESH. To have
76*> every test ratio printed, use THRESH = 0.
77*> \endverbatim
78*>
79*> \param[in] TSTERR
80*> \verbatim
81*> TSTERR is LOGICAL
82*> Flag that indicates whether error exits are to be tested.
83*> \endverbatim
84*>
85*> \param[in] NMAX
86*> \verbatim
87*> NMAX is INTEGER
88*> The leading dimension of the work arrays.
89*> NMAX >= the maximum value of N in NVAL.
90*> \endverbatim
91*>
92*> \param[out] AB
93*> \verbatim
94*> AB is COMPLEX*16 array, dimension (NMAX*NMAX)
95*> \endverbatim
96*>
97*> \param[out] AINV
98*> \verbatim
99*> AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
100*> \endverbatim
101*>
102*> \param[out] B
103*> \verbatim
104*> B is COMPLEX*16 array, dimension (NMAX*NSMAX)
105*> where NSMAX is the largest entry in NSVAL.
106*> \endverbatim
107*>
108*> \param[out] X
109*> \verbatim
110*> X is COMPLEX*16 array, dimension (NMAX*NSMAX)
111*> \endverbatim
112*>
113*> \param[out] XACT
114*> \verbatim
115*> XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
116*> \endverbatim
117*>
118*> \param[out] WORK
119*> \verbatim
120*> WORK is COMPLEX*16 array, dimension
121*> (NMAX*max(3,NSMAX))
122*> \endverbatim
123*>
124*> \param[out] RWORK
125*> \verbatim
126*> RWORK is DOUBLE PRECISION array, dimension
127*> (max(NMAX,2*NSMAX))
128*> \endverbatim
129*>
130*> \param[in] NOUT
131*> \verbatim
132*> NOUT is INTEGER
133*> The unit number for output.
134*> \endverbatim
135*
136* Authors:
137* ========
138*
139*> \author Univ. of Tennessee
140*> \author Univ. of California Berkeley
141*> \author Univ. of Colorado Denver
142*> \author NAG Ltd.
143*
144*> \ingroup complex16_lin
145*
146* =====================================================================
147 SUBROUTINE zchktb( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
148 $ NMAX, AB, AINV, B, X, XACT, WORK, RWORK, NOUT )
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 DOUBLE PRECISION THRESH
158* ..
159* .. Array Arguments ..
160 LOGICAL DOTYPE( * )
161 INTEGER NSVAL( * ), NVAL( * )
162 DOUBLE PRECISION RWORK( * )
163 COMPLEX*16 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 DOUBLE PRECISION ONE, ZERO
177 parameter( one = 1.0d+0, zero = 0.0d+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 DOUBLE PRECISION 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 DOUBLE PRECISION RESULT( NTESTS )
192* ..
193* .. External Functions ..
194 LOGICAL LSAME
195 DOUBLE PRECISION ZLANTB, ZLANTR
196 EXTERNAL lsame, zlantb, zlantr
197* ..
198* .. External Subroutines ..
199 EXTERNAL alaerh, alahd, alasum, zcopy, zerrtr, zget04,
202 $ ztbtrs
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 dcmplx, 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 ) = 'Zomplex 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 zerrtr( 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 ZLATTB to generate a triangular test matrix.
284*
285 srnamt = 'ZLATTB'
286 CALL zlattb( 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 zlaset( 'Full', n, n, dcmplx( zero ),
301 $ dcmplx( one ), ainv, lda )
302 IF( lsame( uplo, 'U' ) ) THEN
303 DO 20 j = 1, n
304 CALL ztbsv( 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 ztbsv( 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 = zlantb( '1', uplo, diag, n, kd, ab, ldab,
318 $ rwork )
319 ainvnm = zlantr( '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 = zlantb( 'I', uplo, diag, n, kd, ab, ldab,
330 $ rwork )
331 ainvnm = zlantr( '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 = 'ZLARHS'
360 CALL zlarhs( path, xtype, uplo, trans, n, n, kd,
361 $ idiag, nrhs, ab, ldab, xact, lda,
362 $ b, lda, iseed, info )
363 xtype = 'C'
364 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
365*
366 srnamt = 'ZTBTRS'
367 CALL ztbtrs( uplo, trans, diag, n, kd, nrhs, ab,
368 $ ldab, x, lda, info )
369*
370* Check error code from ZTBTRS.
371*
372 IF( info.NE.0 )
373 $ CALL alaerh( path, 'ZTBTRS', info, 0,
374 $ uplo // trans // diag, n, n, kd,
375 $ kd, nrhs, imat, nfail, nerrs,
376 $ nout )
377*
378 CALL ztbt02( 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 zget04( 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 = 'ZTBRFS'
393 CALL ztbrfs( 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 ZTBRFS.
399*
400 IF( info.NE.0 )
401 $ CALL alaerh( path, 'ZTBRFS', info, 0,
402 $ uplo // trans // diag, n, n, kd,
403 $ kd, nrhs, imat, nfail, nerrs,
404 $ nout )
405*
406 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
407 $ result( 3 ) )
408 CALL ztbt05( 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 = 'ZTBCON'
441 CALL ztbcon( norm, uplo, diag, n, kd, ab, ldab,
442 $ rcond, work, rwork, info )
443*
444* Check error code from ZTBCON.
445*
446 IF( info.NE.0 )
447 $ CALL alaerh( path, 'ZTBCON', info, 0,
448 $ norm // uplo // diag, n, n, kd, kd,
449 $ -1, imat, nfail, nerrs, nout )
450*
451 CALL ztbt06( 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 ) 'ZTBCON', 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 ZLATBS.
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 ZLATTB to generate a triangular test matrix.
489*
490 srnamt = 'ZLATTB'
491 CALL zlattb( 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 = 'ZLATBS'
498 CALL zcopy( n, x, 1, b, 1 )
499 CALL zlatbs( uplo, trans, diag, 'N', n, kd, ab,
500 $ ldab, b, scale, rwork, info )
501*
502* Check error code from ZLATBS.
503*
504 IF( info.NE.0 )
505 $ CALL alaerh( path, 'ZLATBS', info, 0,
506 $ uplo // trans // diag // 'N', n, n,
507 $ kd, kd, -1, imat, nfail, nerrs,
508 $ nout )
509*
510 CALL ztbt03( 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 zcopy( n, x, 1, b, 1 )
518 CALL zlatbs( uplo, trans, diag, 'Y', n, kd, ab,
519 $ ldab, b, scale, rwork, info )
520*
521* Check error code from ZLATBS.
522*
523 IF( info.NE.0 )
524 $ CALL alaerh( path, 'ZLATBS', info, 0,
525 $ uplo // trans // diag // 'Y', n, n,
526 $ kd, kd, -1, imat, nfail, nerrs,
527 $ nout )
528*
529 CALL ztbt03( 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 )'ZLATBS', 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 )'ZLATBS', 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 ZCHKTB
573*
574 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.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 zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zlatbs(uplo, trans, diag, normin, n, kd, ab, ldab, x, scale, cnorm, info)
ZLATBS solves a triangular banded system of equations.
Definition zlatbs.f:243
subroutine ztbcon(norm, uplo, diag, n, kd, ab, ldab, rcond, work, rwork, info)
ZTBCON
Definition ztbcon.f:143
subroutine ztbrfs(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZTBRFS
Definition ztbrfs.f:188
subroutine ztbsv(uplo, trans, diag, n, k, a, lda, x, incx)
ZTBSV
Definition ztbsv.f:189
subroutine ztbtrs(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, info)
ZTBTRS
Definition ztbtrs.f:146
subroutine zchktb(dotype, nn, nval, nns, nsval, thresh, tsterr, nmax, ab, ainv, b, x, xact, work, rwork, nout)
ZCHKTB
Definition zchktb.f:149
subroutine zerrtr(path, nunit)
ZERRTR
Definition zerrtr.f:54
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zlattb(imat, uplo, trans, diag, iseed, n, kd, ab, ldab, b, work, rwork, info)
ZLATTB
Definition zlattb.f:141
subroutine ztbt02(uplo, trans, diag, n, kd, nrhs, ab, ldab, x, ldx, b, ldb, work, rwork, resid)
ZTBT02
Definition ztbt02.f:159
subroutine ztbt03(uplo, trans, diag, n, kd, nrhs, ab, ldab, scale, cnorm, tscal, x, ldx, b, ldb, work, resid)
ZTBT03
Definition ztbt03.f:177
subroutine ztbt05(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZTBT05
Definition ztbt05.f:189
subroutine ztbt06(rcond, rcondc, uplo, diag, n, kd, ab, ldab, rwork, rat)
ZTBT06
Definition ztbt06.f:126