LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
schksy.f
Go to the documentation of this file.
1 *> \brief \b SCHKSY
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE SCHKSY( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
12 * THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
13 * XACT, WORK, RWORK, IWORK, NOUT )
14 *
15 * .. Scalar Arguments ..
16 * LOGICAL TSTERR
17 * INTEGER NMAX, NN, NNB, NNS, NOUT
18 * REAL THRESH
19 * ..
20 * .. Array Arguments ..
21 * LOGICAL DOTYPE( * )
22 * INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
23 * REAL A( * ), AFAC( * ), AINV( * ), B( * ),
24 * $ RWORK( * ), WORK( * ), X( * ), XACT( * )
25 * ..
26 *
27 *
28 *> \par Purpose:
29 * =============
30 *>
31 *> \verbatim
32 *>
33 *> SCHKSY tests SSYTRF, -TRI2, -TRS, -TRS2, -RFS, and -CON.
34 *> \endverbatim
35 *
36 * Arguments:
37 * ==========
38 *
39 *> \param[in] DOTYPE
40 *> \verbatim
41 *> DOTYPE is LOGICAL array, dimension (NTYPES)
42 *> The matrix types to be used for testing. Matrices of type j
43 *> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
44 *> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
45 *> \endverbatim
46 *>
47 *> \param[in] NN
48 *> \verbatim
49 *> NN is INTEGER
50 *> The number of values of N contained in the vector NVAL.
51 *> \endverbatim
52 *>
53 *> \param[in] NVAL
54 *> \verbatim
55 *> NVAL is INTEGER array, dimension (NN)
56 *> The values of the matrix dimension N.
57 *> \endverbatim
58 *>
59 *> \param[in] NNB
60 *> \verbatim
61 *> NNB is INTEGER
62 *> The number of values of NB contained in the vector NBVAL.
63 *> \endverbatim
64 *>
65 *> \param[in] NBVAL
66 *> \verbatim
67 *> NBVAL is INTEGER array, dimension (NBVAL)
68 *> The values of the blocksize NB.
69 *> \endverbatim
70 *>
71 *> \param[in] NNS
72 *> \verbatim
73 *> NNS is INTEGER
74 *> The number of values of NRHS contained in the vector NSVAL.
75 *> \endverbatim
76 *>
77 *> \param[in] NSVAL
78 *> \verbatim
79 *> NSVAL is INTEGER array, dimension (NNS)
80 *> The values of the number of right hand sides NRHS.
81 *> \endverbatim
82 *>
83 *> \param[in] THRESH
84 *> \verbatim
85 *> THRESH is REAL
86 *> The threshold value for the test ratios. A result is
87 *> included in the output file if RESULT >= THRESH. To have
88 *> every test ratio printed, use THRESH = 0.
89 *> \endverbatim
90 *>
91 *> \param[in] TSTERR
92 *> \verbatim
93 *> TSTERR is LOGICAL
94 *> Flag that indicates whether error exits are to be tested.
95 *> \endverbatim
96 *>
97 *> \param[in] NMAX
98 *> \verbatim
99 *> NMAX is INTEGER
100 *> The maximum value permitted for N, used in dimensioning the
101 *> work arrays.
102 *> \endverbatim
103 *>
104 *> \param[out] A
105 *> \verbatim
106 *> A is REAL array, dimension (NMAX*NMAX)
107 *> \endverbatim
108 *>
109 *> \param[out] AFAC
110 *> \verbatim
111 *> AFAC is REAL array, dimension (NMAX*NMAX)
112 *> \endverbatim
113 *>
114 *> \param[out] AINV
115 *> \verbatim
116 *> AINV is REAL array, dimension (NMAX*NMAX)
117 *> \endverbatim
118 *>
119 *> \param[out] B
120 *> \verbatim
121 *> B is REAL array, dimension (NMAX*NSMAX)
122 *> where NSMAX is the largest entry in NSVAL.
123 *> \endverbatim
124 *>
125 *> \param[out] X
126 *> \verbatim
127 *> X is REAL array, dimension (NMAX*NSMAX)
128 *> \endverbatim
129 *>
130 *> \param[out] XACT
131 *> \verbatim
132 *> XACT is REAL array, dimension (NMAX*NSMAX)
133 *> \endverbatim
134 *>
135 *> \param[out] WORK
136 *> \verbatim
137 *> WORK is REAL array, dimension (NMAX*max(3,NSMAX))
138 *> \endverbatim
139 *>
140 *> \param[out] RWORK
141 *> \verbatim
142 *> RWORK is REAL array, dimension (max(NMAX,2*NSMAX))
143 *> \endverbatim
144 *>
145 *> \param[out] IWORK
146 *> \verbatim
147 *> IWORK is INTEGER array, dimension (2*NMAX)
148 *> \endverbatim
149 *>
150 *> \param[in] NOUT
151 *> \verbatim
152 *> NOUT is INTEGER
153 *> The unit number for output.
154 *> \endverbatim
155 *
156 * Authors:
157 * ========
158 *
159 *> \author Univ. of Tennessee
160 *> \author Univ. of California Berkeley
161 *> \author Univ. of Colorado Denver
162 *> \author NAG Ltd.
163 *
164 *> \date November 2013
165 *
166 *> \ingroup single_lin
167 *
168 * =====================================================================
169  SUBROUTINE schksy( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
170  $ thresh, tsterr, nmax, a, afac, ainv, b, x,
171  $ xact, work, rwork, iwork, nout )
172 *
173 * -- LAPACK test routine (version 3.5.0) --
174 * -- LAPACK is a software package provided by Univ. of Tennessee, --
175 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
176 * November 2013
177 *
178 * .. Scalar Arguments ..
179  LOGICAL TSTERR
180  INTEGER NMAX, NN, NNB, NNS, NOUT
181  REAL THRESH
182 * ..
183 * .. Array Arguments ..
184  LOGICAL DOTYPE( * )
185  INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
186  REAL A( * ), AFAC( * ), AINV( * ), B( * ),
187  $ rwork( * ), work( * ), x( * ), xact( * )
188 * ..
189 *
190 * =====================================================================
191 *
192 * .. Parameters ..
193  REAL ZERO
194  parameter ( zero = 0.0e+0 )
195  INTEGER NTYPES
196  parameter ( ntypes = 10 )
197  INTEGER NTESTS
198  parameter ( ntests = 9 )
199 * ..
200 * .. Local Scalars ..
201  LOGICAL TRFCON, ZEROT
202  CHARACTER DIST, TYPE, UPLO, XTYPE
203  CHARACTER*3 PATH
204  INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
205  $ iuplo, izero, j, k, kl, ku, lda, lwork, mode,
206  $ n, nb, nerrs, nfail, nimat, nrhs, nrun, nt
207  REAL ANORM, CNDNUM, RCOND, RCONDC
208 * ..
209 * .. Local Arrays ..
210  CHARACTER UPLOS( 2 )
211  INTEGER ISEED( 4 ), ISEEDY( 4 )
212  REAL RESULT( ntests )
213 * ..
214 * .. External Functions ..
215  REAL SGET06, SLANSY
216  EXTERNAL sget06, slansy
217 * ..
218 * .. External Subroutines ..
219  EXTERNAL alaerh, alahd, alasum, serrsy, sget04, slacpy,
223 * ..
224 * .. Intrinsic Functions ..
225  INTRINSIC max, min
226 * ..
227 * .. Scalars in Common ..
228  LOGICAL LERR, OK
229  CHARACTER*32 SRNAMT
230  INTEGER INFOT, NUNIT
231 * ..
232 * .. Common blocks ..
233  COMMON / infoc / infot, nunit, ok, lerr
234  COMMON / srnamc / srnamt
235 * ..
236 * .. Data statements ..
237  DATA iseedy / 1988, 1989, 1990, 1991 /
238  DATA uplos / 'U', 'L' /
239 * ..
240 * .. Executable Statements ..
241 *
242 * Initialize constants and the random number seed.
243 *
244  path( 1: 1 ) = 'Single precision'
245  path( 2: 3 ) = 'SY'
246  nrun = 0
247  nfail = 0
248  nerrs = 0
249  DO 10 i = 1, 4
250  iseed( i ) = iseedy( i )
251  10 CONTINUE
252 *
253 * Test the error exits
254 *
255  IF( tsterr )
256  $ CALL serrsy( path, nout )
257  infot = 0
258 *
259 * Set the minimum block size for which the block routine should
260 * be used, which will be later returned by ILAENV
261 *
262  CALL xlaenv( 2, 2 )
263 *
264 * Do for each value of N in NVAL
265 *
266  DO 180 in = 1, nn
267  n = nval( in )
268  lda = max( n, 1 )
269  xtype = 'N'
270  nimat = ntypes
271  IF( n.LE.0 )
272  $ nimat = 1
273 *
274  izero = 0
275 *
276 * Do for each value of matrix type IMAT
277 *
278  DO 170 imat = 1, nimat
279 *
280 * Do the tests only if DOTYPE( IMAT ) is true.
281 *
282  IF( .NOT.dotype( imat ) )
283  $ GO TO 170
284 *
285 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
286 *
287  zerot = imat.GE.3 .AND. imat.LE.6
288  IF( zerot .AND. n.LT.imat-2 )
289  $ GO TO 170
290 *
291 * Do first for UPLO = 'U', then for UPLO = 'L'
292 *
293  DO 160 iuplo = 1, 2
294  uplo = uplos( iuplo )
295 *
296 * Begin generate the test matrix A.
297 *
298 * Set up parameters with SLATB4 for the matrix generator
299 * based on the type of matrix to be generated.
300 *
301  CALL slatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
302  $ cndnum, dist )
303 *
304 * Generate a matrix with SLATMS.
305 *
306  srnamt = 'SLATMS'
307  CALL slatms( n, n, dist, iseed, TYPE, RWORK, MODE,
308  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
309  $ info )
310 *
311 * Check error code from SLATMS and handle error.
312 *
313  IF( info.NE.0 ) THEN
314  CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
315  $ -1, -1, imat, nfail, nerrs, nout )
316 *
317 * Skip all tests for this generated matrix
318 *
319  GO TO 160
320  END IF
321 *
322 * For matrix types 3-6, zero one or more rows and
323 * columns of the matrix to test that INFO is returned
324 * correctly.
325 *
326  IF( zerot ) THEN
327  IF( imat.EQ.3 ) THEN
328  izero = 1
329  ELSE IF( imat.EQ.4 ) THEN
330  izero = n
331  ELSE
332  izero = n / 2 + 1
333  END IF
334 *
335  IF( imat.LT.6 ) THEN
336 *
337 * Set row and column IZERO to zero.
338 *
339  IF( iuplo.EQ.1 ) THEN
340  ioff = ( izero-1 )*lda
341  DO 20 i = 1, izero - 1
342  a( ioff+i ) = zero
343  20 CONTINUE
344  ioff = ioff + izero
345  DO 30 i = izero, n
346  a( ioff ) = zero
347  ioff = ioff + lda
348  30 CONTINUE
349  ELSE
350  ioff = izero
351  DO 40 i = 1, izero - 1
352  a( ioff ) = zero
353  ioff = ioff + lda
354  40 CONTINUE
355  ioff = ioff - izero
356  DO 50 i = izero, n
357  a( ioff+i ) = zero
358  50 CONTINUE
359  END IF
360  ELSE
361  IF( iuplo.EQ.1 ) THEN
362 *
363 * Set the first IZERO rows and columns to zero.
364 *
365  ioff = 0
366  DO 70 j = 1, n
367  i2 = min( j, izero )
368  DO 60 i = 1, i2
369  a( ioff+i ) = zero
370  60 CONTINUE
371  ioff = ioff + lda
372  70 CONTINUE
373  ELSE
374 *
375 * Set the last IZERO rows and columns to zero.
376 *
377  ioff = 0
378  DO 90 j = 1, n
379  i1 = max( j, izero )
380  DO 80 i = i1, n
381  a( ioff+i ) = zero
382  80 CONTINUE
383  ioff = ioff + lda
384  90 CONTINUE
385  END IF
386  END IF
387  ELSE
388  izero = 0
389  END IF
390 *
391 * End generate the test matrix A.
392 *
393 *
394 * Do for each value of NB in NBVAL
395 *
396  DO 150 inb = 1, nnb
397 *
398 * Set the optimal blocksize, which will be later
399 * returned by ILAENV.
400 *
401  nb = nbval( inb )
402  CALL xlaenv( 1, nb )
403 *
404 * Copy the test matrix A into matrix AFAC which
405 * will be factorized in place. This is needed to
406 * preserve the test matrix A for subsequent tests.
407 *
408  CALL slacpy( uplo, n, n, a, lda, afac, lda )
409 *
410 * Compute the L*D*L**T or U*D*U**T factorization of the
411 * matrix. IWORK stores details of the interchanges and
412 * the block structure of D. AINV is a work array for
413 * block factorization, LWORK is the length of AINV.
414 *
415  lwork = max( 2, nb )*lda
416  srnamt = 'SSYTRF'
417  CALL ssytrf( uplo, n, afac, lda, iwork, ainv, lwork,
418  $ info )
419 *
420 * Adjust the expected value of INFO to account for
421 * pivoting.
422 *
423  k = izero
424  IF( k.GT.0 ) THEN
425  100 CONTINUE
426  IF( iwork( k ).LT.0 ) THEN
427  IF( iwork( k ).NE.-k ) THEN
428  k = -iwork( k )
429  GO TO 100
430  END IF
431  ELSE IF( iwork( k ).NE.k ) THEN
432  k = iwork( k )
433  GO TO 100
434  END IF
435  END IF
436 *
437 * Check error code from SSYTRF and handle error.
438 *
439  IF( info.NE.k )
440  $ CALL alaerh( path, 'SSYTRF', info, k, uplo, n, n,
441  $ -1, -1, nb, imat, nfail, nerrs, nout )
442 *
443 * Set the condition estimate flag if the INFO is not 0.
444 *
445  IF( info.NE.0 ) THEN
446  trfcon = .true.
447  ELSE
448  trfcon = .false.
449  END IF
450 *
451 *+ TEST 1
452 * Reconstruct matrix from factors and compute residual.
453 *
454  CALL ssyt01( uplo, n, a, lda, afac, lda, iwork, ainv,
455  $ lda, rwork, result( 1 ) )
456  nt = 1
457 *
458 *+ TEST 2
459 * Form the inverse and compute the residual,
460 * if the factorization was competed without INFO > 0
461 * (i.e. there is no zero rows and columns).
462 * Do it only for the first block size.
463 *
464  IF( inb.EQ.1 .AND. .NOT.trfcon ) THEN
465  CALL slacpy( uplo, n, n, afac, lda, ainv, lda )
466  srnamt = 'SSYTRI2'
467  lwork = (n+nb+1)*(nb+3)
468  CALL ssytri2( uplo, n, ainv, lda, iwork, work,
469  $ lwork, info )
470 *
471 * Check error code from SSYTRI2 and handle error.
472 *
473  IF( info.NE.0 )
474  $ CALL alaerh( path, 'SSYTRI2', info, -1, uplo, n,
475  $ n, -1, -1, -1, imat, nfail, nerrs,
476  $ nout )
477 *
478 * Compute the residual for a symmetric matrix times
479 * its inverse.
480 *
481  CALL spot03( uplo, n, a, lda, ainv, lda, work, lda,
482  $ rwork, rcondc, result( 2 ) )
483  nt = 2
484  END IF
485 *
486 * Print information about the tests that did not pass
487 * the threshold.
488 *
489  DO 110 k = 1, nt
490  IF( result( k ).GE.thresh ) THEN
491  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
492  $ CALL alahd( nout, path )
493  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
494  $ result( k )
495  nfail = nfail + 1
496  END IF
497  110 CONTINUE
498  nrun = nrun + nt
499 *
500 * Skip the other tests if this is not the first block
501 * size.
502 *
503  IF( inb.GT.1 )
504  $ GO TO 150
505 *
506 * Do only the condition estimate if INFO is not 0.
507 *
508  IF( trfcon ) THEN
509  rcondc = zero
510  GO TO 140
511  END IF
512 *
513 * Do for each value of NRHS in NSVAL.
514 *
515  DO 130 irhs = 1, nns
516  nrhs = nsval( irhs )
517 *
518 *+ TEST 3 (Using DSYTRS)
519 * Solve and compute residual for A * X = B.
520 *
521 * Choose a set of NRHS random solution vectors
522 * stored in XACT and set up the right hand side B
523 *
524  srnamt = 'SLARHS'
525  CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
526  $ nrhs, a, lda, xact, lda, b, lda,
527  $ iseed, info )
528  CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
529 *
530  srnamt = 'SSYTRS'
531  CALL ssytrs( uplo, n, nrhs, afac, lda, iwork, x,
532  $ lda, info )
533 *
534 * Check error code from SSYTRS and handle error.
535 *
536  IF( info.NE.0 )
537  $ CALL alaerh( path, 'SSYTRS', info, 0, uplo, n,
538  $ n, -1, -1, nrhs, imat, nfail,
539  $ nerrs, nout )
540 *
541  CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
542 *
543 * Compute the residual for the solution
544 *
545  CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
546  $ lda, rwork, result( 3 ) )
547 *
548 *+ TEST 4 (Using DSYTRS2)
549 * Solve and compute residual for A * X = B.
550 *
551 * Choose a set of NRHS random solution vectors
552 * stored in XACT and set up the right hand side B
553 *
554  srnamt = 'SLARHS'
555  CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
556  $ nrhs, a, lda, xact, lda, b, lda,
557  $ iseed, info )
558  CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
559 *
560  srnamt = 'DSYTRS2'
561  CALL ssytrs2( uplo, n, nrhs, afac, lda, iwork, x,
562  $ lda, work, info )
563 *
564 * Check error code from SSYTRS2 and handle error.
565 *
566  IF( info.NE.0 )
567  $ CALL alaerh( path, 'SSYTRS2', info, 0, uplo, n,
568  $ n, -1, -1, nrhs, imat, nfail,
569  $ nerrs, nout )
570 *
571  CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
572 *
573 * Compute the residual for the solution
574 *
575  CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
576  $ lda, rwork, result( 4 ) )
577 *
578 *+ TEST 5
579 * Check solution from generated exact solution.
580 *
581  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
582  $ result( 5 ) )
583 *
584 *+ TESTS 6, 7, and 8
585 * Use iterative refinement to improve the solution.
586 *
587  srnamt = 'SSYRFS'
588  CALL ssyrfs( uplo, n, nrhs, a, lda, afac, lda,
589  $ iwork, b, lda, x, lda, rwork,
590  $ rwork( nrhs+1 ), work, iwork( n+1 ),
591  $ info )
592 *
593 * Check error code from SSYRFS and handle error.
594 *
595  IF( info.NE.0 )
596  $ CALL alaerh( path, 'SSYRFS', info, 0, uplo, n,
597  $ n, -1, -1, nrhs, imat, nfail,
598  $ nerrs, nout )
599 *
600  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
601  $ result( 6 ) )
602  CALL spot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
603  $ xact, lda, rwork, rwork( nrhs+1 ),
604  $ result( 7 ) )
605 *
606 * Print information about the tests that did not pass
607 * the threshold.
608 *
609  DO 120 k = 3, 8
610  IF( result( k ).GE.thresh ) THEN
611  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
612  $ CALL alahd( nout, path )
613  WRITE( nout, fmt = 9998 )uplo, n, nrhs,
614  $ imat, k, result( k )
615  nfail = nfail + 1
616  END IF
617  120 CONTINUE
618  nrun = nrun + 6
619 *
620 * End do for each value of NRHS in NSVAL.
621 *
622  130 CONTINUE
623 *
624 *+ TEST 9
625 * Get an estimate of RCOND = 1/CNDNUM.
626 *
627  140 CONTINUE
628  anorm = slansy( '1', uplo, n, a, lda, rwork )
629  srnamt = 'SSYCON'
630  CALL ssycon( uplo, n, afac, lda, iwork, anorm, rcond,
631  $ work, iwork( n+1 ), info )
632 *
633 * Check error code from SSYCON and handle error.
634 *
635  IF( info.NE.0 )
636  $ CALL alaerh( path, 'SSYCON', info, 0, uplo, n, n,
637  $ -1, -1, -1, imat, nfail, nerrs, nout )
638 *
639 * Compute the test ratio to compare to values of RCOND
640 *
641  result( 9 ) = sget06( rcond, rcondc )
642 *
643 * Print information about the tests that did not pass
644 * the threshold.
645 *
646  IF( result( 9 ).GE.thresh ) THEN
647  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
648  $ CALL alahd( nout, path )
649  WRITE( nout, fmt = 9997 )uplo, n, imat, 9,
650  $ result( 9 )
651  nfail = nfail + 1
652  END IF
653  nrun = nrun + 1
654  150 CONTINUE
655 *
656  160 CONTINUE
657  170 CONTINUE
658  180 CONTINUE
659 *
660 * Print a summary of the results.
661 *
662  CALL alasum( path, nout, nfail, nrun, nerrs )
663 *
664  9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
665  $ i2, ', test ', i2, ', ratio =', g12.5 )
666  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
667  $ i2, ', test(', i2, ') =', g12.5 )
668  9997 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
669  $ ', test(', i2, ') =', g12.5 )
670  RETURN
671 *
672 * End of SCHKSY
673 *
674  END
subroutine ssytri2(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
SSYTRI2
Definition: ssytri2.f:129
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 slatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
SLATB4
Definition: slatb4.f:122
subroutine ssytrs(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
SSYTRS
Definition: ssytrs.f:122
subroutine ssytrf(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
SSYTRF
Definition: ssytrf.f:184
subroutine ssyconv(UPLO, WAY, N, A, LDA, IPIV, E, INFO)
SSYCONV
Definition: ssyconv.f:116
subroutine slarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
SLARHS
Definition: slarhs.f:206
subroutine ssycon(UPLO, N, A, LDA, IPIV, ANORM, RCOND, WORK, IWORK, INFO)
SSYCON
Definition: ssycon.f:132
subroutine spot05(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
SPOT05
Definition: spot05.f:166
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine sget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
SGET04
Definition: sget04.f:104
subroutine spot03(UPLO, N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
SPOT03
Definition: spot03.f:127
subroutine schksy(DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL, THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK, NOUT)
SCHKSY
Definition: schksy.f:172
subroutine serrsy(PATH, NUNIT)
SERRSY
Definition: serrsy.f:57
subroutine spot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
SPOT02
Definition: spot02.f:129
subroutine ssyrfs(UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
SSYRFS
Definition: ssyrfs.f:193
subroutine ssyt01(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
SSYT01
Definition: ssyt01.f:126
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75
subroutine ssytrs2(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK, INFO)
SSYTRS2
Definition: ssytrs2.f:134