LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ dchksy_aa()

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
November 2017

Definition at line 174 of file dchksy_aa.f.

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