LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dchksy()

subroutine dchksy ( 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

Purpose:
 DCHKSY tests DSYTRF, -TRI2, -TRS, -TRS2, -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]NNB
          NNB is INTEGER
          The number of values of NB contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          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.

Definition at line 167 of file dchksy.f.

170*
171* -- LAPACK test routine --
172* -- LAPACK is a software package provided by Univ. of Tennessee, --
173* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
174*
175* .. Scalar Arguments ..
176 LOGICAL TSTERR
177 INTEGER NMAX, NN, NNB, NNS, NOUT
178 DOUBLE PRECISION THRESH
179* ..
180* .. Array Arguments ..
181 LOGICAL DOTYPE( * )
182 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
183 DOUBLE PRECISION A( * ), AFAC( * ), AINV( * ), B( * ),
184 $ RWORK( * ), WORK( * ), X( * ), XACT( * )
185* ..
186*
187* =====================================================================
188*
189* .. Parameters ..
190 DOUBLE PRECISION ZERO
191 parameter( zero = 0.0d+0 )
192 INTEGER NTYPES
193 parameter( ntypes = 10 )
194 INTEGER NTESTS
195 parameter( ntests = 9 )
196* ..
197* .. Local Scalars ..
198 LOGICAL TRFCON, ZEROT
199 CHARACTER DIST, TYPE, UPLO, XTYPE
200 CHARACTER*3 PATH
201 INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
202 $ IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
203 $ N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
204 DOUBLE PRECISION ANORM, CNDNUM, RCOND, RCONDC
205* ..
206* .. Local Arrays ..
207 CHARACTER UPLOS( 2 )
208 INTEGER ISEED( 4 ), ISEEDY( 4 )
209 DOUBLE PRECISION RESULT( NTESTS )
210* ..
211* .. External Functions ..
212 DOUBLE PRECISION DGET06, DLANSY
213 EXTERNAL dget06, dlansy
214* ..
215* .. External Subroutines ..
216 EXTERNAL alaerh, alahd, alasum, derrsy, dget04, dlacpy,
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC max, min
223* ..
224* .. Scalars in Common ..
225 LOGICAL LERR, OK
226 CHARACTER*32 SRNAMT
227 INTEGER INFOT, NUNIT
228* ..
229* .. Common blocks ..
230 COMMON / infoc / infot, nunit, ok, lerr
231 COMMON / srnamc / srnamt
232* ..
233* .. Data statements ..
234 DATA iseedy / 1988, 1989, 1990, 1991 /
235 DATA uplos / 'U', 'L' /
236* ..
237* .. Executable Statements ..
238*
239* Initialize constants and the random number seed.
240*
241 path( 1: 1 ) = 'Double precision'
242 path( 2: 3 ) = 'SY'
243 nrun = 0
244 nfail = 0
245 nerrs = 0
246 DO 10 i = 1, 4
247 iseed( i ) = iseedy( i )
248 10 CONTINUE
249*
250* Test the error exits
251*
252 IF( tsterr )
253 $ CALL derrsy( path, nout )
254 infot = 0
255*
256* Set the minimum block size for which the block routine should
257* be used, which will be later returned by ILAENV
258*
259 CALL xlaenv( 2, 2 )
260*
261* Do for each value of N in NVAL
262*
263 DO 180 in = 1, nn
264 n = nval( in )
265 lda = max( n, 1 )
266 xtype = 'N'
267 nimat = ntypes
268 IF( n.LE.0 )
269 $ nimat = 1
270*
271 izero = 0
272*
273* Do for each value of matrix type IMAT
274*
275 DO 170 imat = 1, nimat
276*
277* Do the tests only if DOTYPE( IMAT ) is true.
278*
279 IF( .NOT.dotype( imat ) )
280 $ GO TO 170
281*
282* Skip types 3, 4, 5, or 6 if the matrix size is too small.
283*
284 zerot = imat.GE.3 .AND. imat.LE.6
285 IF( zerot .AND. n.LT.imat-2 )
286 $ GO TO 170
287*
288* Do first for UPLO = 'U', then for UPLO = 'L'
289*
290 DO 160 iuplo = 1, 2
291 uplo = uplos( iuplo )
292*
293* Begin generate the test matrix A.
294*
295*
296* Set up parameters with DLATB4 for the matrix generator
297* based on the type of matrix to be generated.
298*
299 CALL dlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
300 $ CNDNUM, DIST )
301*
302* Generate a matrix with DLATMS.
303*
304 srnamt = 'DLATMS'
305 CALL dlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
306 $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
307 $ INFO )
308*
309* Check error code from DLATMS and handle error.
310*
311 IF( info.NE.0 ) THEN
312 CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
313 $ -1, -1, imat, nfail, nerrs, nout )
314*
315* Skip all tests for this generated matrix
316*
317 GO TO 160
318 END IF
319*
320* For matrix types 3-6, zero one or more rows and
321* columns of the matrix to test that INFO is returned
322* correctly.
323*
324 IF( zerot ) THEN
325 IF( imat.EQ.3 ) THEN
326 izero = 1
327 ELSE IF( imat.EQ.4 ) THEN
328 izero = n
329 ELSE
330 izero = n / 2 + 1
331 END IF
332*
333 IF( imat.LT.6 ) THEN
334*
335* Set row and column IZERO to zero.
336*
337 IF( iuplo.EQ.1 ) THEN
338 ioff = ( izero-1 )*lda
339 DO 20 i = 1, izero - 1
340 a( ioff+i ) = zero
341 20 CONTINUE
342 ioff = ioff + izero
343 DO 30 i = izero, n
344 a( ioff ) = zero
345 ioff = ioff + lda
346 30 CONTINUE
347 ELSE
348 ioff = izero
349 DO 40 i = 1, izero - 1
350 a( ioff ) = zero
351 ioff = ioff + lda
352 40 CONTINUE
353 ioff = ioff - izero
354 DO 50 i = izero, n
355 a( ioff+i ) = zero
356 50 CONTINUE
357 END IF
358 ELSE
359 IF( iuplo.EQ.1 ) THEN
360*
361* Set the first IZERO rows and columns to zero.
362*
363 ioff = 0
364 DO 70 j = 1, n
365 i2 = min( j, izero )
366 DO 60 i = 1, i2
367 a( ioff+i ) = zero
368 60 CONTINUE
369 ioff = ioff + lda
370 70 CONTINUE
371 ELSE
372*
373* Set the last IZERO rows and columns to zero.
374*
375 ioff = 0
376 DO 90 j = 1, n
377 i1 = max( j, izero )
378 DO 80 i = i1, n
379 a( ioff+i ) = zero
380 80 CONTINUE
381 ioff = ioff + lda
382 90 CONTINUE
383 END IF
384 END IF
385 ELSE
386 izero = 0
387 END IF
388*
389* End generate the test matrix A.
390*
391* Do for each value of NB in NBVAL
392*
393 DO 150 inb = 1, nnb
394*
395* Set the optimal blocksize, which will be later
396* returned by ILAENV.
397*
398 nb = nbval( inb )
399 CALL xlaenv( 1, nb )
400*
401* Copy the test matrix A into matrix AFAC which
402* will be factorized in place. This is needed to
403* preserve the test matrix A for subsequent tests.
404*
405 CALL dlacpy( uplo, n, n, a, lda, afac, lda )
406*
407* Compute the L*D*L**T or U*D*U**T factorization of the
408* matrix. IWORK stores details of the interchanges and
409* the block structure of D. AINV is a work array for
410* block factorization, LWORK is the length of AINV.
411*
412 lwork = max( 2, nb )*lda
413 srnamt = 'DSYTRF'
414 CALL dsytrf( uplo, n, afac, lda, iwork, ainv, lwork,
415 $ info )
416*
417* Adjust the expected value of INFO to account for
418* pivoting.
419*
420 k = izero
421 IF( k.GT.0 ) THEN
422 100 CONTINUE
423 IF( iwork( k ).LT.0 ) THEN
424 IF( iwork( k ).NE.-k ) THEN
425 k = -iwork( k )
426 GO TO 100
427 END IF
428 ELSE IF( iwork( k ).NE.k ) THEN
429 k = iwork( k )
430 GO TO 100
431 END IF
432 END IF
433*
434* Check error code from DSYTRF and handle error.
435*
436 IF( info.NE.k )
437 $ CALL alaerh( path, 'DSYTRF', info, k, uplo, n, n,
438 $ -1, -1, nb, imat, nfail, nerrs, nout )
439*
440* Set the condition estimate flag if the INFO is not 0.
441*
442 IF( info.NE.0 ) THEN
443 trfcon = .true.
444 ELSE
445 trfcon = .false.
446 END IF
447*
448*+ TEST 1
449* Reconstruct matrix from factors and compute residual.
450*
451 CALL dsyt01( uplo, n, a, lda, afac, lda, iwork, ainv,
452 $ lda, rwork, result( 1 ) )
453 nt = 1
454*
455*+ TEST 2
456* Form the inverse and compute the residual,
457* if the factorization was competed without INFO > 0
458* (i.e. there is no zero rows and columns).
459* Do it only for the first block size.
460*
461 IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
462 CALL dlacpy( uplo, n, n, afac, lda, ainv, lda )
463 srnamt = 'DSYTRI2'
464 lwork = (n+nb+1)*(nb+3)
465 CALL dsytri2( uplo, n, ainv, lda, iwork, work,
466 $ lwork, info )
467*
468* Check error code from DSYTRI2 and handle error.
469*
470 IF( info.NE.0 )
471 $ CALL alaerh( path, 'DSYTRI2', info, -1, uplo, n,
472 $ n, -1, -1, -1, imat, nfail, nerrs,
473 $ nout )
474*
475* Compute the residual for a symmetric matrix times
476* its inverse.
477*
478 CALL dpot03( uplo, n, a, lda, ainv, lda, work, lda,
479 $ rwork, rcondc, result( 2 ) )
480 nt = 2
481 END IF
482*
483* Print information about the tests that did not pass
484* the threshold.
485*
486 DO 110 k = 1, nt
487 IF( result( k ).GE.thresh ) THEN
488 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
489 $ CALL alahd( nout, path )
490 WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
491 $ result( k )
492 nfail = nfail + 1
493 END IF
494 110 CONTINUE
495 nrun = nrun + nt
496*
497* Skip the other tests if this is not the first block
498* size.
499*
500 IF( inb.GT.1 )
501 $ GO TO 150
502*
503* Do only the condition estimate if INFO is not 0.
504*
505 IF( trfcon ) THEN
506 rcondc = zero
507 GO TO 140
508 END IF
509*
510* Do for each value of NRHS in NSVAL.
511*
512 DO 130 irhs = 1, nns
513 nrhs = nsval( irhs )
514*
515*+ TEST 3 ( Using TRS)
516* Solve and compute residual for A * X = B.
517*
518* Choose a set of NRHS random solution vectors
519* stored in XACT and set up the right hand side B
520*
521 srnamt = 'DLARHS'
522 CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
523 $ nrhs, a, lda, xact, lda, b, lda,
524 $ iseed, info )
525 CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
526*
527 srnamt = 'DSYTRS'
528 CALL dsytrs( uplo, n, nrhs, afac, lda, iwork, x,
529 $ lda, info )
530*
531* Check error code from DSYTRS and handle error.
532*
533 IF( info.NE.0 )
534 $ CALL alaerh( path, 'DSYTRS', info, 0, uplo, n,
535 $ n, -1, -1, nrhs, imat, nfail,
536 $ nerrs, nout )
537*
538 CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
539*
540* Compute the residual for the solution
541*
542 CALL dpot02( uplo, n, nrhs, a, lda, x, lda, work,
543 $ lda, rwork, result( 3 ) )
544*
545*+ TEST 4 (Using TRS2)
546*
547* Solve and compute residual for A * X = B.
548*
549* Choose a set of NRHS random solution vectors
550* stored in XACT and set up the right hand side B
551*
552 srnamt = 'DLARHS'
553 CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
554 $ nrhs, a, lda, xact, lda, b, lda,
555 $ iseed, info )
556 CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
557*
558 srnamt = 'DSYTRS2'
559 CALL dsytrs2( uplo, n, nrhs, afac, lda, iwork, x,
560 $ lda, work, info )
561*
562* Check error code from DSYTRS2 and handle error.
563*
564 IF( info.NE.0 )
565 $ CALL alaerh( path, 'DSYTRS2', info, 0, uplo, n,
566 $ n, -1, -1, nrhs, imat, nfail,
567 $ nerrs, nout )
568*
569 CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
570*
571* Compute the residual for the solution
572*
573 CALL dpot02( uplo, n, nrhs, a, lda, x, lda, work,
574 $ lda, rwork, result( 4 ) )
575*
576*+ TEST 5
577* Check solution from generated exact solution.
578*
579 CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
580 $ result( 5 ) )
581*
582*+ TESTS 6, 7, and 8
583* Use iterative refinement to improve the solution.
584*
585 srnamt = 'DSYRFS'
586 CALL dsyrfs( uplo, n, nrhs, a, lda, afac, lda,
587 $ iwork, b, lda, x, lda, rwork,
588 $ rwork( nrhs+1 ), work, iwork( n+1 ),
589 $ info )
590*
591* Check error code from DSYRFS and handle error.
592*
593 IF( info.NE.0 )
594 $ CALL alaerh( path, 'DSYRFS', info, 0, uplo, n,
595 $ n, -1, -1, nrhs, imat, nfail,
596 $ nerrs, nout )
597*
598 CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
599 $ result( 6 ) )
600 CALL dpot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
601 $ xact, lda, rwork, rwork( nrhs+1 ),
602 $ result( 7 ) )
603*
604* Print information about the tests that did not pass
605* the threshold.
606*
607 DO 120 k = 3, 8
608 IF( result( k ).GE.thresh ) THEN
609 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
610 $ CALL alahd( nout, path )
611 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
612 $ imat, k, result( k )
613 nfail = nfail + 1
614 END IF
615 120 CONTINUE
616 nrun = nrun + 6
617*
618* End do for each value of NRHS in NSVAL.
619*
620 130 CONTINUE
621*
622*+ TEST 9
623* Get an estimate of RCOND = 1/CNDNUM.
624*
625 140 CONTINUE
626 anorm = dlansy( '1', uplo, n, a, lda, rwork )
627 srnamt = 'DSYCON'
628 CALL dsycon( uplo, n, afac, lda, iwork, anorm, rcond,
629 $ work, iwork( n+1 ), info )
630*
631* Check error code from DSYCON and handle error.
632*
633 IF( info.NE.0 )
634 $ CALL alaerh( path, 'DSYCON', info, 0, uplo, n, n,
635 $ -1, -1, -1, imat, nfail, nerrs, nout )
636*
637* Compute the test ratio to compare values of RCOND
638*
639 result( 9 ) = dget06( rcond, rcondc )
640*
641* Print information about the tests that did not pass
642* the threshold.
643*
644 IF( result( 9 ).GE.thresh ) THEN
645 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
646 $ CALL alahd( nout, path )
647 WRITE( nout, fmt = 9997 )uplo, n, imat, 9,
648 $ result( 9 )
649 nfail = nfail + 1
650 END IF
651 nrun = nrun + 1
652 150 CONTINUE
653*
654 160 CONTINUE
655 170 CONTINUE
656 180 CONTINUE
657*
658* Print a summary of the results.
659*
660 CALL alasum( path, nout, nfail, nrun, nerrs )
661*
662 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
663 $ i2, ', test ', i2, ', ratio =', g12.5 )
664 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
665 $ i2, ', test(', i2, ') =', g12.5 )
666 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
667 $ ', test(', i2, ') =', g12.5 )
668 RETURN
669*
670* End of DCHKSY
671*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine dlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
DLARHS
Definition dlarhs.f:205
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine derrsy(path, nunit)
DERRSY
Definition derrsy.f:55
subroutine dget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
DGET04
Definition dget04.f:102
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
subroutine dlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
DLATB4
Definition dlatb4.f:120
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
DPOT02
Definition dpot02.f:127
subroutine dpot03(uplo, n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
DPOT03
Definition dpot03.f:125
subroutine dpot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
DPOT05
Definition dpot05.f:164
subroutine dsyt01(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
DSYT01
Definition dsyt01.f:124
subroutine dsycon(uplo, n, a, lda, ipiv, anorm, rcond, work, iwork, info)
DSYCON
Definition dsycon.f:130
subroutine dsyrfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DSYRFS
Definition dsyrfs.f:191
subroutine dsytrf(uplo, n, a, lda, ipiv, work, lwork, info)
DSYTRF
Definition dsytrf.f:182
subroutine dsytri2(uplo, n, a, lda, ipiv, work, lwork, info)
DSYTRI2
Definition dsytri2.f:127
subroutine dsytrs2(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, info)
DSYTRS2
Definition dsytrs2.f:132
subroutine dsytrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
DSYTRS
Definition dsytrs.f:120
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
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,...
Definition dlansy.f:122
Here is the call graph for this function:
Here is the caller graph for this function: