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

CCHKPT

Purpose:
 CCHKPT tests CPTTRF, -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 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.
[out]A
          A is COMPLEX array, dimension (NMAX*2)
[out]D
          D is REAL array, dimension (NMAX*2)
[out]E
          E is COMPLEX array, dimension (NMAX*2)
[out]B
          B is COMPLEX array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL 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 149 of file cchkpt.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: