LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine dchksy_aa ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
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 
)

DCHKSY_AA

Purpose:
 DCHKSY_AA tests DSYTRF_AA, -TRS_AA.
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 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)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AINV
          AINV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[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))
[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
December 2016

Definition at line 174 of file dchksy_aa.f.

174 *
175 * -- LAPACK test routine (version 3.7.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 * December 2016
179 *
180  IMPLICIT NONE
181 *
182 * .. Scalar Arguments ..
183  LOGICAL tsterr
184  INTEGER nn, nnb, nns, nmax, nout
185  DOUBLE PRECISION thresh
186 * ..
187 * .. Array Arguments ..
188  LOGICAL dotype( * )
189  INTEGER iwork( * ), nbval( * ), nsval( * ), nval( * )
190  DOUBLE PRECISION a( * ), afac( * ), ainv( * ), b( * ),
191  $ rwork( * ), work( * ), x( * ), xact( * )
192 * ..
193 *
194 * =====================================================================
195 *
196 * .. Parameters ..
197  DOUBLE PRECISION zero, one
198  parameter ( zero = 0.0d+0, one = 1.0d+0 )
199  INTEGER ntypes
200  parameter ( ntypes = 10 )
201  INTEGER ntests
202  parameter ( ntests = 9 )
203 * ..
204 * .. Local Scalars ..
205  LOGICAL zerot
206  CHARACTER dist, TYPE, uplo, xtype
207  CHARACTER*3 path, matpath
208  INTEGER i, i1, i2, imat, in, inb, info, ioff, irhs,
209  $ iuplo, izero, j, k, kl, ku, lda, lwork, mode,
210  $ n, nb, nerrs, nfail, nimat, nrhs, nrun, nt
211  DOUBLE PRECISION anorm, cndnum
212 * ..
213 * .. Local Arrays ..
214  CHARACTER uplos( 2 )
215  INTEGER iseed( 4 ), iseedy( 4 )
216  DOUBLE PRECISION result( ntests )
217 * ..
218 * .. External Functions ..
219  DOUBLE PRECISION dget06, dlansy
220  EXTERNAL dget06, dlansy
221 * ..
222 * .. External Subroutines ..
223  EXTERNAL alaerh, alahd, alasum, derrsy, dget04, dlacpy,
227 * ..
228 * .. Intrinsic Functions ..
229  INTRINSIC max, min
230 * ..
231 * .. Scalars in Common ..
232  LOGICAL lerr, ok
233  CHARACTER*32 srnamt
234  INTEGER infot, nunit
235 * ..
236 * .. Common blocks ..
237  COMMON / infoc / infot, nunit, ok, lerr
238  COMMON / srnamc / srnamt
239 * ..
240 * .. Data statements ..
241  DATA iseedy / 1988, 1989, 1990, 1991 /
242  DATA uplos / 'U', 'L' /
243 * ..
244 * .. Executable Statements ..
245 *
246 * Initialize constants and the random number seed.
247 *
248 * Test path
249 *
250  path( 1: 1 ) = 'Double precision'
251  path( 2: 3 ) = 'SA'
252 *
253 * Path to generate matrices
254 *
255  matpath( 1: 1 ) = 'Double precision'
256  matpath( 2: 3 ) = 'SY'
257  nrun = 0
258  nfail = 0
259  nerrs = 0
260  DO 10 i = 1, 4
261  iseed( i ) = iseedy( i )
262  10 CONTINUE
263 *
264 * Test the error exits
265 *
266  IF( tsterr )
267  $ CALL derrsy( path, nout )
268  infot = 0
269 *
270 * Set the minimum block size for which the block routine should
271 * be used, which will be later returned by ILAENV
272 *
273  CALL xlaenv( 2, 2 )
274 *
275 * Do for each value of N in NVAL
276 *
277  DO 180 in = 1, nn
278  n = nval( in )
279  IF( n .GT. nmax ) THEN
280  nfail = nfail + 1
281  WRITE(nout, 9995) 'M ', n, nmax
282  GO TO 180
283  END IF
284  lda = max( n, 1 )
285  xtype = 'N'
286  nimat = ntypes
287  IF( n.LE.0 )
288  $ nimat = 1
289 *
290  izero = 0
291 *
292 * Do for each value of matrix type IMAT
293 *
294  DO 170 imat = 1, nimat
295 *
296 * Do the tests only if DOTYPE( IMAT ) is true.
297 *
298  IF( .NOT.dotype( imat ) )
299  $ GO TO 170
300 *
301 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
302 *
303  zerot = imat.GE.3 .AND. imat.LE.6
304  IF( zerot .AND. n.LT.imat-2 )
305  $ GO TO 170
306 *
307 * Do first for UPLO = 'U', then for UPLO = 'L'
308 *
309  DO 160 iuplo = 1, 2
310  uplo = uplos( iuplo )
311 *
312 * Begin generate the test matrix A.
313 *
314 *
315 * Set up parameters with DLATB4 for the matrix generator
316 * based on the type of matrix to be generated.
317 *
318  CALL dlatb4( matpath, imat, n, n, TYPE, kl, ku,
319  $ anorm, mode, cndnum, dist )
320 *
321 * Generate a matrix with DLATMS.
322 *
323  srnamt = 'DLATMS'
324  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
325  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
326  $ info )
327 *
328 * Check error code from DLATMS and handle error.
329 *
330  IF( info.NE.0 ) THEN
331  CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
332  $ -1, -1, imat, nfail, nerrs, nout )
333 *
334 * Skip all tests for this generated matrix
335 *
336  GO TO 160
337  END IF
338 *
339 * For matrix types 3-6, zero one or more rows and
340 * columns of the matrix to test that INFO is returned
341 * correctly.
342 *
343  IF( zerot ) THEN
344  IF( imat.EQ.3 ) THEN
345  izero = 1
346  ELSE IF( imat.EQ.4 ) THEN
347  izero = n
348  ELSE
349  izero = n / 2 + 1
350  END IF
351 *
352  IF( imat.LT.6 ) THEN
353 *
354 * Set row and column IZERO to zero.
355 *
356  IF( iuplo.EQ.1 ) THEN
357  ioff = ( izero-1 )*lda
358  DO 20 i = 1, izero - 1
359  a( ioff+i ) = zero
360  20 CONTINUE
361  ioff = ioff + izero
362  DO 30 i = izero, n
363  a( ioff ) = zero
364  ioff = ioff + lda
365  30 CONTINUE
366  ELSE
367  ioff = izero
368  DO 40 i = 1, izero - 1
369  a( ioff ) = zero
370  ioff = ioff + lda
371  40 CONTINUE
372  ioff = ioff - izero
373  DO 50 i = izero, n
374  a( ioff+i ) = zero
375  50 CONTINUE
376  END IF
377  ELSE
378  IF( iuplo.EQ.1 ) THEN
379 *
380 * Set the first IZERO rows and columns to zero.
381 *
382  ioff = 0
383  DO 70 j = 1, n
384  i2 = min( j, izero )
385  DO 60 i = 1, i2
386  a( ioff+i ) = zero
387  60 CONTINUE
388  ioff = ioff + lda
389  70 CONTINUE
390  izero = 1
391  ELSE
392 *
393 * Set the last IZERO rows and columns to zero.
394 *
395  ioff = 0
396  DO 90 j = 1, n
397  i1 = max( j, izero )
398  DO 80 i = i1, n
399  a( ioff+i ) = zero
400  80 CONTINUE
401  ioff = ioff + lda
402  90 CONTINUE
403  END IF
404  END IF
405  ELSE
406  izero = 0
407  END IF
408 *
409 * End generate the test matrix A.
410 *
411 * Do for each value of NB in NBVAL
412 *
413  DO 150 inb = 1, nnb
414 *
415 * Set the optimal blocksize, which will be later
416 * returned by ILAENV.
417 *
418  nb = nbval( inb )
419  CALL xlaenv( 1, nb )
420 *
421 * Copy the test matrix A into matrix AFAC which
422 * will be factorized in place. This is needed to
423 * preserve the test matrix A for subsequent tests.
424 *
425  CALL dlacpy( uplo, n, n, a, lda, afac, lda )
426 *
427 * Compute the L*D*L**T or U*D*U**T factorization of the
428 * matrix. IWORK stores details of the interchanges and
429 * the block structure of D. AINV is a work array for
430 * block factorization, LWORK is the length of AINV.
431 *
432  srnamt = 'DSYTRF_AA'
433  lwork = max( 1, n*nb + n )
434  CALL dsytrf_aa( uplo, n, afac, lda, iwork, ainv,
435  $ lwork, info )
436 *
437 * Adjust the expected value of INFO to account for
438 * pivoting.
439 *
440  IF( izero.GT.0 ) THEN
441  j = 1
442  k = izero
443  100 CONTINUE
444  IF( j.EQ.k ) THEN
445  k = iwork( j )
446  ELSE IF( iwork( j ).EQ.k ) THEN
447  k = j
448  END IF
449  IF( j.LT.k ) THEN
450  j = j + 1
451  GO TO 100
452  END IF
453  ELSE
454  k = 0
455  END IF
456 *
457 * Check error code from DSYTRF and handle error.
458 *
459  IF( info.NE.k ) THEN
460  CALL alaerh( path, 'DSYTRF_AA', info, k, uplo,
461  $ n, n, -1, -1, nb, imat, nfail, nerrs,
462  $ nout )
463  END IF
464 *
465 *+ TEST 1
466 * Reconstruct matrix from factors and compute residual.
467 *
468  CALL dsyt01_aa( uplo, n, a, lda, afac, lda, iwork,
469  $ ainv, lda, rwork, result( 1 ) )
470  nt = 1
471 *
472 *
473 * Print information about the tests that did not pass
474 * the threshold.
475 *
476  DO 110 k = 1, nt
477  IF( result( k ).GE.thresh ) THEN
478  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
479  $ CALL alahd( nout, path )
480  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
481  $ result( k )
482  nfail = nfail + 1
483  END IF
484  110 CONTINUE
485  nrun = nrun + nt
486 *
487 * Skip solver test if INFO is not 0.
488 *
489  IF( info.NE.0 ) THEN
490  GO TO 140
491  END IF
492 *
493 * Do for each value of NRHS in NSVAL.
494 *
495  DO 130 irhs = 1, nns
496  nrhs = nsval( irhs )
497 *
498 *+ TEST 2 (Using TRS)
499 * Solve and compute residual for A * X = B.
500 *
501 * Choose a set of NRHS random solution vectors
502 * stored in XACT and set up the right hand side B
503 *
504  srnamt = 'DLARHS'
505  CALL dlarhs( matpath, xtype, uplo, ' ', n, n,
506  $ kl, ku, nrhs, a, lda, xact, lda,
507  $ b, lda, iseed, info )
508  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
509 *
510  srnamt = 'DSYTRS_AA'
511  lwork = max( 1, 3*n-2 )
512  CALL dsytrs_aa( uplo, n, nrhs, afac, lda,
513  $ iwork, x, lda, work, lwork,
514  $ info )
515 *
516 * Check error code from DSYTRS and handle error.
517 *
518  IF( info.NE.0 ) THEN
519  CALL alaerh( path, 'DSYTRS_AA', info, 0,
520  $ uplo, n, n, -1, -1, nrhs, imat,
521  $ nfail, nerrs, nout )
522  END IF
523 *
524  CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
525 *
526 * Compute the residual for the solution
527 *
528  CALL dpot02( uplo, n, nrhs, a, lda, x, lda, work,
529  $ lda, rwork, result( 2 ) )
530 *
531 *
532 * Print information about the tests that did not pass
533 * the threshold.
534 *
535  DO 120 k = 2, 2
536  IF( result( k ).GE.thresh ) THEN
537  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
538  $ CALL alahd( nout, path )
539  WRITE( nout, fmt = 9998 )uplo, n, nrhs,
540  $ imat, k, result( k )
541  nfail = nfail + 1
542  END IF
543  120 CONTINUE
544  nrun = nrun + 1
545 *
546 * End do for each value of NRHS in NSVAL.
547 *
548  130 CONTINUE
549  140 CONTINUE
550  150 CONTINUE
551  160 CONTINUE
552  170 CONTINUE
553  180 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, ', NB =', i4, ', type ',
560  $ i2, ', test ', i2, ', ratio =', g12.5 )
561  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
562  $ i2, ', test(', i2, ') =', g12.5 )
563  9995 FORMAT( ' Invalid input value: ', a4, '=', i6, '; must be <=',
564  $ i6 )
565  RETURN
566 *
567 * End of DCHKSY_AA
568 *
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 alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:107
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY 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 matrix.
Definition: dlansy.f:124
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 dsyrfs(UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DSYRFS
Definition: dsyrfs.f:193
subroutine dsytri2(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
DSYTRI2
Definition: dsytri2.f:129
subroutine dsycon(UPLO, N, A, LDA, IPIV, ANORM, RCOND, WORK, IWORK, INFO)
DSYCON
Definition: dsycon.f:132
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine dsyt01_aa(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
DSYT01
Definition: dsyt01_aa.f:128
subroutine dpot03(UPLO, N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
DPOT03
Definition: dpot03.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 dsytrf_aa(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
DSYTRF_AA
Definition: dsytrf_aa.f:138
subroutine derrsy(PATH, NUNIT)
DERRSY
Definition: derrsy.f:57
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dpot05(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DPOT05
Definition: dpot05.f:166
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dsytrs_aa(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK, LWORK, INFO)
DSYTRS_AA
Definition: dsytrs_aa.f:131
subroutine dpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DPOT02
Definition: dpot02.f:129
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: