LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine schkpb ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
real  THRESH,
logical  TSTERR,
integer  NMAX,
real, dimension( * )  A,
real, dimension( * )  AFAC,
real, dimension( * )  AINV,
real, dimension( * )  B,
real, dimension( * )  X,
real, dimension( * )  XACT,
real, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

SCHKPB

Purpose:
 SCHKPB tests SPBTRF, -TRS, -RFS, and -CON.
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 dimension N.
[in]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NBVAL)
          The values of the blocksize NB.
[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 maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]AINV
          AINV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is REAL array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is REAL array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is REAL array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL 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.
Date
November 2011

Definition at line 174 of file schkpb.f.

174 *
175 * -- LAPACK test routine (version 3.4.0) --
176 * -- LAPACK is a software package provided by Univ. of Tennessee, --
177 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178 * November 2011
179 *
180 * .. Scalar Arguments ..
181  LOGICAL tsterr
182  INTEGER nmax, nn, nnb, nns, nout
183  REAL thresh
184 * ..
185 * .. Array Arguments ..
186  LOGICAL dotype( * )
187  INTEGER iwork( * ), nbval( * ), nsval( * ), nval( * )
188  REAL a( * ), afac( * ), ainv( * ), b( * ),
189  $ rwork( * ), work( * ), x( * ), xact( * )
190 * ..
191 *
192 * =====================================================================
193 *
194 * .. Parameters ..
195  REAL one, zero
196  parameter ( one = 1.0e+0, zero = 0.0e+0 )
197  INTEGER ntypes, ntests
198  parameter ( ntypes = 8, ntests = 7 )
199  INTEGER nbw
200  parameter ( nbw = 4 )
201 * ..
202 * .. Local Scalars ..
203  LOGICAL zerot
204  CHARACTER dist, packit, TYPE, uplo, xtype
205  CHARACTER*3 path
206  INTEGER i, i1, i2, ikd, imat, in, inb, info, ioff,
207  $ irhs, iuplo, iw, izero, k, kd, kl, koff, ku,
208  $ lda, ldab, mode, n, nb, nerrs, nfail, nimat,
209  $ nkd, nrhs, nrun
210  REAL ainvnm, anorm, cndnum, rcond, rcondc
211 * ..
212 * .. Local Arrays ..
213  INTEGER iseed( 4 ), iseedy( 4 ), kdval( nbw )
214  REAL result( ntests )
215 * ..
216 * .. External Functions ..
217  REAL sget06, slange, slansb
218  EXTERNAL sget06, slange, slansb
219 * ..
220 * .. External Subroutines ..
221  EXTERNAL alaerh, alahd, alasum, scopy, serrpo, sget04,
224  $ sswap, xlaenv
225 * ..
226 * .. Intrinsic Functions ..
227  INTRINSIC max, min
228 * ..
229 * .. Scalars in Common ..
230  LOGICAL lerr, ok
231  CHARACTER*32 srnamt
232  INTEGER infot, nunit
233 * ..
234 * .. Common blocks ..
235  COMMON / infoc / infot, nunit, ok, lerr
236  COMMON / srnamc / srnamt
237 * ..
238 * .. Data statements ..
239  DATA iseedy / 1988, 1989, 1990, 1991 /
240 * ..
241 * .. Executable Statements ..
242 *
243 * Initialize constants and the random number seed.
244 *
245  path( 1: 1 ) = 'Single precision'
246  path( 2: 3 ) = 'PB'
247  nrun = 0
248  nfail = 0
249  nerrs = 0
250  DO 10 i = 1, 4
251  iseed( i ) = iseedy( i )
252  10 CONTINUE
253 *
254 * Test the error exits
255 *
256  IF( tsterr )
257  $ CALL serrpo( path, nout )
258  infot = 0
259  CALL xlaenv( 2, 2 )
260  kdval( 1 ) = 0
261 *
262 * Do for each value of N in NVAL
263 *
264  DO 90 in = 1, nn
265  n = nval( in )
266  lda = max( n, 1 )
267  xtype = 'N'
268 *
269 * Set limits on the number of loop iterations.
270 *
271  nkd = max( 1, min( n, 4 ) )
272  nimat = ntypes
273  IF( n.EQ.0 )
274  $ nimat = 1
275 *
276  kdval( 2 ) = n + ( n+1 ) / 4
277  kdval( 3 ) = ( 3*n-1 ) / 4
278  kdval( 4 ) = ( n+1 ) / 4
279 *
280  DO 80 ikd = 1, nkd
281 *
282 * Do for KD = 0, (5*N+1)/4, (3N-1)/4, and (N+1)/4. This order
283 * makes it easier to skip redundant values for small values
284 * of N.
285 *
286  kd = kdval( ikd )
287  ldab = kd + 1
288 *
289 * Do first for UPLO = 'U', then for UPLO = 'L'
290 *
291  DO 70 iuplo = 1, 2
292  koff = 1
293  IF( iuplo.EQ.1 ) THEN
294  uplo = 'U'
295  koff = max( 1, kd+2-n )
296  packit = 'Q'
297  ELSE
298  uplo = 'L'
299  packit = 'B'
300  END IF
301 *
302  DO 60 imat = 1, nimat
303 *
304 * Do the tests only if DOTYPE( IMAT ) is true.
305 *
306  IF( .NOT.dotype( imat ) )
307  $ GO TO 60
308 *
309 * Skip types 2, 3, or 4 if the matrix size is too small.
310 *
311  zerot = imat.GE.2 .AND. imat.LE.4
312  IF( zerot .AND. n.LT.imat-1 )
313  $ GO TO 60
314 *
315  IF( .NOT.zerot .OR. .NOT.dotype( 1 ) ) THEN
316 *
317 * Set up parameters with SLATB4 and generate a test
318 * matrix with SLATMS.
319 *
320  CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm,
321  $ mode, cndnum, dist )
322 *
323  srnamt = 'SLATMS'
324  CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
325  $ cndnum, anorm, kd, kd, packit,
326  $ a( koff ), ldab, work, info )
327 *
328 * Check error code from SLATMS.
329 *
330  IF( info.NE.0 ) THEN
331  CALL alaerh( path, 'SLATMS', info, 0, uplo, n,
332  $ n, kd, kd, -1, imat, nfail, nerrs,
333  $ nout )
334  GO TO 60
335  END IF
336  ELSE IF( izero.GT.0 ) THEN
337 *
338 * Use the same matrix for types 3 and 4 as for type
339 * 2 by copying back the zeroed out column,
340 *
341  iw = 2*lda + 1
342  IF( iuplo.EQ.1 ) THEN
343  ioff = ( izero-1 )*ldab + kd + 1
344  CALL scopy( izero-i1, work( iw ), 1,
345  $ a( ioff-izero+i1 ), 1 )
346  iw = iw + izero - i1
347  CALL scopy( i2-izero+1, work( iw ), 1,
348  $ a( ioff ), max( ldab-1, 1 ) )
349  ELSE
350  ioff = ( i1-1 )*ldab + 1
351  CALL scopy( izero-i1, work( iw ), 1,
352  $ a( ioff+izero-i1 ),
353  $ max( ldab-1, 1 ) )
354  ioff = ( izero-1 )*ldab + 1
355  iw = iw + izero - i1
356  CALL scopy( i2-izero+1, work( iw ), 1,
357  $ a( ioff ), 1 )
358  END IF
359  END IF
360 *
361 * For types 2-4, zero one row and column of the matrix
362 * to test that INFO is returned correctly.
363 *
364  izero = 0
365  IF( zerot ) THEN
366  IF( imat.EQ.2 ) THEN
367  izero = 1
368  ELSE IF( imat.EQ.3 ) THEN
369  izero = n
370  ELSE
371  izero = n / 2 + 1
372  END IF
373 *
374 * Save the zeroed out row and column in WORK(*,3)
375 *
376  iw = 2*lda
377  DO 20 i = 1, min( 2*kd+1, n )
378  work( iw+i ) = zero
379  20 CONTINUE
380  iw = iw + 1
381  i1 = max( izero-kd, 1 )
382  i2 = min( izero+kd, n )
383 *
384  IF( iuplo.EQ.1 ) THEN
385  ioff = ( izero-1 )*ldab + kd + 1
386  CALL sswap( izero-i1, a( ioff-izero+i1 ), 1,
387  $ work( iw ), 1 )
388  iw = iw + izero - i1
389  CALL sswap( i2-izero+1, a( ioff ),
390  $ max( ldab-1, 1 ), work( iw ), 1 )
391  ELSE
392  ioff = ( i1-1 )*ldab + 1
393  CALL sswap( izero-i1, a( ioff+izero-i1 ),
394  $ max( ldab-1, 1 ), work( iw ), 1 )
395  ioff = ( izero-1 )*ldab + 1
396  iw = iw + izero - i1
397  CALL sswap( i2-izero+1, a( ioff ), 1,
398  $ work( iw ), 1 )
399  END IF
400  END IF
401 *
402 * Do for each value of NB in NBVAL
403 *
404  DO 50 inb = 1, nnb
405  nb = nbval( inb )
406  CALL xlaenv( 1, nb )
407 *
408 * Compute the L*L' or U'*U factorization of the band
409 * matrix.
410 *
411  CALL slacpy( 'Full', kd+1, n, a, ldab, afac, ldab )
412  srnamt = 'SPBTRF'
413  CALL spbtrf( uplo, n, kd, afac, ldab, info )
414 *
415 * Check error code from SPBTRF.
416 *
417  IF( info.NE.izero ) THEN
418  CALL alaerh( path, 'SPBTRF', info, izero, uplo,
419  $ n, n, kd, kd, nb, imat, nfail,
420  $ nerrs, nout )
421  GO TO 50
422  END IF
423 *
424 * Skip the tests if INFO is not 0.
425 *
426  IF( info.NE.0 )
427  $ GO TO 50
428 *
429 *+ TEST 1
430 * Reconstruct matrix from factors and compute
431 * residual.
432 *
433  CALL slacpy( 'Full', kd+1, n, afac, ldab, ainv,
434  $ ldab )
435  CALL spbt01( uplo, n, kd, a, ldab, ainv, ldab,
436  $ rwork, result( 1 ) )
437 *
438 * Print the test ratio if it is .GE. THRESH.
439 *
440  IF( result( 1 ).GE.thresh ) THEN
441  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
442  $ CALL alahd( nout, path )
443  WRITE( nout, fmt = 9999 )uplo, n, kd, nb, imat,
444  $ 1, result( 1 )
445  nfail = nfail + 1
446  END IF
447  nrun = nrun + 1
448 *
449 * Only do other tests if this is the first blocksize.
450 *
451  IF( inb.GT.1 )
452  $ GO TO 50
453 *
454 * Form the inverse of A so we can get a good estimate
455 * of RCONDC = 1/(norm(A) * norm(inv(A))).
456 *
457  CALL slaset( 'Full', n, n, zero, one, ainv, lda )
458  srnamt = 'SPBTRS'
459  CALL spbtrs( uplo, n, kd, n, afac, ldab, ainv, lda,
460  $ info )
461 *
462 * Compute RCONDC = 1/(norm(A) * norm(inv(A))).
463 *
464  anorm = slansb( '1', uplo, n, kd, a, ldab, rwork )
465  ainvnm = slange( '1', n, n, ainv, lda, rwork )
466  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
467  rcondc = one
468  ELSE
469  rcondc = ( one / anorm ) / ainvnm
470  END IF
471 *
472  DO 40 irhs = 1, nns
473  nrhs = nsval( irhs )
474 *
475 *+ TEST 2
476 * Solve and compute residual for A * X = B.
477 *
478  srnamt = 'SLARHS'
479  CALL slarhs( path, xtype, uplo, ' ', n, n, kd,
480  $ kd, nrhs, a, ldab, xact, lda, b,
481  $ lda, iseed, info )
482  CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
483 *
484  srnamt = 'SPBTRS'
485  CALL spbtrs( uplo, n, kd, nrhs, afac, ldab, x,
486  $ lda, info )
487 *
488 * Check error code from SPBTRS.
489 *
490  IF( info.NE.0 )
491  $ CALL alaerh( path, 'SPBTRS', info, 0, uplo,
492  $ n, n, kd, kd, nrhs, imat, nfail,
493  $ nerrs, nout )
494 *
495  CALL slacpy( 'Full', n, nrhs, b, lda, work,
496  $ lda )
497  CALL spbt02( uplo, n, kd, nrhs, a, ldab, x, lda,
498  $ work, lda, rwork, result( 2 ) )
499 *
500 *+ TEST 3
501 * Check solution from generated exact solution.
502 *
503  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
504  $ result( 3 ) )
505 *
506 *+ TESTS 4, 5, and 6
507 * Use iterative refinement to improve the solution.
508 *
509  srnamt = 'SPBRFS'
510  CALL spbrfs( uplo, n, kd, nrhs, a, ldab, afac,
511  $ ldab, b, lda, x, lda, rwork,
512  $ rwork( nrhs+1 ), work, iwork,
513  $ info )
514 *
515 * Check error code from SPBRFS.
516 *
517  IF( info.NE.0 )
518  $ CALL alaerh( path, 'SPBRFS', info, 0, uplo,
519  $ n, n, kd, kd, nrhs, imat, nfail,
520  $ nerrs, nout )
521 *
522  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
523  $ result( 4 ) )
524  CALL spbt05( uplo, n, kd, nrhs, a, ldab, b, lda,
525  $ x, lda, xact, lda, rwork,
526  $ rwork( nrhs+1 ), result( 5 ) )
527 *
528 * Print information about the tests that did not
529 * pass the threshold.
530 *
531  DO 30 k = 2, 6
532  IF( result( k ).GE.thresh ) THEN
533  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
534  $ CALL alahd( nout, path )
535  WRITE( nout, fmt = 9998 )uplo, n, kd,
536  $ nrhs, imat, k, result( k )
537  nfail = nfail + 1
538  END IF
539  30 CONTINUE
540  nrun = nrun + 5
541  40 CONTINUE
542 *
543 *+ TEST 7
544 * Get an estimate of RCOND = 1/CNDNUM.
545 *
546  srnamt = 'SPBCON'
547  CALL spbcon( uplo, n, kd, afac, ldab, anorm, rcond,
548  $ work, iwork, info )
549 *
550 * Check error code from SPBCON.
551 *
552  IF( info.NE.0 )
553  $ CALL alaerh( path, 'SPBCON', info, 0, uplo, n,
554  $ n, kd, kd, -1, imat, nfail, nerrs,
555  $ nout )
556 *
557  result( 7 ) = sget06( rcond, rcondc )
558 *
559 * Print the test ratio if it is .GE. THRESH.
560 *
561  IF( result( 7 ).GE.thresh ) THEN
562  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
563  $ CALL alahd( nout, path )
564  WRITE( nout, fmt = 9997 )uplo, n, kd, imat, 7,
565  $ result( 7 )
566  nfail = nfail + 1
567  END IF
568  nrun = nrun + 1
569  50 CONTINUE
570  60 CONTINUE
571  70 CONTINUE
572  80 CONTINUE
573  90 CONTINUE
574 *
575 * Print a summary of the results.
576 *
577  CALL alasum( path, nout, nfail, nrun, nerrs )
578 *
579  9999 FORMAT( ' UPLO=''', a1, ''', N=', i5, ', KD=', i5, ', NB=', i4,
580  $ ', type ', i2, ', test ', i2, ', ratio= ', g12.5 )
581  9998 FORMAT( ' UPLO=''', a1, ''', N=', i5, ', KD=', i5, ', NRHS=', i3,
582  $ ', type ', i2, ', test(', i2, ') = ', g12.5 )
583  9997 FORMAT( ' UPLO=''', a1, ''', N=', i5, ', KD=', i5, ',', 10x,
584  $ ' type ', i2, ', test(', i2, ') = ', g12.5 )
585  RETURN
586 *
587 * End of SCHKPB
588 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine slatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
SLATB4
Definition: slatb4.f:122
subroutine spbt02(UPLO, N, KD, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
SPBT02
Definition: spbt02.f:138
subroutine slarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
SLARHS
Definition: slarhs.f:206
subroutine spbcon(UPLO, N, KD, AB, LDAB, ANORM, RCOND, WORK, IWORK, INFO)
SPBCON
Definition: spbcon.f:134
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine spbtrs(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
SPBTRS
Definition: spbtrs.f:123
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
subroutine spbrfs(UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
SPBRFS
Definition: spbrfs.f:191
real function slansb(NORM, UPLO, N, K, AB, LDAB, WORK)
SLANSB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric band matrix.
Definition: slansb.f:131
real function slange(NORM, M, N, A, LDA, WORK)
SLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: slange.f:116
subroutine sget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
SGET04
Definition: sget04.f:104
subroutine spbtrf(UPLO, N, KD, AB, LDAB, INFO)
SPBTRF
Definition: spbtrf.f:144
subroutine spbt05(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
SPBT05
Definition: spbt05.f:173
subroutine spbt01(UPLO, N, KD, A, LDA, AFAC, LDAFAC, RWORK, RESID)
SPBT01
Definition: spbt01.f:121
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP
Definition: sswap.f:53
subroutine serrpo(PATH, NUNIT)
SERRPO
Definition: serrpo.f:57
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75

Here is the call graph for this function:

Here is the caller graph for this function: