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

◆ schksy()

subroutine schksy ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
integer  nns,
integer, dimension( * )  nsval,
real  thresh,
logical  tsterr,
integer  nmax,
real, dimension( * )  a,
real, dimension( * )  afac,
real, dimension( * )  ainv,
real, dimension( * )  b,
real, dimension( * )  x,
real, dimension( * )  xact,
real, dimension( * )  work,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

SCHKSY

Purpose:
 SCHKSY tests SSYTRF, -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 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.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]AINV
          AINV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is REAL array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is REAL array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is REAL array, dimension (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL 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 schksy.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 REAL THRESH
179* ..
180* .. Array Arguments ..
181 LOGICAL DOTYPE( * )
182 INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
183 REAL A( * ), AFAC( * ), AINV( * ), B( * ),
184 $ RWORK( * ), WORK( * ), X( * ), XACT( * )
185* ..
186*
187* =====================================================================
188*
189* .. Parameters ..
190 REAL ZERO
191 parameter( zero = 0.0e+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 REAL ANORM, CNDNUM, RCOND, RCONDC
205* ..
206* .. Local Arrays ..
207 CHARACTER UPLOS( 2 )
208 INTEGER ISEED( 4 ), ISEEDY( 4 )
209 REAL RESULT( NTESTS )
210* ..
211* .. External Functions ..
212 REAL SGET06, SLANSY
213 EXTERNAL sget06, slansy
214* ..
215* .. External Subroutines ..
216 EXTERNAL alaerh, alahd, alasum, serrsy, sget04, slacpy,
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 ) = 'Single 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 serrsy( 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* Set up parameters with SLATB4 for the matrix generator
296* based on the type of matrix to be generated.
297*
298 CALL slatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
299 $ CNDNUM, DIST )
300*
301* Generate a matrix with SLATMS.
302*
303 srnamt = 'SLATMS'
304 CALL slatms( n, n, dist, iseed, TYPE, RWORK, MODE,
305 $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
306 $ INFO )
307*
308* Check error code from SLATMS and handle error.
309*
310 IF( info.NE.0 ) THEN
311 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
312 $ -1, -1, imat, nfail, nerrs, nout )
313*
314* Skip all tests for this generated matrix
315*
316 GO TO 160
317 END IF
318*
319* For matrix types 3-6, zero one or more rows and
320* columns of the matrix to test that INFO is returned
321* correctly.
322*
323 IF( zerot ) THEN
324 IF( imat.EQ.3 ) THEN
325 izero = 1
326 ELSE IF( imat.EQ.4 ) THEN
327 izero = n
328 ELSE
329 izero = n / 2 + 1
330 END IF
331*
332 IF( imat.LT.6 ) THEN
333*
334* Set row and column IZERO to zero.
335*
336 IF( iuplo.EQ.1 ) THEN
337 ioff = ( izero-1 )*lda
338 DO 20 i = 1, izero - 1
339 a( ioff+i ) = zero
340 20 CONTINUE
341 ioff = ioff + izero
342 DO 30 i = izero, n
343 a( ioff ) = zero
344 ioff = ioff + lda
345 30 CONTINUE
346 ELSE
347 ioff = izero
348 DO 40 i = 1, izero - 1
349 a( ioff ) = zero
350 ioff = ioff + lda
351 40 CONTINUE
352 ioff = ioff - izero
353 DO 50 i = izero, n
354 a( ioff+i ) = zero
355 50 CONTINUE
356 END IF
357 ELSE
358 IF( iuplo.EQ.1 ) THEN
359*
360* Set the first IZERO rows and columns to zero.
361*
362 ioff = 0
363 DO 70 j = 1, n
364 i2 = min( j, izero )
365 DO 60 i = 1, i2
366 a( ioff+i ) = zero
367 60 CONTINUE
368 ioff = ioff + lda
369 70 CONTINUE
370 ELSE
371*
372* Set the last IZERO rows and columns to zero.
373*
374 ioff = 0
375 DO 90 j = 1, n
376 i1 = max( j, izero )
377 DO 80 i = i1, n
378 a( ioff+i ) = zero
379 80 CONTINUE
380 ioff = ioff + lda
381 90 CONTINUE
382 END IF
383 END IF
384 ELSE
385 izero = 0
386 END IF
387*
388* End generate the test matrix A.
389*
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 slacpy( 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 = 'SSYTRF'
414 CALL ssytrf( 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 SSYTRF and handle error.
435*
436 IF( info.NE.k )
437 $ CALL alaerh( path, 'SSYTRF', 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 ssyt01( 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 slacpy( uplo, n, n, afac, lda, ainv, lda )
463 srnamt = 'SSYTRI2'
464 lwork = (n+nb+1)*(nb+3)
465 CALL ssytri2( uplo, n, ainv, lda, iwork, work,
466 $ lwork, info )
467*
468* Check error code from SSYTRI2 and handle error.
469*
470 IF( info.NE.0 )
471 $ CALL alaerh( path, 'SSYTRI2', 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 spot03( 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 DSYTRS)
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 = 'SLARHS'
522 CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
523 $ nrhs, a, lda, xact, lda, b, lda,
524 $ iseed, info )
525 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
526*
527 srnamt = 'SSYTRS'
528 CALL ssytrs( uplo, n, nrhs, afac, lda, iwork, x,
529 $ lda, info )
530*
531* Check error code from SSYTRS and handle error.
532*
533 IF( info.NE.0 )
534 $ CALL alaerh( path, 'SSYTRS', info, 0, uplo, n,
535 $ n, -1, -1, nrhs, imat, nfail,
536 $ nerrs, nout )
537*
538 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
539*
540* Compute the residual for the solution
541*
542 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
543 $ lda, rwork, result( 3 ) )
544*
545*+ TEST 4 (Using DSYTRS2)
546* Solve and compute residual for A * X = B.
547*
548* Choose a set of NRHS random solution vectors
549* stored in XACT and set up the right hand side B
550*
551 srnamt = 'SLARHS'
552 CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
553 $ nrhs, a, lda, xact, lda, b, lda,
554 $ iseed, info )
555 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
556*
557 srnamt = 'DSYTRS2'
558 CALL ssytrs2( uplo, n, nrhs, afac, lda, iwork, x,
559 $ lda, work, info )
560*
561* Check error code from SSYTRS2 and handle error.
562*
563 IF( info.NE.0 )
564 $ CALL alaerh( path, 'SSYTRS2', info, 0, uplo, n,
565 $ n, -1, -1, nrhs, imat, nfail,
566 $ nerrs, nout )
567*
568 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
569*
570* Compute the residual for the solution
571*
572 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
573 $ lda, rwork, result( 4 ) )
574*
575*+ TEST 5
576* Check solution from generated exact solution.
577*
578 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
579 $ result( 5 ) )
580*
581*+ TESTS 6, 7, and 8
582* Use iterative refinement to improve the solution.
583*
584 srnamt = 'SSYRFS'
585 CALL ssyrfs( uplo, n, nrhs, a, lda, afac, lda,
586 $ iwork, b, lda, x, lda, rwork,
587 $ rwork( nrhs+1 ), work, iwork( n+1 ),
588 $ info )
589*
590* Check error code from SSYRFS and handle error.
591*
592 IF( info.NE.0 )
593 $ CALL alaerh( path, 'SSYRFS', info, 0, uplo, n,
594 $ n, -1, -1, nrhs, imat, nfail,
595 $ nerrs, nout )
596*
597 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
598 $ result( 6 ) )
599 CALL spot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
600 $ xact, lda, rwork, rwork( nrhs+1 ),
601 $ result( 7 ) )
602*
603* Print information about the tests that did not pass
604* the threshold.
605*
606 DO 120 k = 3, 8
607 IF( result( k ).GE.thresh ) THEN
608 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
609 $ CALL alahd( nout, path )
610 WRITE( nout, fmt = 9998 )uplo, n, nrhs,
611 $ imat, k, result( k )
612 nfail = nfail + 1
613 END IF
614 120 CONTINUE
615 nrun = nrun + 6
616*
617* End do for each value of NRHS in NSVAL.
618*
619 130 CONTINUE
620*
621*+ TEST 9
622* Get an estimate of RCOND = 1/CNDNUM.
623*
624 140 CONTINUE
625 anorm = slansy( '1', uplo, n, a, lda, rwork )
626 srnamt = 'SSYCON'
627 CALL ssycon( uplo, n, afac, lda, iwork, anorm, rcond,
628 $ work, iwork( n+1 ), info )
629*
630* Check error code from SSYCON and handle error.
631*
632 IF( info.NE.0 )
633 $ CALL alaerh( path, 'SSYCON', info, 0, uplo, n, n,
634 $ -1, -1, -1, imat, nfail, nerrs, nout )
635*
636* Compute the test ratio to compare to values of RCOND
637*
638 result( 9 ) = sget06( rcond, rcondc )
639*
640* Print information about the tests that did not pass
641* the threshold.
642*
643 IF( result( 9 ).GE.thresh ) THEN
644 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
645 $ CALL alahd( nout, path )
646 WRITE( nout, fmt = 9997 )uplo, n, imat, 9,
647 $ result( 9 )
648 nfail = nfail + 1
649 END IF
650 nrun = nrun + 1
651 150 CONTINUE
652*
653 160 CONTINUE
654 170 CONTINUE
655 180 CONTINUE
656*
657* Print a summary of the results.
658*
659 CALL alasum( path, nout, nfail, nrun, nerrs )
660*
661 9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
662 $ i2, ', test ', i2, ', ratio =', g12.5 )
663 9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
664 $ i2, ', test(', i2, ') =', g12.5 )
665 9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
666 $ ', test(', i2, ') =', g12.5 )
667 RETURN
668*
669* End of SCHKSY
670*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine slarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
SLARHS
Definition slarhs.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 ssycon(uplo, n, a, lda, ipiv, anorm, rcond, work, iwork, info)
SSYCON
Definition ssycon.f:130
subroutine ssyrfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
SSYRFS
Definition ssyrfs.f:191
subroutine ssytrf(uplo, n, a, lda, ipiv, work, lwork, info)
SSYTRF
Definition ssytrf.f:182
subroutine ssytri2(uplo, n, a, lda, ipiv, work, lwork, info)
SSYTRI2
Definition ssytri2.f:127
subroutine ssytrs2(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, info)
SSYTRS2
Definition ssytrs2.f:132
subroutine ssytrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
SSYTRS
Definition ssytrs.f:120
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
real function slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:122
subroutine ssyconv(uplo, way, n, a, lda, ipiv, e, info)
SSYCONV
Definition ssyconv.f:114
subroutine serrsy(path, nunit)
SERRSY
Definition serrsy.f:55
subroutine sget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
SGET04
Definition sget04.f:102
real function sget06(rcond, rcondc)
SGET06
Definition sget06.f:55
subroutine slatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
SLATB4
Definition slatb4.f:120
subroutine slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321
subroutine spot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
SPOT02
Definition spot02.f:127
subroutine spot03(uplo, n, a, lda, ainv, ldainv, work, ldwork, rwork, rcond, resid)
SPOT03
Definition spot03.f:125
subroutine spot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
SPOT05
Definition spot05.f:164
subroutine ssyt01(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
SSYT01
Definition ssyt01.f:124
Here is the call graph for this function:
Here is the caller graph for this function: