LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvsyx.f
Go to the documentation of this file.
1*> \brief \b SDRVSYX
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 SDRVSY( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12* A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK,
13* NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NOUT, NRHS
18* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NVAL( * )
23* REAL A( * ), AFAC( * ), AINV( * ), B( * ),
24* $ RWORK( * ), WORK( * ), X( * ), XACT( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> SDRVSY tests the driver routines SSYSV, -SVX, and -SVXX
34*>
35*> Note that this file is used only when the XBLAS are available,
36*> otherwise sdrvsy.f defines this subroutine.
37*> \endverbatim
38*
39* Arguments:
40* ==========
41*
42*> \param[in] DOTYPE
43*> \verbatim
44*> DOTYPE is LOGICAL array, dimension (NTYPES)
45*> The matrix types to be used for testing. Matrices of type j
46*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
47*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
48*> \endverbatim
49*>
50*> \param[in] NN
51*> \verbatim
52*> NN is INTEGER
53*> The number of values of N contained in the vector NVAL.
54*> \endverbatim
55*>
56*> \param[in] NVAL
57*> \verbatim
58*> NVAL is INTEGER array, dimension (NN)
59*> The values of the matrix dimension N.
60*> \endverbatim
61*>
62*> \param[in] NRHS
63*> \verbatim
64*> NRHS is INTEGER
65*> The number of right hand side vectors to be generated for
66*> each linear system.
67*> \endverbatim
68*>
69*> \param[in] THRESH
70*> \verbatim
71*> THRESH is REAL
72*> The threshold value for the test ratios. A result is
73*> included in the output file if RESULT >= THRESH. To have
74*> every test ratio printed, use THRESH = 0.
75*> \endverbatim
76*>
77*> \param[in] TSTERR
78*> \verbatim
79*> TSTERR is LOGICAL
80*> Flag that indicates whether error exits are to be tested.
81*> \endverbatim
82*>
83*> \param[in] NMAX
84*> \verbatim
85*> NMAX is INTEGER
86*> The maximum value permitted for N, used in dimensioning the
87*> work arrays.
88*> \endverbatim
89*>
90*> \param[out] A
91*> \verbatim
92*> A is REAL array, dimension (NMAX*NMAX)
93*> \endverbatim
94*>
95*> \param[out] AFAC
96*> \verbatim
97*> AFAC is REAL array, dimension (NMAX*NMAX)
98*> \endverbatim
99*>
100*> \param[out] AINV
101*> \verbatim
102*> AINV is REAL array, dimension (NMAX*NMAX)
103*> \endverbatim
104*>
105*> \param[out] B
106*> \verbatim
107*> B is REAL array, dimension (NMAX*NRHS)
108*> \endverbatim
109*>
110*> \param[out] X
111*> \verbatim
112*> X is REAL array, dimension (NMAX*NRHS)
113*> \endverbatim
114*>
115*> \param[out] XACT
116*> \verbatim
117*> XACT is REAL array, dimension (NMAX*NRHS)
118*> \endverbatim
119*>
120*> \param[out] WORK
121*> \verbatim
122*> WORK is REAL array, dimension
123*> (NMAX*max(2,NRHS))
124*> \endverbatim
125*>
126*> \param[out] RWORK
127*> \verbatim
128*> RWORK is REAL array, dimension (NMAX+2*NRHS)
129*> \endverbatim
130*>
131*> \param[out] IWORK
132*> \verbatim
133*> IWORK is INTEGER array, dimension (2*NMAX)
134*> \endverbatim
135*>
136*> \param[in] NOUT
137*> \verbatim
138*> NOUT is INTEGER
139*> The unit number for output.
140*> \endverbatim
141*
142* Authors:
143* ========
144*
145*> \author Univ. of Tennessee
146*> \author Univ. of California Berkeley
147*> \author Univ. of Colorado Denver
148*> \author NAG Ltd.
149*
150*> \ingroup single_lin
151*
152* =====================================================================
153 SUBROUTINE sdrvsy( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
154 $ A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK,
155 $ NOUT )
156*
157* -- LAPACK test routine --
158* -- LAPACK is a software package provided by Univ. of Tennessee, --
159* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
160*
161* .. Scalar Arguments ..
162 LOGICAL TSTERR
163 INTEGER NMAX, NN, NOUT, NRHS
164 REAL THRESH
165* ..
166* .. Array Arguments ..
167 LOGICAL DOTYPE( * )
168 INTEGER IWORK( * ), NVAL( * )
169 REAL A( * ), AFAC( * ), AINV( * ), B( * ),
170 $ rwork( * ), work( * ), x( * ), xact( * )
171* ..
172*
173* =====================================================================
174*
175* .. Parameters ..
176 REAL ONE, ZERO
177 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
178 INTEGER NTYPES, NTESTS
179 parameter( ntypes = 10, ntests = 6 )
180 INTEGER NFACT
181 parameter( nfact = 2 )
182* ..
183* .. Local Scalars ..
184 LOGICAL ZEROT
185 CHARACTER DIST, EQUED, FACT, TYPE, UPLO, XTYPE
186 CHARACTER*3 PATH
187 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
188 $ izero, j, k, k1, kl, ku, lda, lwork, mode, n,
189 $ nb, nbmin, nerrs, nfail, nimat, nrun, nt,
190 $ n_err_bnds
191 REAL AINVNM, ANORM, CNDNUM, RCOND, RCONDC,
192 $ RPVGRW_SVXX
193* ..
194* .. Local Arrays ..
195 CHARACTER FACTS( NFACT ), UPLOS( 2 )
196 INTEGER ISEED( 4 ), ISEEDY( 4 )
197 REAL RESULT( NTESTS ), BERR( NRHS ),
198 $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
199* ..
200* .. External Functions ..
201 REAL SGET06, SLANSY
202 EXTERNAL SGET06, SLANSY
203* ..
204* .. External Subroutines ..
205 EXTERNAL aladhd, alaerh, alasvm, serrvx, sget04, slacpy,
208 $ ssysvxx
209* ..
210* .. Scalars in Common ..
211 LOGICAL LERR, OK
212 CHARACTER*32 SRNAMT
213 INTEGER INFOT, NUNIT
214* ..
215* .. Common blocks ..
216 COMMON / infoc / infot, nunit, ok, lerr
217 COMMON / srnamc / srnamt
218* ..
219* .. Intrinsic Functions ..
220 INTRINSIC max, min
221* ..
222* .. Data statements ..
223 DATA iseedy / 1988, 1989, 1990, 1991 /
224 DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
225* ..
226* .. Executable Statements ..
227*
228* Initialize constants and the random number seed.
229*
230 path( 1: 1 ) = 'Single precision'
231 path( 2: 3 ) = 'SY'
232 nrun = 0
233 nfail = 0
234 nerrs = 0
235 DO 10 i = 1, 4
236 iseed( i ) = iseedy( i )
237 10 CONTINUE
238 lwork = max( 2*nmax, nmax*nrhs )
239*
240* Test the error exits
241*
242 IF( tsterr )
243 $ CALL serrvx( path, nout )
244 infot = 0
245*
246* Set the block size and minimum block size for testing.
247*
248 nb = 1
249 nbmin = 2
250 CALL xlaenv( 1, nb )
251 CALL xlaenv( 2, nbmin )
252*
253* Do for each value of N in NVAL
254*
255 DO 180 in = 1, nn
256 n = nval( in )
257 lda = max( n, 1 )
258 xtype = 'N'
259 nimat = ntypes
260 IF( n.LE.0 )
261 $ nimat = 1
262*
263 DO 170 imat = 1, nimat
264*
265* Do the tests only if DOTYPE( IMAT ) is true.
266*
267 IF( .NOT.dotype( imat ) )
268 $ GO TO 170
269*
270* Skip types 3, 4, 5, or 6 if the matrix size is too small.
271*
272 zerot = imat.GE.3 .AND. imat.LE.6
273 IF( zerot .AND. n.LT.imat-2 )
274 $ GO TO 170
275*
276* Do first for UPLO = 'U', then for UPLO = 'L'
277*
278 DO 160 iuplo = 1, 2
279 uplo = uplos( iuplo )
280*
281* Set up parameters with SLATB4 and generate a test matrix
282* with SLATMS.
283*
284 CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
285 $ cndnum, dist )
286*
287 srnamt = 'SLATMS'
288 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
289 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
290 $ info )
291*
292* Check error code from SLATMS.
293*
294 IF( info.NE.0 ) THEN
295 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
296 $ -1, -1, imat, nfail, nerrs, nout )
297 GO TO 160
298 END IF
299*
300* For types 3-6, zero one or more rows and columns of the
301* matrix to test that INFO is returned correctly.
302*
303 IF( zerot ) THEN
304 IF( imat.EQ.3 ) THEN
305 izero = 1
306 ELSE IF( imat.EQ.4 ) THEN
307 izero = n
308 ELSE
309 izero = n / 2 + 1
310 END IF
311*
312 IF( imat.LT.6 ) THEN
313*
314* Set row and column IZERO to zero.
315*
316 IF( iuplo.EQ.1 ) THEN
317 ioff = ( izero-1 )*lda
318 DO 20 i = 1, izero - 1
319 a( ioff+i ) = zero
320 20 CONTINUE
321 ioff = ioff + izero
322 DO 30 i = izero, n
323 a( ioff ) = zero
324 ioff = ioff + lda
325 30 CONTINUE
326 ELSE
327 ioff = izero
328 DO 40 i = 1, izero - 1
329 a( ioff ) = zero
330 ioff = ioff + lda
331 40 CONTINUE
332 ioff = ioff - izero
333 DO 50 i = izero, n
334 a( ioff+i ) = zero
335 50 CONTINUE
336 END IF
337 ELSE
338 ioff = 0
339 IF( iuplo.EQ.1 ) THEN
340*
341* Set the first IZERO rows and columns to zero.
342*
343 DO 70 j = 1, n
344 i2 = min( j, izero )
345 DO 60 i = 1, i2
346 a( ioff+i ) = zero
347 60 CONTINUE
348 ioff = ioff + lda
349 70 CONTINUE
350 ELSE
351*
352* Set the last IZERO rows and columns to zero.
353*
354 DO 90 j = 1, n
355 i1 = max( j, izero )
356 DO 80 i = i1, n
357 a( ioff+i ) = zero
358 80 CONTINUE
359 ioff = ioff + lda
360 90 CONTINUE
361 END IF
362 END IF
363 ELSE
364 izero = 0
365 END IF
366*
367 DO 150 ifact = 1, nfact
368*
369* Do first for FACT = 'F', then for other values.
370*
371 fact = facts( ifact )
372*
373* Compute the condition number for comparison with
374* the value returned by SSYSVX.
375*
376 IF( zerot ) THEN
377 IF( ifact.EQ.1 )
378 $ GO TO 150
379 rcondc = zero
380*
381 ELSE IF( ifact.EQ.1 ) THEN
382*
383* Compute the 1-norm of A.
384*
385 anorm = slansy( '1', uplo, n, a, lda, rwork )
386*
387* Factor the matrix A.
388*
389 CALL slacpy( uplo, n, n, a, lda, afac, lda )
390 CALL ssytrf( uplo, n, afac, lda, iwork, work,
391 $ lwork, info )
392*
393* Compute inv(A) and take its norm.
394*
395 CALL slacpy( uplo, n, n, afac, lda, ainv, lda )
396 lwork = (n+nb+1)*(nb+3)
397 CALL ssytri2( uplo, n, ainv, lda, iwork, work,
398 $ lwork, info )
399 ainvnm = slansy( '1', uplo, n, ainv, lda, rwork )
400*
401* Compute the 1-norm condition number of A.
402*
403 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
404 rcondc = one
405 ELSE
406 rcondc = ( one / anorm ) / ainvnm
407 END IF
408 END IF
409*
410* Form an exact solution and set the right hand side.
411*
412 srnamt = 'SLARHS'
413 CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
414 $ nrhs, a, lda, xact, lda, b, lda, iseed,
415 $ info )
416 xtype = 'C'
417*
418* --- Test SSYSV ---
419*
420 IF( ifact.EQ.2 ) THEN
421 CALL slacpy( uplo, n, n, a, lda, afac, lda )
422 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
423*
424* Factor the matrix and solve the system using SSYSV.
425*
426 srnamt = 'SSYSV '
427 CALL ssysv( uplo, n, nrhs, afac, lda, iwork, x,
428 $ lda, work, lwork, info )
429*
430* Adjust the expected value of INFO to account for
431* pivoting.
432*
433 k = izero
434 IF( k.GT.0 ) THEN
435 100 CONTINUE
436 IF( iwork( k ).LT.0 ) THEN
437 IF( iwork( k ).NE.-k ) THEN
438 k = -iwork( k )
439 GO TO 100
440 END IF
441 ELSE IF( iwork( k ).NE.k ) THEN
442 k = iwork( k )
443 GO TO 100
444 END IF
445 END IF
446*
447* Check error code from SSYSV .
448*
449 IF( info.NE.k ) THEN
450 CALL alaerh( path, 'SSYSV ', info, k, uplo, n,
451 $ n, -1, -1, nrhs, imat, nfail,
452 $ nerrs, nout )
453 GO TO 120
454 ELSE IF( info.NE.0 ) THEN
455 GO TO 120
456 END IF
457*
458* Reconstruct matrix from factors and compute
459* residual.
460*
461 CALL ssyt01( uplo, n, a, lda, afac, lda, iwork,
462 $ ainv, lda, rwork, result( 1 ) )
463*
464* Compute residual of the computed solution.
465*
466 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
467 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
468 $ lda, rwork, result( 2 ) )
469*
470* Check solution from generated exact solution.
471*
472 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
473 $ result( 3 ) )
474 nt = 3
475*
476* Print information about the tests that did not pass
477* the threshold.
478*
479 DO 110 k = 1, nt
480 IF( result( k ).GE.thresh ) THEN
481 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
482 $ CALL aladhd( nout, path )
483 WRITE( nout, fmt = 9999 )'SSYSV ', uplo, n,
484 $ imat, k, result( k )
485 nfail = nfail + 1
486 END IF
487 110 CONTINUE
488 nrun = nrun + nt
489 120 CONTINUE
490 END IF
491*
492* --- Test SSYSVX ---
493*
494 IF( ifact.EQ.2 )
495 $ CALL slaset( uplo, n, n, zero, zero, afac, lda )
496 CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
497*
498* Solve the system and compute the condition number and
499* error bounds using SSYSVX.
500*
501 srnamt = 'SSYSVX'
502 CALL ssysvx( fact, uplo, n, nrhs, a, lda, afac, lda,
503 $ iwork, b, lda, x, lda, rcond, rwork,
504 $ rwork( nrhs+1 ), work, lwork,
505 $ iwork( n+1 ), info )
506*
507* Adjust the expected value of INFO to account for
508* pivoting.
509*
510 k = izero
511 IF( k.GT.0 ) THEN
512 130 CONTINUE
513 IF( iwork( k ).LT.0 ) THEN
514 IF( iwork( k ).NE.-k ) THEN
515 k = -iwork( k )
516 GO TO 130
517 END IF
518 ELSE IF( iwork( k ).NE.k ) THEN
519 k = iwork( k )
520 GO TO 130
521 END IF
522 END IF
523*
524* Check the error code from SSYSVX.
525*
526 IF( info.NE.k ) THEN
527 CALL alaerh( path, 'SSYSVX', info, k, fact // uplo,
528 $ n, n, -1, -1, nrhs, imat, nfail,
529 $ nerrs, nout )
530 GO TO 150
531 END IF
532*
533 IF( info.EQ.0 ) THEN
534 IF( ifact.GE.2 ) THEN
535*
536* Reconstruct matrix from factors and compute
537* residual.
538*
539 CALL ssyt01( uplo, n, a, lda, afac, lda, iwork,
540 $ ainv, lda, rwork( 2*nrhs+1 ),
541 $ result( 1 ) )
542 k1 = 1
543 ELSE
544 k1 = 2
545 END IF
546*
547* Compute residual of the computed solution.
548*
549 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
550 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
551 $ lda, rwork( 2*nrhs+1 ), result( 2 ) )
552*
553* Check solution from generated exact solution.
554*
555 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
556 $ result( 3 ) )
557*
558* Check the error bounds from iterative refinement.
559*
560 CALL spot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
561 $ xact, lda, rwork, rwork( nrhs+1 ),
562 $ result( 4 ) )
563 ELSE
564 k1 = 6
565 END IF
566*
567* Compare RCOND from SSYSVX with the computed value
568* in RCONDC.
569*
570 result( 6 ) = sget06( rcond, rcondc )
571*
572* Print information about the tests that did not pass
573* the threshold.
574*
575 DO 140 k = k1, 6
576 IF( result( k ).GE.thresh ) THEN
577 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
578 $ CALL aladhd( nout, path )
579 WRITE( nout, fmt = 9998 )'SSYSVX', fact, uplo,
580 $ n, imat, k, result( k )
581 nfail = nfail + 1
582 END IF
583 140 CONTINUE
584 nrun = nrun + 7 - k1
585*
586* --- Test SSYSVXX ---
587*
588* Restore the matrices A and B.
589*
590 IF( ifact.EQ.2 )
591 $ CALL slaset( uplo, n, n, zero, zero, afac, lda )
592 CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
593*
594* Solve the system and compute the condition number
595* and error bounds using SSYSVXX.
596*
597 srnamt = 'SSYSVXX'
598 n_err_bnds = 3
599 equed = 'N'
600 CALL ssysvxx( fact, uplo, n, nrhs, a, lda, afac,
601 $ lda, iwork, equed, work( n+1 ), b, lda, x,
602 $ lda, rcond, rpvgrw_svxx, berr, n_err_bnds,
603 $ errbnds_n, errbnds_c, 0, zero, work,
604 $ iwork( n+1 ), info )
605*
606* Adjust the expected value of INFO to account for
607* pivoting.
608*
609 k = izero
610 IF( k.GT.0 ) THEN
611 135 CONTINUE
612 IF( iwork( k ).LT.0 ) THEN
613 IF( iwork( k ).NE.-k ) THEN
614 k = -iwork( k )
615 GO TO 135
616 END IF
617 ELSE IF( iwork( k ).NE.k ) THEN
618 k = iwork( k )
619 GO TO 135
620 END IF
621 END IF
622*
623* Check the error code from SSYSVXX.
624*
625 IF( info.NE.k .AND. info.LE.n ) THEN
626 CALL alaerh( path, 'SSYSVXX', info, k,
627 $ fact // uplo, n, n, -1, -1, nrhs, imat, nfail,
628 $ nerrs, nout )
629 GO TO 150
630 END IF
631*
632 IF( info.EQ.0 ) THEN
633 IF( ifact.GE.2 ) THEN
634*
635* Reconstruct matrix from factors and compute
636* residual.
637*
638 CALL ssyt01( uplo, n, a, lda, afac, lda, iwork,
639 $ ainv, lda, rwork(2*nrhs+1),
640 $ result( 1 ) )
641 k1 = 1
642 ELSE
643 k1 = 2
644 END IF
645*
646* Compute residual of the computed solution.
647*
648 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
649 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
650 $ lda, rwork( 2*nrhs+1 ), result( 2 ) )
651*
652* Check solution from generated exact solution.
653*
654 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
655 $ result( 3 ) )
656*
657* Check the error bounds from iterative refinement.
658*
659 CALL spot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
660 $ xact, lda, rwork, rwork( nrhs+1 ),
661 $ result( 4 ) )
662 ELSE
663 k1 = 6
664 END IF
665*
666* Compare RCOND from SSYSVXX with the computed value
667* in RCONDC.
668*
669 result( 6 ) = sget06( rcond, rcondc )
670*
671* Print information about the tests that did not pass
672* the threshold.
673*
674 DO 85 k = k1, 6
675 IF( result( k ).GE.thresh ) THEN
676 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
677 $ CALL aladhd( nout, path )
678 WRITE( nout, fmt = 9998 )'SSYSVXX',
679 $ fact, uplo, n, imat, k,
680 $ result( k )
681 nfail = nfail + 1
682 END IF
683 85 CONTINUE
684 nrun = nrun + 7 - k1
685*
686 150 CONTINUE
687*
688 160 CONTINUE
689 170 CONTINUE
690 180 CONTINUE
691*
692* Print a summary of the results.
693*
694 CALL alasvm( path, nout, nfail, nrun, nerrs )
695*
696
697* Test Error Bounds from SSYSVXX
698
699 CALL sebchvxx(thresh, path)
700
701 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
702 $ ', test ', i2, ', ratio =', g12.5 )
703 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N =', i5,
704 $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
705 RETURN
706*
707* End of SDRVSYX
708*
709 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.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 aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine ssysv(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
SSYSV computes the solution to system of linear equations A * X = B for SY matrices
Definition ssysv.f:171
subroutine ssysvx(fact, uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, lwork, iwork, info)
SSYSVX computes the solution to system of linear equations A * X = B for SY matrices
Definition ssysvx.f:284
subroutine ssysvxx(fact, uplo, n, nrhs, a, lda, af, ldaf, ipiv, equed, s, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
SSYSVXX
Definition ssysvxx.f:508
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 slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sdrvsy(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
SDRVSY
Definition sdrvsy.f:152
subroutine sebchvxx(thresh, path)
SEBCHVXX
Definition sebchvxx.f:96
subroutine serrvx(path, nunit)
SERRVX
Definition serrvx.f:55
subroutine sget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
SGET04
Definition sget04.f:102
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 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