LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cchkpb.f
Go to the documentation of this file.
1*> \brief \b CCHKPB
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 CCHKPB( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
12* THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
13* XACT, WORK, RWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NNB, NNS, NOUT
18* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER NBVAL( * ), NSVAL( * ), NVAL( * )
23* REAL RWORK( * )
24* COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> CCHKPB tests CPBTRF, -TRS, -RFS, and -CON.
35*> \endverbatim
36*
37* Arguments:
38* ==========
39*
40*> \param[in] DOTYPE
41*> \verbatim
42*> DOTYPE is LOGICAL array, dimension (NTYPES)
43*> The matrix types to be used for testing. Matrices of type j
44*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
45*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
46*> \endverbatim
47*>
48*> \param[in] NN
49*> \verbatim
50*> NN is INTEGER
51*> The number of values of N contained in the vector NVAL.
52*> \endverbatim
53*>
54*> \param[in] NVAL
55*> \verbatim
56*> NVAL is INTEGER array, dimension (NN)
57*> The values of the matrix dimension N.
58*> \endverbatim
59*>
60*> \param[in] NNB
61*> \verbatim
62*> NNB is INTEGER
63*> The number of values of NB contained in the vector NBVAL.
64*> \endverbatim
65*>
66*> \param[in] NBVAL
67*> \verbatim
68*> NBVAL is INTEGER array, dimension (NNB)
69*> The values of the blocksize NB.
70*> \endverbatim
71*>
72*> \param[in] NNS
73*> \verbatim
74*> NNS is INTEGER
75*> The number of values of NRHS contained in the vector NSVAL.
76*> \endverbatim
77*>
78*> \param[in] NSVAL
79*> \verbatim
80*> NSVAL is INTEGER array, dimension (NNS)
81*> The values of the number of right hand sides NRHS.
82*> \endverbatim
83*>
84*> \param[in] THRESH
85*> \verbatim
86*> THRESH is REAL
87*> The threshold value for the test ratios. A result is
88*> included in the output file if RESULT >= THRESH. To have
89*> every test ratio printed, use THRESH = 0.
90*> \endverbatim
91*>
92*> \param[in] TSTERR
93*> \verbatim
94*> TSTERR is LOGICAL
95*> Flag that indicates whether error exits are to be tested.
96*> \endverbatim
97*>
98*> \param[in] NMAX
99*> \verbatim
100*> NMAX is INTEGER
101*> The maximum value permitted for N, used in dimensioning the
102*> work arrays.
103*> \endverbatim
104*>
105*> \param[out] A
106*> \verbatim
107*> A is REAL array, dimension (NMAX*NMAX)
108*> \endverbatim
109*>
110*> \param[out] AFAC
111*> \verbatim
112*> AFAC is REAL array, dimension (NMAX*NMAX)
113*> \endverbatim
114*>
115*> \param[out] AINV
116*> \verbatim
117*> AINV is REAL array, dimension (NMAX*NMAX)
118*> \endverbatim
119*>
120*> \param[out] B
121*> \verbatim
122*> B is REAL array, dimension (NMAX*NSMAX)
123*> where NSMAX is the largest entry in NSVAL.
124*> \endverbatim
125*>
126*> \param[out] X
127*> \verbatim
128*> X is REAL array, dimension (NMAX*NSMAX)
129*> \endverbatim
130*>
131*> \param[out] XACT
132*> \verbatim
133*> XACT is REAL array, dimension (NMAX*NSMAX)
134*> \endverbatim
135*>
136*> \param[out] WORK
137*> \verbatim
138*> WORK is REAL array, dimension
139*> (NMAX*max(3,NSMAX))
140*> \endverbatim
141*>
142*> \param[out] RWORK
143*> \verbatim
144*> RWORK is REAL array, dimension
145*> (max(NMAX,2*NSMAX))
146*> \endverbatim
147*>
148*> \param[in] NOUT
149*> \verbatim
150*> NOUT is INTEGER
151*> The unit number for output.
152*> \endverbatim
153*
154* Authors:
155* ========
156*
157*> \author Univ. of Tennessee
158*> \author Univ. of California Berkeley
159*> \author Univ. of Colorado Denver
160*> \author NAG Ltd.
161*
162*> \ingroup complex_lin
163*
164* =====================================================================
165 SUBROUTINE cchkpb( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
166 $ THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
167 $ XACT, WORK, RWORK, NOUT )
168*
169* -- LAPACK test routine --
170* -- LAPACK is a software package provided by Univ. of Tennessee, --
171* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
172*
173* .. Scalar Arguments ..
174 LOGICAL TSTERR
175 INTEGER NMAX, NN, NNB, NNS, NOUT
176 REAL THRESH
177* ..
178* .. Array Arguments ..
179 LOGICAL DOTYPE( * )
180 INTEGER NBVAL( * ), NSVAL( * ), NVAL( * )
181 REAL RWORK( * )
182 COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ),
183 $ work( * ), x( * ), xact( * )
184* ..
185*
186* =====================================================================
187*
188* .. Parameters ..
189 REAL ONE, ZERO
190 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
191 INTEGER NTYPES, NTESTS
192 parameter( ntypes = 8, ntests = 7 )
193 INTEGER NBW
194 parameter( nbw = 4 )
195* ..
196* .. Local Scalars ..
197 LOGICAL ZEROT
198 CHARACTER DIST, PACKIT, TYPE, UPLO, XTYPE
199 CHARACTER*3 PATH
200 INTEGER I, I1, I2, IKD, IMAT, IN, INB, INFO, IOFF,
201 $ irhs, iuplo, iw, izero, k, kd, kl, koff, ku,
202 $ lda, ldab, mode, n, nb, nerrs, nfail, nimat,
203 $ nkd, nrhs, nrun
204 REAL AINVNM, ANORM, CNDNUM, RCOND, RCONDC
205* ..
206* .. Local Arrays ..
207 INTEGER ISEED( 4 ), ISEEDY( 4 ), KDVAL( NBW )
208 REAL RESULT( NTESTS )
209* ..
210* .. External Functions ..
211 REAL CLANGE, CLANHB, SGET06
212 EXTERNAL CLANGE, CLANHB, SGET06
213* ..
214* .. External Subroutines ..
215 EXTERNAL alaerh, alahd, alasum, ccopy, cerrpo, cget04,
219* ..
220* .. Intrinsic Functions ..
221 INTRINSIC cmplx, max, min
222* ..
223* .. Scalars in Common ..
224 LOGICAL LERR, OK
225 CHARACTER*32 SRNAMT
226 INTEGER INFOT, NUNIT
227* ..
228* .. Common blocks ..
229 COMMON / infoc / infot, nunit, ok, lerr
230 COMMON / srnamc / srnamt
231* ..
232* .. Data statements ..
233 DATA iseedy / 1988, 1989, 1990, 1991 /
234* ..
235* .. Executable Statements ..
236*
237* Initialize constants and the random number seed.
238*
239 path( 1: 1 ) = 'Complex precision'
240 path( 2: 3 ) = 'PB'
241 nrun = 0
242 nfail = 0
243 nerrs = 0
244 DO 10 i = 1, 4
245 iseed( i ) = iseedy( i )
246 10 CONTINUE
247*
248* Test the error exits
249*
250 IF( tsterr )
251 $ CALL cerrpo( path, nout )
252 infot = 0
253 kdval( 1 ) = 0
254*
255* Do for each value of N in NVAL
256*
257 DO 90 in = 1, nn
258 n = nval( in )
259 lda = max( n, 1 )
260 xtype = 'N'
261*
262* Set limits on the number of loop iterations.
263*
264 nkd = max( 1, min( n, 4 ) )
265 nimat = ntypes
266 IF( n.EQ.0 )
267 $ nimat = 1
268*
269 kdval( 2 ) = n + ( n+1 ) / 4
270 kdval( 3 ) = ( 3*n-1 ) / 4
271 kdval( 4 ) = ( n+1 ) / 4
272*
273 DO 80 ikd = 1, nkd
274*
275* Do for KD = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This order
276* makes it easier to skip redundant values for small values
277* of N.
278*
279 kd = kdval( ikd )
280 ldab = kd + 1
281*
282* Do first for UPLO = 'U', then for UPLO = 'L'
283*
284 DO 70 iuplo = 1, 2
285 koff = 1
286 IF( iuplo.EQ.1 ) THEN
287 uplo = 'U'
288 koff = max( 1, kd+2-n )
289 packit = 'Q'
290 ELSE
291 uplo = 'L'
292 packit = 'B'
293 END IF
294*
295 DO 60 imat = 1, nimat
296*
297* Do the tests only if DOTYPE( IMAT ) is true.
298*
299 IF( .NOT.dotype( imat ) )
300 $ GO TO 60
301*
302* Skip types 2, 3, or 4 if the matrix size is too small.
303*
304 zerot = imat.GE.2 .AND. imat.LE.4
305 IF( zerot .AND. n.LT.imat-1 )
306 $ GO TO 60
307*
308 IF( .NOT.zerot .OR. .NOT.dotype( 1 ) ) THEN
309*
310* Set up parameters with CLATB4 and generate a test
311* matrix with CLATMS.
312*
313 CALL clatb4( path, imat, n, n, TYPE, kl, ku, anorm,
314 $ mode, cndnum, dist )
315*
316 srnamt = 'CLATMS'
317 CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
318 $ cndnum, anorm, kd, kd, packit,
319 $ a( koff ), ldab, work, info )
320*
321* Check error code from CLATMS.
322*
323 IF( info.NE.0 ) THEN
324 CALL alaerh( path, 'CLATMS', info, 0, uplo, n,
325 $ n, kd, kd, -1, imat, nfail, nerrs,
326 $ nout )
327 GO TO 60
328 END IF
329 ELSE IF( izero.GT.0 ) THEN
330*
331* Use the same matrix for types 3 and 4 as for type
332* 2 by copying back the zeroed out column,
333*
334 iw = 2*lda + 1
335 IF( iuplo.EQ.1 ) THEN
336 ioff = ( izero-1 )*ldab + kd + 1
337 CALL ccopy( izero-i1, work( iw ), 1,
338 $ a( ioff-izero+i1 ), 1 )
339 iw = iw + izero - i1
340 CALL ccopy( i2-izero+1, work( iw ), 1,
341 $ a( ioff ), max( ldab-1, 1 ) )
342 ELSE
343 ioff = ( i1-1 )*ldab + 1
344 CALL ccopy( izero-i1, work( iw ), 1,
345 $ a( ioff+izero-i1 ),
346 $ max( ldab-1, 1 ) )
347 ioff = ( izero-1 )*ldab + 1
348 iw = iw + izero - i1
349 CALL ccopy( i2-izero+1, work( iw ), 1,
350 $ a( ioff ), 1 )
351 END IF
352 END IF
353*
354* For types 2-4, zero one row and column of the matrix
355* to test that INFO is returned correctly.
356*
357 izero = 0
358 IF( zerot ) THEN
359 IF( imat.EQ.2 ) THEN
360 izero = 1
361 ELSE IF( imat.EQ.3 ) THEN
362 izero = n
363 ELSE
364 izero = n / 2 + 1
365 END IF
366*
367* Save the zeroed out row and column in WORK(*,3)
368*
369 iw = 2*lda
370 DO 20 i = 1, min( 2*kd+1, n )
371 work( iw+i ) = zero
372 20 CONTINUE
373 iw = iw + 1
374 i1 = max( izero-kd, 1 )
375 i2 = min( izero+kd, n )
376*
377 IF( iuplo.EQ.1 ) THEN
378 ioff = ( izero-1 )*ldab + kd + 1
379 CALL cswap( izero-i1, a( ioff-izero+i1 ), 1,
380 $ work( iw ), 1 )
381 iw = iw + izero - i1
382 CALL cswap( i2-izero+1, a( ioff ),
383 $ max( ldab-1, 1 ), work( iw ), 1 )
384 ELSE
385 ioff = ( i1-1 )*ldab + 1
386 CALL cswap( izero-i1, a( ioff+izero-i1 ),
387 $ max( ldab-1, 1 ), work( iw ), 1 )
388 ioff = ( izero-1 )*ldab + 1
389 iw = iw + izero - i1
390 CALL cswap( i2-izero+1, a( ioff ), 1,
391 $ work( iw ), 1 )
392 END IF
393 END IF
394*
395* Set the imaginary part of the diagonals.
396*
397 IF( iuplo.EQ.1 ) THEN
398 CALL claipd( n, a( kd+1 ), ldab, 0 )
399 ELSE
400 CALL claipd( n, a( 1 ), ldab, 0 )
401 END IF
402*
403* Do for each value of NB in NBVAL
404*
405 DO 50 inb = 1, nnb
406 nb = nbval( inb )
407 CALL xlaenv( 1, nb )
408*
409* Compute the L*L' or U'*U factorization of the band
410* matrix.
411*
412 CALL clacpy( 'Full', kd+1, n, a, ldab, afac, ldab )
413 srnamt = 'CPBTRF'
414 CALL cpbtrf( uplo, n, kd, afac, ldab, info )
415*
416* Check error code from CPBTRF.
417*
418 IF( info.NE.izero ) THEN
419 CALL alaerh( path, 'CPBTRF', info, izero, uplo,
420 $ n, n, kd, kd, nb, imat, nfail,
421 $ nerrs, nout )
422 GO TO 50
423 END IF
424*
425* Skip the tests if INFO is not 0.
426*
427 IF( info.NE.0 )
428 $ GO TO 50
429*
430*+ TEST 1
431* Reconstruct matrix from factors and compute
432* residual.
433*
434 CALL clacpy( 'Full', kd+1, n, afac, ldab, ainv,
435 $ ldab )
436 CALL cpbt01( uplo, n, kd, a, ldab, ainv, ldab,
437 $ rwork, result( 1 ) )
438*
439* Print the test ratio if it is .GE. THRESH.
440*
441 IF( result( 1 ).GE.thresh ) THEN
442 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
443 $ CALL alahd( nout, path )
444 WRITE( nout, fmt = 9999 )uplo, n, kd, nb, imat,
445 $ 1, result( 1 )
446 nfail = nfail + 1
447 END IF
448 nrun = nrun + 1
449*
450* Only do other tests if this is the first blocksize.
451*
452 IF( inb.GT.1 )
453 $ GO TO 50
454*
455* Form the inverse of A so we can get a good estimate
456* of RCONDC = 1/(norm(A) * norm(inv(A))).
457*
458 CALL claset( 'Full', n, n, cmplx( zero ),
459 $ cmplx( one ), ainv, lda )
460 srnamt = 'CPBTRS'
461 CALL cpbtrs( uplo, n, kd, n, afac, ldab, ainv, lda,
462 $ info )
463*
464* Compute RCONDC = 1/(norm(A) * norm(inv(A))).
465*
466 anorm = clanhb( '1', uplo, n, kd, a, ldab, rwork )
467 ainvnm = clange( '1', n, n, ainv, lda, rwork )
468 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
469 rcondc = one
470 ELSE
471 rcondc = ( one / anorm ) / ainvnm
472 END IF
473*
474 DO 40 irhs = 1, nns
475 nrhs = nsval( irhs )
476*
477*+ TEST 2
478* Solve and compute residual for A * X = B.
479*
480 srnamt = 'CLARHS'
481 CALL clarhs( path, xtype, uplo, ' ', n, n, kd,
482 $ kd, nrhs, a, ldab, xact, lda, b,
483 $ lda, iseed, info )
484 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
485*
486 srnamt = 'CPBTRS'
487 CALL cpbtrs( uplo, n, kd, nrhs, afac, ldab, x,
488 $ lda, info )
489*
490* Check error code from CPBTRS.
491*
492 IF( info.NE.0 )
493 $ CALL alaerh( path, 'CPBTRS', info, 0, uplo,
494 $ n, n, kd, kd, nrhs, imat, nfail,
495 $ nerrs, nout )
496*
497 CALL clacpy( 'Full', n, nrhs, b, lda, work,
498 $ lda )
499 CALL cpbt02( uplo, n, kd, nrhs, a, ldab, x, lda,
500 $ work, lda, rwork, result( 2 ) )
501*
502*+ TEST 3
503* Check solution from generated exact solution.
504*
505 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
506 $ result( 3 ) )
507*
508*+ TESTS 4, 5, and 6
509* Use iterative refinement to improve the solution.
510*
511 srnamt = 'CPBRFS'
512 CALL cpbrfs( uplo, n, kd, nrhs, a, ldab, afac,
513 $ ldab, b, lda, x, lda, rwork,
514 $ rwork( nrhs+1 ), work,
515 $ rwork( 2*nrhs+1 ), info )
516*
517* Check error code from CPBRFS.
518*
519 IF( info.NE.0 )
520 $ CALL alaerh( path, 'CPBRFS', info, 0, uplo,
521 $ n, n, kd, kd, nrhs, imat, nfail,
522 $ nerrs, nout )
523*
524 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
525 $ result( 4 ) )
526 CALL cpbt05( uplo, n, kd, nrhs, a, ldab, b, lda,
527 $ x, lda, xact, lda, rwork,
528 $ rwork( nrhs+1 ), result( 5 ) )
529*
530* Print information about the tests that did not
531* pass the threshold.
532*
533 DO 30 k = 2, 6
534 IF( result( k ).GE.thresh ) THEN
535 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
536 $ CALL alahd( nout, path )
537 WRITE( nout, fmt = 9998 )uplo, n, kd,
538 $ nrhs, imat, k, result( k )
539 nfail = nfail + 1
540 END IF
541 30 CONTINUE
542 nrun = nrun + 5
543 40 CONTINUE
544*
545*+ TEST 7
546* Get an estimate of RCOND = 1/CNDNUM.
547*
548 srnamt = 'CPBCON'
549 CALL cpbcon( uplo, n, kd, afac, ldab, anorm, rcond,
550 $ work, rwork, info )
551*
552* Check error code from CPBCON.
553*
554 IF( info.NE.0 )
555 $ CALL alaerh( path, 'CPBCON', info, 0, uplo, n,
556 $ n, kd, kd, -1, imat, nfail, nerrs,
557 $ nout )
558*
559 result( 7 ) = sget06( rcond, rcondc )
560*
561* Print the test ratio if it is .GE. THRESH.
562*
563 IF( result( 7 ).GE.thresh ) THEN
564 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
565 $ CALL alahd( nout, path )
566 WRITE( nout, fmt = 9997 )uplo, n, kd, imat, 7,
567 $ result( 7 )
568 nfail = nfail + 1
569 END IF
570 nrun = nrun + 1
571 50 CONTINUE
572 60 CONTINUE
573 70 CONTINUE
574 80 CONTINUE
575 90 CONTINUE
576*
577* Print a summary of the results.
578*
579 CALL alasum( path, nout, nfail, nrun, nerrs )
580*
581 9999 FORMAT( ' UPLO=''', a1, ''', N=', i5, ', KD=', i5, ', NB=', i4,
582 $ ', type ', i2, ', test ', i2, ', ratio= ', g12.5 )
583 9998 FORMAT( ' UPLO=''', a1, ''', N=', i5, ', KD=', i5, ', NRHS=', i3,
584 $ ', type ', i2, ', test(', i2, ') = ', g12.5 )
585 9997 FORMAT( ' UPLO=''', a1, ''', N=', i5, ', KD=', i5, ',', 10x,
586 $ ' type ', i2, ', test(', i2, ') = ', g12.5 )
587 RETURN
588*
589* End of CCHKPB
590*
591 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 cchkpb(dotype, nn, nval, nnb, nbval, nns, nsval, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, nout)
CCHKPB
Definition cchkpb.f:168
subroutine cerrpo(path, nunit)
CERRPO
Definition cerrpo.f:55
subroutine cget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
CGET04
Definition cget04.f:102
subroutine claipd(n, a, inda, vinda)
CLAIPD
Definition claipd.f:83
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 cpbt01(uplo, n, kd, a, lda, afac, ldafac, rwork, resid)
CPBT01
Definition cpbt01.f:120
subroutine cpbt02(uplo, n, kd, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
CPBT02
Definition cpbt02.f:136
subroutine cpbt05(uplo, n, kd, nrhs, ab, ldab, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
CPBT05
Definition cpbt05.f:171
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
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 cpbcon(uplo, n, kd, ab, ldab, anorm, rcond, work, rwork, info)
CPBCON
Definition cpbcon.f:133
subroutine cpbrfs(uplo, n, kd, nrhs, ab, ldab, afb, ldafb, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CPBRFS
Definition cpbrfs.f:189
subroutine cpbtrf(uplo, n, kd, ab, ldab, info)
CPBTRF
Definition cpbtrf.f:142
subroutine cpbtrs(uplo, n, kd, nrhs, ab, ldab, b, ldb, info)
CPBTRS
Definition cpbtrs.f:121
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81