LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dchkpt ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
double precision  THRESH,
logical  TSTERR,
double precision, dimension( * )  A,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
double precision, dimension( * )  B,
double precision, dimension( * )  X,
double precision, dimension( * )  XACT,
double precision, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer  NOUT 
)

DCHKPT

Purpose:
 DCHKPT tests DPTTRF, -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.
[out]A
          A is DOUBLE PRECISION array, dimension (NMAX*2)
[out]D
          D is DOUBLE PRECISION array, dimension (NMAX*2)
[out]E
          E is DOUBLE PRECISION array, dimension (NMAX*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]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 148 of file dchkpt.f.

148 *
149 * -- LAPACK test routine (version 3.4.0) --
150 * -- LAPACK is a software package provided by Univ. of Tennessee, --
151 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
152 * November 2011
153 *
154 * .. Scalar Arguments ..
155  LOGICAL tsterr
156  INTEGER nn, nns, nout
157  DOUBLE PRECISION thresh
158 * ..
159 * .. Array Arguments ..
160  LOGICAL dotype( * )
161  INTEGER nsval( * ), nval( * )
162  DOUBLE PRECISION a( * ), b( * ), d( * ), e( * ), rwork( * ),
163  $ work( * ), x( * ), xact( * )
164 * ..
165 *
166 * =====================================================================
167 *
168 * .. Parameters ..
169  DOUBLE PRECISION one, zero
170  parameter ( one = 1.0d+0, zero = 0.0d+0 )
171  INTEGER ntypes
172  parameter ( ntypes = 12 )
173  INTEGER ntests
174  parameter ( ntests = 7 )
175 * ..
176 * .. Local Scalars ..
177  LOGICAL zerot
178  CHARACTER dist, type
179  CHARACTER*3 path
180  INTEGER i, ia, imat, in, info, irhs, ix, izero, j, k,
181  $ kl, ku, lda, mode, n, nerrs, nfail, nimat,
182  $ nrhs, nrun
183  DOUBLE PRECISION ainvnm, anorm, cond, dmax, rcond, rcondc
184 * ..
185 * .. Local Arrays ..
186  INTEGER iseed( 4 ), iseedy( 4 )
187  DOUBLE PRECISION result( ntests ), z( 3 )
188 * ..
189 * .. External Functions ..
190  INTEGER idamax
191  DOUBLE PRECISION dasum, dget06, dlanst
192  EXTERNAL idamax, dasum, dget06, dlanst
193 * ..
194 * .. External Subroutines ..
195  EXTERNAL alaerh, alahd, alasum, dcopy, derrgt, dget04,
198  $ dscal
199 * ..
200 * .. Intrinsic Functions ..
201  INTRINSIC abs, max
202 * ..
203 * .. Scalars in Common ..
204  LOGICAL lerr, ok
205  CHARACTER*32 srnamt
206  INTEGER infot, nunit
207 * ..
208 * .. Common blocks ..
209  COMMON / infoc / infot, nunit, ok, lerr
210  COMMON / srnamc / srnamt
211 * ..
212 * .. Data statements ..
213  DATA iseedy / 0, 0, 0, 1 /
214 * ..
215 * .. Executable Statements ..
216 *
217  path( 1: 1 ) = 'Double precision'
218  path( 2: 3 ) = 'PT'
219  nrun = 0
220  nfail = 0
221  nerrs = 0
222  DO 10 i = 1, 4
223  iseed( i ) = iseedy( i )
224  10 CONTINUE
225 *
226 * Test the error exits
227 *
228  IF( tsterr )
229  $ CALL derrgt( path, nout )
230  infot = 0
231 *
232  DO 110 in = 1, nn
233 *
234 * Do for each value of N in NVAL.
235 *
236  n = nval( in )
237  lda = max( 1, n )
238  nimat = ntypes
239  IF( n.LE.0 )
240  $ nimat = 1
241 *
242  DO 100 imat = 1, nimat
243 *
244 * Do the tests only if DOTYPE( IMAT ) is true.
245 *
246  IF( n.GT.0 .AND. .NOT.dotype( imat ) )
247  $ GO TO 100
248 *
249 * Set up parameters with DLATB4.
250 *
251  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
252  $ cond, dist )
253 *
254  zerot = imat.GE.8 .AND. imat.LE.10
255  IF( imat.LE.6 ) THEN
256 *
257 * Type 1-6: generate a symmetric tridiagonal matrix of
258 * known condition number in lower triangular band storage.
259 *
260  srnamt = 'DLATMS'
261  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode, cond,
262  $ anorm, kl, ku, 'B', a, 2, work, info )
263 *
264 * Check the error code from DLATMS.
265 *
266  IF( info.NE.0 ) THEN
267  CALL alaerh( path, 'DLATMS', info, 0, ' ', n, n, kl,
268  $ ku, -1, imat, nfail, nerrs, nout )
269  GO TO 100
270  END IF
271  izero = 0
272 *
273 * Copy the matrix to D and E.
274 *
275  ia = 1
276  DO 20 i = 1, n - 1
277  d( i ) = a( ia )
278  e( i ) = a( ia+1 )
279  ia = ia + 2
280  20 CONTINUE
281  IF( n.GT.0 )
282  $ d( n ) = a( ia )
283  ELSE
284 *
285 * Type 7-12: generate a diagonally dominant matrix with
286 * unknown condition number in the vectors D and E.
287 *
288  IF( .NOT.zerot .OR. .NOT.dotype( 7 ) ) THEN
289 *
290 * Let D and E have values from [-1,1].
291 *
292  CALL dlarnv( 2, iseed, n, d )
293  CALL dlarnv( 2, iseed, n-1, e )
294 *
295 * Make the tridiagonal matrix diagonally dominant.
296 *
297  IF( n.EQ.1 ) THEN
298  d( 1 ) = abs( d( 1 ) )
299  ELSE
300  d( 1 ) = abs( d( 1 ) ) + abs( e( 1 ) )
301  d( n ) = abs( d( n ) ) + abs( e( n-1 ) )
302  DO 30 i = 2, n - 1
303  d( i ) = abs( d( i ) ) + abs( e( i ) ) +
304  $ abs( e( i-1 ) )
305  30 CONTINUE
306  END IF
307 *
308 * Scale D and E so the maximum element is ANORM.
309 *
310  ix = idamax( n, d, 1 )
311  dmax = d( ix )
312  CALL dscal( n, anorm / dmax, d, 1 )
313  CALL dscal( n-1, anorm / dmax, e, 1 )
314 *
315  ELSE IF( izero.GT.0 ) THEN
316 *
317 * Reuse the last matrix by copying back the zeroed out
318 * elements.
319 *
320  IF( izero.EQ.1 ) THEN
321  d( 1 ) = z( 2 )
322  IF( n.GT.1 )
323  $ e( 1 ) = z( 3 )
324  ELSE IF( izero.EQ.n ) THEN
325  e( n-1 ) = z( 1 )
326  d( n ) = z( 2 )
327  ELSE
328  e( izero-1 ) = z( 1 )
329  d( izero ) = z( 2 )
330  e( izero ) = z( 3 )
331  END IF
332  END IF
333 *
334 * For types 8-10, set one row and column of the matrix to
335 * zero.
336 *
337  izero = 0
338  IF( imat.EQ.8 ) THEN
339  izero = 1
340  z( 2 ) = d( 1 )
341  d( 1 ) = zero
342  IF( n.GT.1 ) THEN
343  z( 3 ) = e( 1 )
344  e( 1 ) = zero
345  END IF
346  ELSE IF( imat.EQ.9 ) THEN
347  izero = n
348  IF( n.GT.1 ) THEN
349  z( 1 ) = e( n-1 )
350  e( n-1 ) = zero
351  END IF
352  z( 2 ) = d( n )
353  d( n ) = zero
354  ELSE IF( imat.EQ.10 ) THEN
355  izero = ( n+1 ) / 2
356  IF( izero.GT.1 ) THEN
357  z( 1 ) = e( izero-1 )
358  e( izero-1 ) = zero
359  z( 3 ) = e( izero )
360  e( izero ) = zero
361  END IF
362  z( 2 ) = d( izero )
363  d( izero ) = zero
364  END IF
365  END IF
366 *
367  CALL dcopy( n, d, 1, d( n+1 ), 1 )
368  IF( n.GT.1 )
369  $ CALL dcopy( n-1, e, 1, e( n+1 ), 1 )
370 *
371 *+ TEST 1
372 * Factor A as L*D*L' and compute the ratio
373 * norm(L*D*L' - A) / (n * norm(A) * EPS )
374 *
375  CALL dpttrf( n, d( n+1 ), e( n+1 ), info )
376 *
377 * Check error code from DPTTRF.
378 *
379  IF( info.NE.izero ) THEN
380  CALL alaerh( path, 'DPTTRF', info, izero, ' ', n, n, -1,
381  $ -1, -1, imat, nfail, nerrs, nout )
382  GO TO 100
383  END IF
384 *
385  IF( info.GT.0 ) THEN
386  rcondc = zero
387  GO TO 90
388  END IF
389 *
390  CALL dptt01( n, d, e, d( n+1 ), e( n+1 ), work,
391  $ result( 1 ) )
392 *
393 * Print the test ratio if greater than or equal to THRESH.
394 *
395  IF( result( 1 ).GE.thresh ) THEN
396  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
397  $ CALL alahd( nout, path )
398  WRITE( nout, fmt = 9999 )n, imat, 1, result( 1 )
399  nfail = nfail + 1
400  END IF
401  nrun = nrun + 1
402 *
403 * Compute RCONDC = 1 / (norm(A) * norm(inv(A))
404 *
405 * Compute norm(A).
406 *
407  anorm = dlanst( '1', n, d, e )
408 *
409 * Use DPTTRS to solve for one column at a time of inv(A),
410 * computing the maximum column sum as we go.
411 *
412  ainvnm = zero
413  DO 50 i = 1, n
414  DO 40 j = 1, n
415  x( j ) = zero
416  40 CONTINUE
417  x( i ) = one
418  CALL dpttrs( n, 1, d( n+1 ), e( n+1 ), x, lda, info )
419  ainvnm = max( ainvnm, dasum( n, x, 1 ) )
420  50 CONTINUE
421  rcondc = one / max( one, anorm*ainvnm )
422 *
423  DO 80 irhs = 1, nns
424  nrhs = nsval( irhs )
425 *
426 * Generate NRHS random solution vectors.
427 *
428  ix = 1
429  DO 60 j = 1, nrhs
430  CALL dlarnv( 2, iseed, n, xact( ix ) )
431  ix = ix + lda
432  60 CONTINUE
433 *
434 * Set the right hand side.
435 *
436  CALL dlaptm( n, nrhs, one, d, e, xact, lda, zero, b,
437  $ lda )
438 *
439 *+ TEST 2
440 * Solve A*x = b and compute the residual.
441 *
442  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
443  CALL dpttrs( n, nrhs, d( n+1 ), e( n+1 ), x, lda, info )
444 *
445 * Check error code from DPTTRS.
446 *
447  IF( info.NE.0 )
448  $ CALL alaerh( path, 'DPTTRS', info, 0, ' ', n, n, -1,
449  $ -1, nrhs, imat, nfail, nerrs, nout )
450 *
451  CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
452  CALL dptt02( n, nrhs, d, e, x, lda, work, lda,
453  $ result( 2 ) )
454 *
455 *+ TEST 3
456 * Check solution from generated exact solution.
457 *
458  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
459  $ result( 3 ) )
460 *
461 *+ TESTS 4, 5, and 6
462 * Use iterative refinement to improve the solution.
463 *
464  srnamt = 'DPTRFS'
465  CALL dptrfs( n, nrhs, d, e, d( n+1 ), e( n+1 ), b, lda,
466  $ x, lda, rwork, rwork( nrhs+1 ), work, info )
467 *
468 * Check error code from DPTRFS.
469 *
470  IF( info.NE.0 )
471  $ CALL alaerh( path, 'DPTRFS', info, 0, ' ', n, n, -1,
472  $ -1, nrhs, imat, nfail, nerrs, nout )
473 *
474  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
475  $ result( 4 ) )
476  CALL dptt05( n, nrhs, d, e, b, lda, x, lda, xact, lda,
477  $ rwork, rwork( nrhs+1 ), result( 5 ) )
478 *
479 * Print information about the tests that did not pass the
480 * threshold.
481 *
482  DO 70 k = 2, 6
483  IF( result( k ).GE.thresh ) THEN
484  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
485  $ CALL alahd( nout, path )
486  WRITE( nout, fmt = 9998 )n, nrhs, imat, k,
487  $ result( k )
488  nfail = nfail + 1
489  END IF
490  70 CONTINUE
491  nrun = nrun + 5
492  80 CONTINUE
493 *
494 *+ TEST 7
495 * Estimate the reciprocal of the condition number of the
496 * matrix.
497 *
498  90 CONTINUE
499  srnamt = 'DPTCON'
500  CALL dptcon( n, d( n+1 ), e( n+1 ), anorm, rcond, rwork,
501  $ info )
502 *
503 * Check error code from DPTCON.
504 *
505  IF( info.NE.0 )
506  $ CALL alaerh( path, 'DPTCON', info, 0, ' ', n, n, -1, -1,
507  $ -1, imat, nfail, nerrs, nout )
508 *
509  result( 7 ) = dget06( rcond, rcondc )
510 *
511 * Print the test ratio if greater than or equal to THRESH.
512 *
513  IF( result( 7 ).GE.thresh ) THEN
514  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
515  $ CALL alahd( nout, path )
516  WRITE( nout, fmt = 9999 )n, imat, 7, result( 7 )
517  nfail = nfail + 1
518  END IF
519  nrun = nrun + 1
520  100 CONTINUE
521  110 CONTINUE
522 *
523 * Print a summary of the results.
524 *
525  CALL alasum( path, nout, nfail, nrun, nerrs )
526 *
527  9999 FORMAT( ' N =', i5, ', type ', i2, ', test ', i2, ', ratio = ',
528  $ g12.5 )
529  9998 FORMAT( ' N =', i5, ', NRHS=', i3, ', type ', i2, ', test(', i2,
530  $ ') = ', g12.5 )
531  RETURN
532 *
533 * End of DCHKPT
534 *
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 dptt01(N, D, E, DF, EF, WORK, RESID)
DPTT01
Definition: dptt01.f:93
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
subroutine dptt02(N, NRHS, D, E, X, LDX, B, LDB, RESID)
DPTT02
Definition: dptt02.f:106
subroutine dpttrf(N, D, E, INFO)
DPTTRF
Definition: dpttrf.f:93
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 dlaptm(N, NRHS, ALPHA, D, E, X, LDX, BETA, B, LDB)
DLAPTM
Definition: dlaptm.f:118
subroutine dptcon(N, D, E, ANORM, RCOND, WORK, INFO)
DPTCON
Definition: dptcon.f:120
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 dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
double precision function dasum(N, DX, INCX)
DASUM
Definition: dasum.f:53
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine derrgt(PATH, NUNIT)
DERRGT
Definition: derrgt.f:57
subroutine dptrfs(N, NRHS, D, E, DF, EF, B, LDB, X, LDX, FERR, BERR, WORK, INFO)
DPTRFS
Definition: dptrfs.f:165
double precision function dlanst(NORM, N, D, E)
DLANST returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric tridiagonal matrix.
Definition: dlanst.f:102
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dlarnv(IDIST, ISEED, N, X)
DLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: dlarnv.f:99
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75
subroutine dptt05(N, NRHS, D, E, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DPTT05
Definition: dptt05.f:152
subroutine dpttrs(N, NRHS, D, E, B, LDB, INFO)
DPTTRS
Definition: dpttrs.f:111

Here is the call graph for this function:

Here is the caller graph for this function: