LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dchksp ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
double precision, dimension( * )  A,
double precision, dimension( * )  AFAC,
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 
)

DCHKSP

Purpose:
 DCHKSP tests DSPTRF, -TRI, -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]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 maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is DOUBLE PRECISION array, dimension
                      (NMAX*(NMAX+1)/2)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension
                      (NMAX*(NMAX+1)/2)
[out]AINV
          AINV is DOUBLE PRECISION array, dimension
                      (NMAX*(NMAX+1)/2)
[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(2,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array,
                                 dimension (NMAX+2*NSMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (2*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 165 of file dchksp.f.

165 *
166 * -- LAPACK test routine (version 3.4.0) --
167 * -- LAPACK is a software package provided by Univ. of Tennessee, --
168 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
169 * November 2011
170 *
171 * .. Scalar Arguments ..
172  LOGICAL tsterr
173  INTEGER nmax, nn, nns, nout
174  DOUBLE PRECISION thresh
175 * ..
176 * .. Array Arguments ..
177  LOGICAL dotype( * )
178  INTEGER iwork( * ), nsval( * ), nval( * )
179  DOUBLE PRECISION a( * ), afac( * ), ainv( * ), b( * ),
180  $ rwork( * ), work( * ), x( * ), xact( * )
181 * ..
182 *
183 * =====================================================================
184 *
185 * .. Parameters ..
186  DOUBLE PRECISION zero
187  parameter ( zero = 0.0d+0 )
188  INTEGER ntypes
189  parameter ( ntypes = 10 )
190  INTEGER ntests
191  parameter ( ntests = 8 )
192 * ..
193 * .. Local Scalars ..
194  LOGICAL trfcon, zerot
195  CHARACTER dist, packit, TYPE, uplo, xtype
196  CHARACTER*3 path
197  INTEGER i, i1, i2, imat, in, info, ioff, irhs, iuplo,
198  $ izero, j, k, kl, ku, lda, mode, n, nerrs,
199  $ nfail, nimat, npp, nrhs, nrun, nt
200  DOUBLE PRECISION anorm, cndnum, rcond, rcondc
201 * ..
202 * .. Local Arrays ..
203  CHARACTER uplos( 2 )
204  INTEGER iseed( 4 ), iseedy( 4 )
205  DOUBLE PRECISION result( ntests )
206 * ..
207 * .. External Functions ..
208  LOGICAL lsame
209  DOUBLE PRECISION dget06, dlansp
210  EXTERNAL lsame, dget06, dlansp
211 * ..
212 * .. External Subroutines ..
213  EXTERNAL alaerh, alahd, alasum, dcopy, derrsy, dget04,
216  $ dsptrs
217 * ..
218 * .. Intrinsic Functions ..
219  INTRINSIC max, min
220 * ..
221 * .. Scalars in Common ..
222  LOGICAL lerr, ok
223  CHARACTER*32 srnamt
224  INTEGER infot, nunit
225 * ..
226 * .. Common blocks ..
227  COMMON / infoc / infot, nunit, ok, lerr
228  COMMON / srnamc / srnamt
229 * ..
230 * .. Data statements ..
231  DATA iseedy / 1988, 1989, 1990, 1991 /
232  DATA uplos / 'U', 'L' /
233 * ..
234 * .. Executable Statements ..
235 *
236 * Initialize constants and the random number seed.
237 *
238  path( 1: 1 ) = 'Double precision'
239  path( 2: 3 ) = 'SP'
240  nrun = 0
241  nfail = 0
242  nerrs = 0
243  DO 10 i = 1, 4
244  iseed( i ) = iseedy( i )
245  10 CONTINUE
246 *
247 * Test the error exits
248 *
249  IF( tsterr )
250  $ CALL derrsy( path, nout )
251  infot = 0
252 *
253 * Do for each value of N in NVAL
254 *
255  DO 170 in = 1, nn
256  n = nval( in )
257  lda = max( n, 1 )
258  xtype = 'N'
259  nimat = ntypes
260  IF( n.LE.0 )
261  $ nimat = 1
262 *
263  izero = 0
264  DO 160 imat = 1, nimat
265 *
266 * Do the tests only if DOTYPE( IMAT ) is true.
267 *
268  IF( .NOT.dotype( imat ) )
269  $ GO TO 160
270 *
271 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
272 *
273  zerot = imat.GE.3 .AND. imat.LE.6
274  IF( zerot .AND. n.LT.imat-2 )
275  $ GO TO 160
276 *
277 * Do first for UPLO = 'U', then for UPLO = 'L'
278 *
279  DO 150 iuplo = 1, 2
280  uplo = uplos( iuplo )
281  IF( lsame( uplo, 'U' ) ) THEN
282  packit = 'C'
283  ELSE
284  packit = 'R'
285  END IF
286 *
287 * Set up parameters with DLATB4 and generate a test matrix
288 * with DLATMS.
289 *
290  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
291  $ cndnum, dist )
292 *
293  srnamt = 'DLATMS'
294  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
295  $ cndnum, anorm, kl, ku, packit, a, lda, work,
296  $ info )
297 *
298 * Check error code from DLATMS.
299 *
300  IF( info.NE.0 ) THEN
301  CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
302  $ -1, -1, imat, nfail, nerrs, nout )
303  GO TO 150
304  END IF
305 *
306 * For types 3-6, zero one or more rows and columns of
307 * the matrix to test that INFO is returned correctly.
308 *
309  IF( zerot ) THEN
310  IF( imat.EQ.3 ) THEN
311  izero = 1
312  ELSE IF( imat.EQ.4 ) THEN
313  izero = n
314  ELSE
315  izero = n / 2 + 1
316  END IF
317 *
318  IF( imat.LT.6 ) THEN
319 *
320 * Set row and column IZERO to zero.
321 *
322  IF( iuplo.EQ.1 ) THEN
323  ioff = ( izero-1 )*izero / 2
324  DO 20 i = 1, izero - 1
325  a( ioff+i ) = zero
326  20 CONTINUE
327  ioff = ioff + izero
328  DO 30 i = izero, n
329  a( ioff ) = zero
330  ioff = ioff + i
331  30 CONTINUE
332  ELSE
333  ioff = izero
334  DO 40 i = 1, izero - 1
335  a( ioff ) = zero
336  ioff = ioff + n - i
337  40 CONTINUE
338  ioff = ioff - izero
339  DO 50 i = izero, n
340  a( ioff+i ) = zero
341  50 CONTINUE
342  END IF
343  ELSE
344  ioff = 0
345  IF( iuplo.EQ.1 ) THEN
346 *
347 * Set the first IZERO rows and columns to zero.
348 *
349  DO 70 j = 1, n
350  i2 = min( j, izero )
351  DO 60 i = 1, i2
352  a( ioff+i ) = zero
353  60 CONTINUE
354  ioff = ioff + j
355  70 CONTINUE
356  ELSE
357 *
358 * Set the last IZERO rows and columns to zero.
359 *
360  DO 90 j = 1, n
361  i1 = max( j, izero )
362  DO 80 i = i1, n
363  a( ioff+i ) = zero
364  80 CONTINUE
365  ioff = ioff + n - j
366  90 CONTINUE
367  END IF
368  END IF
369  ELSE
370  izero = 0
371  END IF
372 *
373 * Compute the L*D*L' or U*D*U' factorization of the matrix.
374 *
375  npp = n*( n+1 ) / 2
376  CALL dcopy( npp, a, 1, afac, 1 )
377  srnamt = 'DSPTRF'
378  CALL dsptrf( uplo, n, afac, iwork, info )
379 *
380 * Adjust the expected value of INFO to account for
381 * pivoting.
382 *
383  k = izero
384  IF( k.GT.0 ) THEN
385  100 CONTINUE
386  IF( iwork( k ).LT.0 ) THEN
387  IF( iwork( k ).NE.-k ) THEN
388  k = -iwork( k )
389  GO TO 100
390  END IF
391  ELSE IF( iwork( k ).NE.k ) THEN
392  k = iwork( k )
393  GO TO 100
394  END IF
395  END IF
396 *
397 * Check error code from DSPTRF.
398 *
399  IF( info.NE.k )
400  $ CALL alaerh( path, 'DSPTRF', info, k, uplo, n, n, -1,
401  $ -1, -1, imat, nfail, nerrs, nout )
402  IF( info.NE.0 ) THEN
403  trfcon = .true.
404  ELSE
405  trfcon = .false.
406  END IF
407 *
408 *+ TEST 1
409 * Reconstruct matrix from factors and compute residual.
410 *
411  CALL dspt01( uplo, n, a, afac, iwork, ainv, lda, rwork,
412  $ result( 1 ) )
413  nt = 1
414 *
415 *+ TEST 2
416 * Form the inverse and compute the residual.
417 *
418  IF( .NOT.trfcon ) THEN
419  CALL dcopy( npp, afac, 1, ainv, 1 )
420  srnamt = 'DSPTRI'
421  CALL dsptri( uplo, n, ainv, iwork, work, info )
422 *
423 * Check error code from DSPTRI.
424 *
425  IF( info.NE.0 )
426  $ CALL alaerh( path, 'DSPTRI', info, 0, uplo, n, n,
427  $ -1, -1, -1, imat, nfail, nerrs, nout )
428 *
429  CALL dppt03( uplo, n, a, ainv, work, lda, rwork,
430  $ rcondc, result( 2 ) )
431  nt = 2
432  END IF
433 *
434 * Print information about the tests that did not pass
435 * the threshold.
436 *
437  DO 110 k = 1, nt
438  IF( result( k ).GE.thresh ) THEN
439  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
440  $ CALL alahd( nout, path )
441  WRITE( nout, fmt = 9999 )uplo, n, imat, k,
442  $ result( k )
443  nfail = nfail + 1
444  END IF
445  110 CONTINUE
446  nrun = nrun + nt
447 *
448 * Do only the condition estimate if INFO is not 0.
449 *
450  IF( trfcon ) THEN
451  rcondc = zero
452  GO TO 140
453  END IF
454 *
455  DO 130 irhs = 1, nns
456  nrhs = nsval( irhs )
457 *
458 *+ TEST 3
459 * Solve and compute residual for A * X = B.
460 *
461  srnamt = 'DLARHS'
462  CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
463  $ nrhs, a, lda, xact, lda, b, lda, iseed,
464  $ info )
465  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
466 *
467  srnamt = 'DSPTRS'
468  CALL dsptrs( uplo, n, nrhs, afac, iwork, x, lda,
469  $ info )
470 *
471 * Check error code from DSPTRS.
472 *
473  IF( info.NE.0 )
474  $ CALL alaerh( path, 'DSPTRS', info, 0, uplo, n, n,
475  $ -1, -1, nrhs, imat, nfail, nerrs,
476  $ nout )
477 *
478  CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
479  CALL dppt02( uplo, n, nrhs, a, x, lda, work, lda,
480  $ rwork, result( 3 ) )
481 *
482 *+ TEST 4
483 * Check solution from generated exact solution.
484 *
485  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
486  $ result( 4 ) )
487 *
488 *+ TESTS 5, 6, and 7
489 * Use iterative refinement to improve the solution.
490 *
491  srnamt = 'DSPRFS'
492  CALL dsprfs( uplo, n, nrhs, a, afac, iwork, b, lda, x,
493  $ lda, rwork, rwork( nrhs+1 ), work,
494  $ iwork( n+1 ), info )
495 *
496 * Check error code from DSPRFS.
497 *
498  IF( info.NE.0 )
499  $ CALL alaerh( path, 'DSPRFS', info, 0, uplo, n, n,
500  $ -1, -1, nrhs, imat, nfail, nerrs,
501  $ nout )
502 *
503  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
504  $ result( 5 ) )
505  CALL dppt05( uplo, n, nrhs, a, b, lda, x, lda, xact,
506  $ lda, rwork, rwork( nrhs+1 ),
507  $ result( 6 ) )
508 *
509 * Print information about the tests that did not pass
510 * the threshold.
511 *
512  DO 120 k = 3, 7
513  IF( result( k ).GE.thresh ) THEN
514  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
515  $ CALL alahd( nout, path )
516  WRITE( nout, fmt = 9998 )uplo, n, nrhs, imat,
517  $ k, result( k )
518  nfail = nfail + 1
519  END IF
520  120 CONTINUE
521  nrun = nrun + 5
522  130 CONTINUE
523 *
524 *+ TEST 8
525 * Get an estimate of RCOND = 1/CNDNUM.
526 *
527  140 CONTINUE
528  anorm = dlansp( '1', uplo, n, a, rwork )
529  srnamt = 'DSPCON'
530  CALL dspcon( uplo, n, afac, iwork, anorm, rcond, work,
531  $ iwork( n+1 ), info )
532 *
533 * Check error code from DSPCON.
534 *
535  IF( info.NE.0 )
536  $ CALL alaerh( path, 'DSPCON', info, 0, uplo, n, n, -1,
537  $ -1, -1, imat, nfail, nerrs, nout )
538 *
539  result( 8 ) = dget06( rcond, rcondc )
540 *
541 * Print the test ratio if it is .GE. THRESH.
542 *
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 = 9999 )uplo, n, imat, 8,
547  $ result( 8 )
548  nfail = nfail + 1
549  END IF
550  nrun = nrun + 1
551  150 CONTINUE
552  160 CONTINUE
553  170 CONTINUE
554 *
555 * Print a summary of the results.
556 *
557  CALL alasum( path, nout, nfail, nrun, nerrs )
558 *
559  9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', type ', i2, ', test ',
560  $ i2, ', ratio =', g12.5 )
561  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
562  $ i2, ', test(', i2, ') =', g12.5 )
563  RETURN
564 *
565 * End of DCHKSP
566 *
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 dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dsptrs(UPLO, N, NRHS, AP, IPIV, B, LDB, INFO)
DSPTRS
Definition: dsptrs.f:117
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dspt01(UPLO, N, A, AFAC, IPIV, C, LDC, RWORK, RESID)
DSPT01
Definition: dspt01.f:112
subroutine dsprfs(UPLO, N, NRHS, AP, AFP, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DSPRFS
Definition: dsprfs.f:181
subroutine dppt02(UPLO, N, NRHS, A, X, LDX, B, LDB, RWORK, RESID)
DPPT02
Definition: dppt02.f:124
subroutine dsptri(UPLO, N, AP, IPIV, WORK, INFO)
DSPTRI
Definition: dsptri.f:111
subroutine dppt03(UPLO, N, A, AINV, WORK, LDWORK, RWORK, RCOND, RESID)
DPPT03
Definition: dppt03.f:112
subroutine dspcon(UPLO, N, AP, IPIV, ANORM, RCOND, WORK, IWORK, INFO)
DSPCON
Definition: dspcon.f:127
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine derrsy(PATH, NUNIT)
DERRSY
Definition: derrsy.f:57
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dppt05(UPLO, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DPPT05
Definition: dppt05.f:158
subroutine dsptrf(UPLO, N, AP, IPIV, INFO)
DSPTRF
Definition: dsptrf.f:161
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
double precision function dlansp(NORM, UPLO, N, AP, WORK)
DLANSP 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 matrix supplied in packed form.
Definition: dlansp.f:116
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
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: