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

DCHKTP

Purpose:
 DCHKTP tests DTPTRI, -TRS, -RFS, and -CON, and DLATPS
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
          maximumm value of N in NVAL.
[out]AP
          AP is DOUBLE PRECISION array, dimension
                      (NMAX*(NMAX+1)/2)
[out]AINVP
          AINVP 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(3,NSMAX))
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(NMAX,2*NSMAX))
[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 159 of file dchktp.f.

159 *
160 * -- LAPACK test routine (version 3.4.0) --
161 * -- LAPACK is a software package provided by Univ. of Tennessee, --
162 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163 * November 2011
164 *
165 * .. Scalar Arguments ..
166  LOGICAL tsterr
167  INTEGER nmax, nn, nns, nout
168  DOUBLE PRECISION thresh
169 * ..
170 * .. Array Arguments ..
171  LOGICAL dotype( * )
172  INTEGER iwork( * ), nsval( * ), nval( * )
173  DOUBLE PRECISION ainvp( * ), ap( * ), b( * ), rwork( * ),
174  $ work( * ), x( * ), xact( * )
175 * ..
176 *
177 * =====================================================================
178 *
179 * .. Parameters ..
180  INTEGER ntype1, ntypes
181  parameter ( ntype1 = 10, ntypes = 18 )
182  INTEGER ntests
183  parameter ( ntests = 9 )
184  INTEGER ntran
185  parameter ( ntran = 3 )
186  DOUBLE PRECISION one, zero
187  parameter ( one = 1.0d+0, zero = 0.0d+0 )
188 * ..
189 * .. Local Scalars ..
190  CHARACTER diag, norm, trans, uplo, xtype
191  CHARACTER*3 path
192  INTEGER i, idiag, imat, in, info, irhs, itran, iuplo,
193  $ k, lap, lda, n, nerrs, nfail, nrhs, nrun
194  DOUBLE PRECISION ainvnm, anorm, rcond, rcondc, rcondi, rcondo,
195  $ scale
196 * ..
197 * .. Local Arrays ..
198  CHARACTER transs( ntran ), uplos( 2 )
199  INTEGER iseed( 4 ), iseedy( 4 )
200  DOUBLE PRECISION result( ntests )
201 * ..
202 * .. External Functions ..
203  LOGICAL lsame
204  DOUBLE PRECISION dlantp
205  EXTERNAL lsame, dlantp
206 * ..
207 * .. External Subroutines ..
208  EXTERNAL alaerh, alahd, alasum, dcopy, derrtr, dget04,
211  $ dtptrs
212 * ..
213 * .. Scalars in Common ..
214  LOGICAL lerr, ok
215  CHARACTER*32 srnamt
216  INTEGER infot, iounit
217 * ..
218 * .. Common blocks ..
219  COMMON / infoc / infot, iounit, ok, lerr
220  COMMON / srnamc / srnamt
221 * ..
222 * .. Intrinsic Functions ..
223  INTRINSIC max
224 * ..
225 * .. Data statements ..
226  DATA iseedy / 1988, 1989, 1990, 1991 /
227  DATA uplos / 'U', 'L' / , transs / 'N', 'T', 'C' /
228 * ..
229 * .. Executable Statements ..
230 *
231 * Initialize constants and the random number seed.
232 *
233  path( 1: 1 ) = 'Double precision'
234  path( 2: 3 ) = 'TP'
235  nrun = 0
236  nfail = 0
237  nerrs = 0
238  DO 10 i = 1, 4
239  iseed( i ) = iseedy( i )
240  10 CONTINUE
241 *
242 * Test the error exits
243 *
244  IF( tsterr )
245  $ CALL derrtr( path, nout )
246  infot = 0
247 *
248  DO 110 in = 1, nn
249 *
250 * Do for each value of N in NVAL
251 *
252  n = nval( in )
253  lda = max( 1, n )
254  lap = lda*( lda+1 ) / 2
255  xtype = 'N'
256 *
257  DO 70 imat = 1, ntype1
258 *
259 * Do the tests only if DOTYPE( IMAT ) is true.
260 *
261  IF( .NOT.dotype( imat ) )
262  $ GO TO 70
263 *
264  DO 60 iuplo = 1, 2
265 *
266 * Do first for UPLO = 'U', then for UPLO = 'L'
267 *
268  uplo = uplos( iuplo )
269 *
270 * Call DLATTP to generate a triangular test matrix.
271 *
272  srnamt = 'DLATTP'
273  CALL dlattp( imat, uplo, 'No transpose', diag, iseed, n,
274  $ ap, x, work, info )
275 *
276 * Set IDIAG = 1 for non-unit matrices, 2 for unit.
277 *
278  IF( lsame( diag, 'N' ) ) THEN
279  idiag = 1
280  ELSE
281  idiag = 2
282  END IF
283 *
284 *+ TEST 1
285 * Form the inverse of A.
286 *
287  IF( n.GT.0 )
288  $ CALL dcopy( lap, ap, 1, ainvp, 1 )
289  srnamt = 'DTPTRI'
290  CALL dtptri( uplo, diag, n, ainvp, info )
291 *
292 * Check error code from DTPTRI.
293 *
294  IF( info.NE.0 )
295  $ CALL alaerh( path, 'DTPTRI', info, 0, uplo // diag, n,
296  $ n, -1, -1, -1, imat, nfail, nerrs, nout )
297 *
298 * Compute the infinity-norm condition number of A.
299 *
300  anorm = dlantp( 'I', uplo, diag, n, ap, rwork )
301  ainvnm = dlantp( 'I', uplo, diag, n, ainvp, rwork )
302  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
303  rcondi = one
304  ELSE
305  rcondi = ( one / anorm ) / ainvnm
306  END IF
307 *
308 * Compute the residual for the triangular matrix times its
309 * inverse. Also compute the 1-norm condition number of A.
310 *
311  CALL dtpt01( uplo, diag, n, ap, ainvp, rcondo, rwork,
312  $ result( 1 ) )
313 *
314 * Print the test ratio if it is .GE. THRESH.
315 *
316  IF( result( 1 ).GE.thresh ) THEN
317  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
318  $ CALL alahd( nout, path )
319  WRITE( nout, fmt = 9999 )uplo, diag, n, imat, 1,
320  $ result( 1 )
321  nfail = nfail + 1
322  END IF
323  nrun = nrun + 1
324 *
325  DO 40 irhs = 1, nns
326  nrhs = nsval( irhs )
327  xtype = 'N'
328 *
329  DO 30 itran = 1, ntran
330 *
331 * Do for op(A) = A, A**T, or A**H.
332 *
333  trans = transs( itran )
334  IF( itran.EQ.1 ) THEN
335  norm = 'O'
336  rcondc = rcondo
337  ELSE
338  norm = 'I'
339  rcondc = rcondi
340  END IF
341 *
342 *+ TEST 2
343 * Solve and compute residual for op(A)*x = b.
344 *
345  srnamt = 'DLARHS'
346  CALL dlarhs( path, xtype, uplo, trans, n, n, 0,
347  $ idiag, nrhs, ap, lap, xact, lda, b,
348  $ lda, iseed, info )
349  xtype = 'C'
350  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
351 *
352  srnamt = 'DTPTRS'
353  CALL dtptrs( uplo, trans, diag, n, nrhs, ap, x,
354  $ lda, info )
355 *
356 * Check error code from DTPTRS.
357 *
358  IF( info.NE.0 )
359  $ CALL alaerh( path, 'DTPTRS', info, 0,
360  $ uplo // trans // diag, n, n, -1,
361  $ -1, -1, imat, nfail, nerrs, nout )
362 *
363  CALL dtpt02( uplo, trans, diag, n, nrhs, ap, x,
364  $ lda, b, lda, work, result( 2 ) )
365 *
366 *+ TEST 3
367 * Check solution from generated exact solution.
368 *
369  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
370  $ result( 3 ) )
371 *
372 *+ TESTS 4, 5, and 6
373 * Use iterative refinement to improve the solution and
374 * compute error bounds.
375 *
376  srnamt = 'DTPRFS'
377  CALL dtprfs( uplo, trans, diag, n, nrhs, ap, b,
378  $ lda, x, lda, rwork, rwork( nrhs+1 ),
379  $ work, iwork, info )
380 *
381 * Check error code from DTPRFS.
382 *
383  IF( info.NE.0 )
384  $ CALL alaerh( path, 'DTPRFS', info, 0,
385  $ uplo // trans // diag, n, n, -1,
386  $ -1, nrhs, imat, nfail, nerrs,
387  $ nout )
388 *
389  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
390  $ result( 4 ) )
391  CALL dtpt05( uplo, trans, diag, n, nrhs, ap, b,
392  $ lda, x, lda, xact, lda, rwork,
393  $ rwork( nrhs+1 ), result( 5 ) )
394 *
395 * Print information about the tests that did not pass
396 * the threshold.
397 *
398  DO 20 k = 2, 6
399  IF( result( k ).GE.thresh ) THEN
400  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
401  $ CALL alahd( nout, path )
402  WRITE( nout, fmt = 9998 )uplo, trans, diag,
403  $ n, nrhs, imat, k, result( k )
404  nfail = nfail + 1
405  END IF
406  20 CONTINUE
407  nrun = nrun + 5
408  30 CONTINUE
409  40 CONTINUE
410 *
411 *+ TEST 7
412 * Get an estimate of RCOND = 1/CNDNUM.
413 *
414  DO 50 itran = 1, 2
415  IF( itran.EQ.1 ) THEN
416  norm = 'O'
417  rcondc = rcondo
418  ELSE
419  norm = 'I'
420  rcondc = rcondi
421  END IF
422 *
423  srnamt = 'DTPCON'
424  CALL dtpcon( norm, uplo, diag, n, ap, rcond, work,
425  $ iwork, info )
426 *
427 * Check error code from DTPCON.
428 *
429  IF( info.NE.0 )
430  $ CALL alaerh( path, 'DTPCON', info, 0,
431  $ norm // uplo // diag, n, n, -1, -1,
432  $ -1, imat, nfail, nerrs, nout )
433 *
434  CALL dtpt06( rcond, rcondc, uplo, diag, n, ap, rwork,
435  $ result( 7 ) )
436 *
437 * Print the test ratio if it is .GE. THRESH.
438 *
439  IF( result( 7 ).GE.thresh ) THEN
440  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
441  $ CALL alahd( nout, path )
442  WRITE( nout, fmt = 9997 ) 'DTPCON', norm, uplo,
443  $ diag, n, imat, 7, result( 7 )
444  nfail = nfail + 1
445  END IF
446  nrun = nrun + 1
447  50 CONTINUE
448  60 CONTINUE
449  70 CONTINUE
450 *
451 * Use pathological test matrices to test DLATPS.
452 *
453  DO 100 imat = ntype1 + 1, ntypes
454 *
455 * Do the tests only if DOTYPE( IMAT ) is true.
456 *
457  IF( .NOT.dotype( imat ) )
458  $ GO TO 100
459 *
460  DO 90 iuplo = 1, 2
461 *
462 * Do first for UPLO = 'U', then for UPLO = 'L'
463 *
464  uplo = uplos( iuplo )
465  DO 80 itran = 1, ntran
466 *
467 * Do for op(A) = A, A**T, or A**H.
468 *
469  trans = transs( itran )
470 *
471 * Call DLATTP to generate a triangular test matrix.
472 *
473  srnamt = 'DLATTP'
474  CALL dlattp( imat, uplo, trans, diag, iseed, n, ap, x,
475  $ work, info )
476 *
477 *+ TEST 8
478 * Solve the system op(A)*x = b.
479 *
480  srnamt = 'DLATPS'
481  CALL dcopy( n, x, 1, b, 1 )
482  CALL dlatps( uplo, trans, diag, 'N', n, ap, b, scale,
483  $ rwork, info )
484 *
485 * Check error code from DLATPS.
486 *
487  IF( info.NE.0 )
488  $ CALL alaerh( path, 'DLATPS', info, 0,
489  $ uplo // trans // diag // 'N', n, n,
490  $ -1, -1, -1, imat, nfail, nerrs, nout )
491 *
492  CALL dtpt03( uplo, trans, diag, n, 1, ap, scale,
493  $ rwork, one, b, lda, x, lda, work,
494  $ result( 8 ) )
495 *
496 *+ TEST 9
497 * Solve op(A)*x = b again with NORMIN = 'Y'.
498 *
499  CALL dcopy( n, x, 1, b( n+1 ), 1 )
500  CALL dlatps( uplo, trans, diag, 'Y', n, ap, b( n+1 ),
501  $ scale, rwork, info )
502 *
503 * Check error code from DLATPS.
504 *
505  IF( info.NE.0 )
506  $ CALL alaerh( path, 'DLATPS', info, 0,
507  $ uplo // trans // diag // 'Y', n, n,
508  $ -1, -1, -1, imat, nfail, nerrs, nout )
509 *
510  CALL dtpt03( uplo, trans, diag, n, 1, ap, scale,
511  $ rwork, one, b( n+1 ), lda, x, lda, work,
512  $ result( 9 ) )
513 *
514 * Print information about the tests that did not pass
515 * the threshold.
516 *
517  IF( result( 8 ).GE.thresh ) THEN
518  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
519  $ CALL alahd( nout, path )
520  WRITE( nout, fmt = 9996 )'DLATPS', uplo, trans,
521  $ diag, 'N', n, imat, 8, result( 8 )
522  nfail = nfail + 1
523  END IF
524  IF( result( 9 ).GE.thresh ) THEN
525  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
526  $ CALL alahd( nout, path )
527  WRITE( nout, fmt = 9996 )'DLATPS', uplo, trans,
528  $ diag, 'Y', n, imat, 9, result( 9 )
529  nfail = nfail + 1
530  END IF
531  nrun = nrun + 2
532  80 CONTINUE
533  90 CONTINUE
534  100 CONTINUE
535  110 CONTINUE
536 *
537 * Print a summary of the results.
538 *
539  CALL alasum( path, nout, nfail, nrun, nerrs )
540 *
541  9999 FORMAT( ' UPLO=''', a1, ''', DIAG=''', a1, ''', N=', i5,
542  $ ', type ', i2, ', test(', i2, ')= ', g12.5 )
543  9998 FORMAT( ' UPLO=''', a1, ''', TRANS=''', a1, ''', DIAG=''', a1,
544  $ ''', N=', i5, ''', NRHS=', i5, ', type ', i2, ', test(',
545  $ i2, ')= ', g12.5 )
546  9997 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ''', a1, ''',',
547  $ i5, ', ... ), type ', i2, ', test(', i2, ')=', g12.5 )
548  9996 FORMAT( 1x, a, '( ''', a1, ''', ''', a1, ''', ''', a1, ''', ''',
549  $ a1, ''',', i5, ', ... ), type ', i2, ', test(', i2, ')=',
550  $ g12.5 )
551  RETURN
552 *
553 * End of DCHKTP
554 *
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 dlattp(IMAT, UPLO, TRANS, DIAG, ISEED, N, A, B, WORK, INFO)
DLATTP
Definition: dlattp.f:127
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dtpt06(RCOND, RCONDC, UPLO, DIAG, N, AP, WORK, RAT)
DTPT06
Definition: dtpt06.f:113
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
double precision function dlantp(NORM, UPLO, DIAG, N, AP, WORK)
DLANTP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a triangular matrix supplied in packed form.
Definition: dlantp.f:126
subroutine dtpt02(UPLO, TRANS, DIAG, N, NRHS, AP, X, LDX, B, LDB, WORK, RESID)
DTPT02
Definition: dtpt02.f:143
subroutine dtpt05(UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DTPT05
Definition: dtpt05.f:176
subroutine dtptrs(UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, INFO)
DTPTRS
Definition: dtptrs.f:132
subroutine dtpt01(UPLO, DIAG, N, AP, AINVP, RCOND, WORK, RESID)
DTPT01
Definition: dtpt01.f:110
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine derrtr(PATH, NUNIT)
DERRTR
Definition: derrtr.f:57
subroutine dlatps(UPLO, TRANS, DIAG, NORMIN, N, AP, X, SCALE, CNORM, INFO)
DLATPS solves a triangular system of equations with the matrix held in packed storage.
Definition: dlatps.f:231
subroutine dtpt03(UPLO, TRANS, DIAG, N, NRHS, AP, SCALE, CNORM, TSCAL, X, LDX, B, LDB, WORK, RESID)
DTPT03
Definition: dtpt03.f:163
subroutine dtptri(UPLO, DIAG, N, AP, INFO)
DTPTRI
Definition: dtptri.f:119
subroutine dtprfs(UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DTPRFS
Definition: dtprfs.f:177
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dtpcon(NORM, UPLO, DIAG, N, AP, RCOND, WORK, IWORK, INFO)
DTPCON
Definition: dtpcon.f:132
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: