LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvsy.f
Go to the documentation of this file.
1*> \brief \b SDRVSY
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 and -SVX.
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] NRHS
60*> \verbatim
61*> NRHS is INTEGER
62*> The number of right hand side vectors to be generated for
63*> each linear system.
64*> \endverbatim
65*>
66*> \param[in] THRESH
67*> \verbatim
68*> THRESH is REAL
69*> The threshold value for the test ratios. A result is
70*> included in the output file if RESULT >= THRESH. To have
71*> every test ratio printed, use THRESH = 0.
72*> \endverbatim
73*>
74*> \param[in] TSTERR
75*> \verbatim
76*> TSTERR is LOGICAL
77*> Flag that indicates whether error exits are to be tested.
78*> \endverbatim
79*>
80*> \param[in] NMAX
81*> \verbatim
82*> NMAX is INTEGER
83*> The maximum value permitted for N, used in dimensioning the
84*> work arrays.
85*> \endverbatim
86*>
87*> \param[out] A
88*> \verbatim
89*> A is REAL array, dimension (NMAX*NMAX)
90*> \endverbatim
91*>
92*> \param[out] AFAC
93*> \verbatim
94*> AFAC is REAL array, dimension (NMAX*NMAX)
95*> \endverbatim
96*>
97*> \param[out] AINV
98*> \verbatim
99*> AINV is REAL array, dimension (NMAX*NMAX)
100*> \endverbatim
101*>
102*> \param[out] B
103*> \verbatim
104*> B is REAL array, dimension (NMAX*NRHS)
105*> \endverbatim
106*>
107*> \param[out] X
108*> \verbatim
109*> X is REAL array, dimension (NMAX*NRHS)
110*> \endverbatim
111*>
112*> \param[out] XACT
113*> \verbatim
114*> XACT is REAL array, dimension (NMAX*NRHS)
115*> \endverbatim
116*>
117*> \param[out] WORK
118*> \verbatim
119*> WORK is REAL array, dimension (NMAX*max(2,NRHS))
120*> \endverbatim
121*>
122*> \param[out] RWORK
123*> \verbatim
124*> RWORK is REAL array, dimension (NMAX+2*NRHS)
125*> \endverbatim
126*>
127*> \param[out] IWORK
128*> \verbatim
129*> IWORK is INTEGER array, dimension (2*NMAX)
130*> \endverbatim
131*>
132*> \param[in] NOUT
133*> \verbatim
134*> NOUT is INTEGER
135*> The unit number for output.
136*> \endverbatim
137*
138* Authors:
139* ========
140*
141*> \author Univ. of Tennessee
142*> \author Univ. of California Berkeley
143*> \author Univ. of Colorado Denver
144*> \author NAG Ltd.
145*
146*> \ingroup single_lin
147*
148* =====================================================================
149 SUBROUTINE sdrvsy( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
150 $ A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK,
151 $ NOUT )
152*
153* -- LAPACK test routine --
154* -- LAPACK is a software package provided by Univ. of Tennessee, --
155* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
156*
157* .. Scalar Arguments ..
158 LOGICAL TSTERR
159 INTEGER NMAX, NN, NOUT, NRHS
160 REAL THRESH
161* ..
162* .. Array Arguments ..
163 LOGICAL DOTYPE( * )
164 INTEGER IWORK( * ), NVAL( * )
165 REAL A( * ), AFAC( * ), AINV( * ), B( * ),
166 $ rwork( * ), work( * ), x( * ), xact( * )
167* ..
168*
169* =====================================================================
170*
171* .. Parameters ..
172 REAL ONE, ZERO
173 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
174 INTEGER NTYPES, NTESTS
175 parameter( ntypes = 10, ntests = 6 )
176 INTEGER NFACT
177 parameter( nfact = 2 )
178* ..
179* .. Local Scalars ..
180 LOGICAL ZEROT
181 CHARACTER DIST, FACT, TYPE, UPLO, XTYPE
182 CHARACTER*3 PATH
183 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
184 $ izero, j, k, k1, kl, ku, lda, lwork, mode, n,
185 $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
186 REAL AINVNM, ANORM, CNDNUM, RCOND, RCONDC
187* ..
188* .. Local Arrays ..
189 CHARACTER FACTS( NFACT ), UPLOS( 2 )
190 INTEGER ISEED( 4 ), ISEEDY( 4 )
191 REAL RESULT( NTESTS )
192* ..
193* .. External Functions ..
194 REAL SGET06, SLANSY
195 EXTERNAL SGET06, SLANSY
196* ..
197* .. External Subroutines ..
198 EXTERNAL aladhd, alaerh, alasvm, serrvx, sget04, slacpy,
201* ..
202* .. Scalars in Common ..
203 LOGICAL LERR, OK
204 CHARACTER*32 SRNAMT
205 INTEGER INFOT, NUNIT
206* ..
207* .. Common blocks ..
208 COMMON / infoc / infot, nunit, ok, lerr
209 COMMON / srnamc / srnamt
210* ..
211* .. Intrinsic Functions ..
212 INTRINSIC max, min
213* ..
214* .. Data statements ..
215 DATA iseedy / 1988, 1989, 1990, 1991 /
216 DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
217* ..
218* .. Executable Statements ..
219*
220* Initialize constants and the random number seed.
221*
222 path( 1: 1 ) = 'Single precision'
223 path( 2: 3 ) = 'SY'
224 nrun = 0
225 nfail = 0
226 nerrs = 0
227 DO 10 i = 1, 4
228 iseed( i ) = iseedy( i )
229 10 CONTINUE
230 lwork = max( 2*nmax, nmax*nrhs )
231*
232* Test the error exits
233*
234 IF( tsterr )
235 $ CALL serrvx( path, nout )
236 infot = 0
237*
238* Set the block size and minimum block size for testing.
239*
240 nb = 1
241 nbmin = 2
242 CALL xlaenv( 1, nb )
243 CALL xlaenv( 2, nbmin )
244*
245* Do for each value of N in NVAL
246*
247 DO 180 in = 1, nn
248 n = nval( in )
249 lda = max( n, 1 )
250 xtype = 'N'
251 nimat = ntypes
252 IF( n.LE.0 )
253 $ nimat = 1
254*
255 DO 170 imat = 1, nimat
256*
257* Do the tests only if DOTYPE( IMAT ) is true.
258*
259 IF( .NOT.dotype( imat ) )
260 $ GO TO 170
261*
262* Skip types 3, 4, 5, or 6 if the matrix size is too small.
263*
264 zerot = imat.GE.3 .AND. imat.LE.6
265 IF( zerot .AND. n.LT.imat-2 )
266 $ GO TO 170
267*
268* Do first for UPLO = 'U', then for UPLO = 'L'
269*
270 DO 160 iuplo = 1, 2
271 uplo = uplos( iuplo )
272*
273* Set up parameters with SLATB4 and generate a test matrix
274* with SLATMS.
275*
276 CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
277 $ cndnum, dist )
278*
279 srnamt = 'SLATMS'
280 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
281 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
282 $ info )
283*
284* Check error code from SLATMS.
285*
286 IF( info.NE.0 ) THEN
287 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
288 $ -1, -1, imat, nfail, nerrs, nout )
289 GO TO 160
290 END IF
291*
292* For types 3-6, zero one or more rows and columns of the
293* matrix to test that INFO is returned correctly.
294*
295 IF( zerot ) THEN
296 IF( imat.EQ.3 ) THEN
297 izero = 1
298 ELSE IF( imat.EQ.4 ) THEN
299 izero = n
300 ELSE
301 izero = n / 2 + 1
302 END IF
303*
304 IF( imat.LT.6 ) THEN
305*
306* Set row and column IZERO to zero.
307*
308 IF( iuplo.EQ.1 ) THEN
309 ioff = ( izero-1 )*lda
310 DO 20 i = 1, izero - 1
311 a( ioff+i ) = zero
312 20 CONTINUE
313 ioff = ioff + izero
314 DO 30 i = izero, n
315 a( ioff ) = zero
316 ioff = ioff + lda
317 30 CONTINUE
318 ELSE
319 ioff = izero
320 DO 40 i = 1, izero - 1
321 a( ioff ) = zero
322 ioff = ioff + lda
323 40 CONTINUE
324 ioff = ioff - izero
325 DO 50 i = izero, n
326 a( ioff+i ) = zero
327 50 CONTINUE
328 END IF
329 ELSE
330 ioff = 0
331 IF( iuplo.EQ.1 ) THEN
332*
333* Set the first IZERO rows and columns to zero.
334*
335 DO 70 j = 1, n
336 i2 = min( j, izero )
337 DO 60 i = 1, i2
338 a( ioff+i ) = zero
339 60 CONTINUE
340 ioff = ioff + lda
341 70 CONTINUE
342 ELSE
343*
344* Set the last IZERO rows and columns to zero.
345*
346 DO 90 j = 1, n
347 i1 = max( j, izero )
348 DO 80 i = i1, n
349 a( ioff+i ) = zero
350 80 CONTINUE
351 ioff = ioff + lda
352 90 CONTINUE
353 END IF
354 END IF
355 ELSE
356 izero = 0
357 END IF
358*
359 DO 150 ifact = 1, nfact
360*
361* Do first for FACT = 'F', then for other values.
362*
363 fact = facts( ifact )
364*
365* Compute the condition number for comparison with
366* the value returned by SSYSVX.
367*
368 IF( zerot ) THEN
369 IF( ifact.EQ.1 )
370 $ GO TO 150
371 rcondc = zero
372*
373 ELSE IF( ifact.EQ.1 ) THEN
374*
375* Compute the 1-norm of A.
376*
377 anorm = slansy( '1', uplo, n, a, lda, rwork )
378*
379* Factor the matrix A.
380*
381 CALL slacpy( uplo, n, n, a, lda, afac, lda )
382 CALL ssytrf( uplo, n, afac, lda, iwork, work,
383 $ lwork, info )
384*
385* Compute inv(A) and take its norm.
386*
387 CALL slacpy( uplo, n, n, afac, lda, ainv, lda )
388 lwork = (n+nb+1)*(nb+3)
389 CALL ssytri2( uplo, n, ainv, lda, iwork, work,
390 $ lwork, info )
391 ainvnm = slansy( '1', uplo, n, ainv, lda, rwork )
392*
393* Compute the 1-norm condition number of A.
394*
395 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
396 rcondc = one
397 ELSE
398 rcondc = ( one / anorm ) / ainvnm
399 END IF
400 END IF
401*
402* Form an exact solution and set the right hand side.
403*
404 srnamt = 'SLARHS'
405 CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
406 $ nrhs, a, lda, xact, lda, b, lda, iseed,
407 $ info )
408 xtype = 'C'
409*
410* --- Test SSYSV ---
411*
412 IF( ifact.EQ.2 ) THEN
413 CALL slacpy( uplo, n, n, a, lda, afac, lda )
414 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
415*
416* Factor the matrix and solve the system using SSYSV.
417*
418 srnamt = 'SSYSV '
419 CALL ssysv( uplo, n, nrhs, afac, lda, iwork, x,
420 $ lda, work, lwork, info )
421*
422* Adjust the expected value of INFO to account for
423* pivoting.
424*
425 k = izero
426 IF( k.GT.0 ) THEN
427 100 CONTINUE
428 IF( iwork( k ).LT.0 ) THEN
429 IF( iwork( k ).NE.-k ) THEN
430 k = -iwork( k )
431 GO TO 100
432 END IF
433 ELSE IF( iwork( k ).NE.k ) THEN
434 k = iwork( k )
435 GO TO 100
436 END IF
437 END IF
438*
439* Check error code from SSYSV .
440*
441 IF( info.NE.k ) THEN
442 CALL alaerh( path, 'SSYSV ', info, k, uplo, n,
443 $ n, -1, -1, nrhs, imat, nfail,
444 $ nerrs, nout )
445 GO TO 120
446 ELSE IF( info.NE.0 ) THEN
447 GO TO 120
448 END IF
449*
450* Reconstruct matrix from factors and compute
451* residual.
452*
453 CALL ssyt01( uplo, n, a, lda, afac, lda, iwork,
454 $ ainv, lda, rwork, result( 1 ) )
455*
456* Compute residual of the computed solution.
457*
458 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
459 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
460 $ lda, rwork, result( 2 ) )
461*
462* Check solution from generated exact solution.
463*
464 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
465 $ result( 3 ) )
466 nt = 3
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 aladhd( nout, path )
475 WRITE( nout, fmt = 9999 )'SSYSV ', uplo, n,
476 $ imat, k, result( k )
477 nfail = nfail + 1
478 END IF
479 110 CONTINUE
480 nrun = nrun + nt
481 120 CONTINUE
482 END IF
483*
484* --- Test SSYSVX ---
485*
486 IF( ifact.EQ.2 )
487 $ CALL slaset( uplo, n, n, zero, zero, afac, lda )
488 CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
489*
490* Solve the system and compute the condition number and
491* error bounds using SSYSVX.
492*
493 srnamt = 'SSYSVX'
494 CALL ssysvx( fact, uplo, n, nrhs, a, lda, afac, lda,
495 $ iwork, b, lda, x, lda, rcond, rwork,
496 $ rwork( nrhs+1 ), work, lwork,
497 $ iwork( n+1 ), info )
498*
499* Adjust the expected value of INFO to account for
500* pivoting.
501*
502 k = izero
503 IF( k.GT.0 ) THEN
504 130 CONTINUE
505 IF( iwork( k ).LT.0 ) THEN
506 IF( iwork( k ).NE.-k ) THEN
507 k = -iwork( k )
508 GO TO 130
509 END IF
510 ELSE IF( iwork( k ).NE.k ) THEN
511 k = iwork( k )
512 GO TO 130
513 END IF
514 END IF
515*
516* Check the error code from SSYSVX.
517*
518 IF( info.NE.k ) THEN
519 CALL alaerh( path, 'SSYSVX', info, k, fact // uplo,
520 $ n, n, -1, -1, nrhs, imat, nfail,
521 $ nerrs, nout )
522 GO TO 150
523 END IF
524*
525 IF( info.EQ.0 ) THEN
526 IF( ifact.GE.2 ) THEN
527*
528* Reconstruct matrix from factors and compute
529* residual.
530*
531 CALL ssyt01( uplo, n, a, lda, afac, lda, iwork,
532 $ ainv, lda, rwork( 2*nrhs+1 ),
533 $ result( 1 ) )
534 k1 = 1
535 ELSE
536 k1 = 2
537 END IF
538*
539* Compute residual of the computed solution.
540*
541 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
542 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
543 $ lda, rwork( 2*nrhs+1 ), result( 2 ) )
544*
545* Check solution from generated exact solution.
546*
547 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
548 $ result( 3 ) )
549*
550* Check the error bounds from iterative refinement.
551*
552 CALL spot05( uplo, n, nrhs, a, lda, b, lda, x, lda,
553 $ xact, lda, rwork, rwork( nrhs+1 ),
554 $ result( 4 ) )
555 ELSE
556 k1 = 6
557 END IF
558*
559* Compare RCOND from SSYSVX with the computed value
560* in RCONDC.
561*
562 result( 6 ) = sget06( rcond, rcondc )
563*
564* Print information about the tests that did not pass
565* the threshold.
566*
567 DO 140 k = k1, 6
568 IF( result( k ).GE.thresh ) THEN
569 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
570 $ CALL aladhd( nout, path )
571 WRITE( nout, fmt = 9998 )'SSYSVX', fact, uplo,
572 $ n, imat, k, result( k )
573 nfail = nfail + 1
574 END IF
575 140 CONTINUE
576 nrun = nrun + 7 - k1
577*
578 150 CONTINUE
579*
580 160 CONTINUE
581 170 CONTINUE
582 180 CONTINUE
583*
584* Print a summary of the results.
585*
586 CALL alasvm( path, nout, nfail, nrun, nerrs )
587*
588 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
589 $ ', test ', i2, ', ratio =', g12.5 )
590 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N =', i5,
591 $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
592 RETURN
593*
594* End of SDRVSY
595*
596 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 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 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