LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvsp.f
Go to the documentation of this file.
1*> \brief \b SDRVSP
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 SDRVSP( 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*> SDRVSP tests the driver routines SSPSV 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
90*> (NMAX*(NMAX+1)/2)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*> AFAC is REAL array, dimension
96*> (NMAX*(NMAX+1)/2)
97*> \endverbatim
98*>
99*> \param[out] AINV
100*> \verbatim
101*> AINV is REAL array, dimension
102*> (NMAX*(NMAX+1)/2)
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 sdrvsp( 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, FACT, PACKIT, 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 $ nerrs, nfail, nimat, npp, nrun, nt
190 REAL AINVNM, ANORM, CNDNUM, RCOND, RCONDC
191* ..
192* .. Local Arrays ..
193 CHARACTER FACTS( NFACT )
194 INTEGER ISEED( 4 ), ISEEDY( 4 )
195 REAL RESULT( NTESTS )
196* ..
197* .. External Functions ..
198 REAL SGET06, SLANSP
199 EXTERNAL SGET06, SLANSP
200* ..
201* .. External Subroutines ..
202 EXTERNAL aladhd, alaerh, alasvm, scopy, serrvx, sget04,
205* ..
206* .. Scalars in Common ..
207 LOGICAL LERR, OK
208 CHARACTER*32 SRNAMT
209 INTEGER INFOT, NUNIT
210* ..
211* .. Common blocks ..
212 COMMON / infoc / infot, nunit, ok, lerr
213 COMMON / srnamc / srnamt
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC max, min
217* ..
218* .. Data statements ..
219 DATA iseedy / 1988, 1989, 1990, 1991 /
220 DATA facts / 'F', 'N' /
221* ..
222* .. Executable Statements ..
223*
224* Initialize constants and the random number seed.
225*
226 path( 1: 1 ) = 'Single precision'
227 path( 2: 3 ) = 'SP'
228 nrun = 0
229 nfail = 0
230 nerrs = 0
231 DO 10 i = 1, 4
232 iseed( i ) = iseedy( i )
233 10 CONTINUE
234 lwork = max( 2*nmax, nmax*nrhs )
235*
236* Test the error exits
237*
238 IF( tsterr )
239 $ CALL serrvx( path, nout )
240 infot = 0
241*
242* Do for each value of N in NVAL
243*
244 DO 180 in = 1, nn
245 n = nval( in )
246 lda = max( n, 1 )
247 npp = n*( n+1 ) / 2
248 xtype = 'N'
249 nimat = ntypes
250 IF( n.LE.0 )
251 $ nimat = 1
252*
253 DO 170 imat = 1, nimat
254*
255* Do the tests only if DOTYPE( IMAT ) is true.
256*
257 IF( .NOT.dotype( imat ) )
258 $ GO TO 170
259*
260* Skip types 3, 4, 5, or 6 if the matrix size is too small.
261*
262 zerot = imat.GE.3 .AND. imat.LE.6
263 IF( zerot .AND. n.LT.imat-2 )
264 $ GO TO 170
265*
266* Do first for UPLO = 'U', then for UPLO = 'L'
267*
268 DO 160 iuplo = 1, 2
269 IF( iuplo.EQ.1 ) THEN
270 uplo = 'U'
271 packit = 'C'
272 ELSE
273 uplo = 'L'
274 packit = 'R'
275 END IF
276*
277* Set up parameters with SLATB4 and generate a test matrix
278* with SLATMS.
279*
280 CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
281 $ cndnum, dist )
282*
283 srnamt = 'SLATMS'
284 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
285 $ cndnum, anorm, kl, ku, packit, a, lda, work,
286 $ info )
287*
288* Check error code from SLATMS.
289*
290 IF( info.NE.0 ) THEN
291 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
292 $ -1, -1, imat, nfail, nerrs, nout )
293 GO TO 160
294 END IF
295*
296* For types 3-6, zero one or more rows and columns of the
297* matrix to test that INFO is returned correctly.
298*
299 IF( zerot ) THEN
300 IF( imat.EQ.3 ) THEN
301 izero = 1
302 ELSE IF( imat.EQ.4 ) THEN
303 izero = n
304 ELSE
305 izero = n / 2 + 1
306 END IF
307*
308 IF( imat.LT.6 ) THEN
309*
310* Set row and column IZERO to zero.
311*
312 IF( iuplo.EQ.1 ) THEN
313 ioff = ( izero-1 )*izero / 2
314 DO 20 i = 1, izero - 1
315 a( ioff+i ) = zero
316 20 CONTINUE
317 ioff = ioff + izero
318 DO 30 i = izero, n
319 a( ioff ) = zero
320 ioff = ioff + i
321 30 CONTINUE
322 ELSE
323 ioff = izero
324 DO 40 i = 1, izero - 1
325 a( ioff ) = zero
326 ioff = ioff + n - i
327 40 CONTINUE
328 ioff = ioff - izero
329 DO 50 i = izero, n
330 a( ioff+i ) = zero
331 50 CONTINUE
332 END IF
333 ELSE
334 ioff = 0
335 IF( iuplo.EQ.1 ) THEN
336*
337* Set the first IZERO rows and columns to zero.
338*
339 DO 70 j = 1, n
340 i2 = min( j, izero )
341 DO 60 i = 1, i2
342 a( ioff+i ) = zero
343 60 CONTINUE
344 ioff = ioff + j
345 70 CONTINUE
346 ELSE
347*
348* Set the last IZERO rows and columns to zero.
349*
350 DO 90 j = 1, n
351 i1 = max( j, izero )
352 DO 80 i = i1, n
353 a( ioff+i ) = zero
354 80 CONTINUE
355 ioff = ioff + n - j
356 90 CONTINUE
357 END IF
358 END IF
359 ELSE
360 izero = 0
361 END IF
362*
363 DO 150 ifact = 1, nfact
364*
365* Do first for FACT = 'F', then for other values.
366*
367 fact = facts( ifact )
368*
369* Compute the condition number for comparison with
370* the value returned by SSPSVX.
371*
372 IF( zerot ) THEN
373 IF( ifact.EQ.1 )
374 $ GO TO 150
375 rcondc = zero
376*
377 ELSE IF( ifact.EQ.1 ) THEN
378*
379* Compute the 1-norm of A.
380*
381 anorm = slansp( '1', uplo, n, a, rwork )
382*
383* Factor the matrix A.
384*
385 CALL scopy( npp, a, 1, afac, 1 )
386 CALL ssptrf( uplo, n, afac, iwork, info )
387*
388* Compute inv(A) and take its norm.
389*
390 CALL scopy( npp, afac, 1, ainv, 1 )
391 CALL ssptri( uplo, n, ainv, iwork, work, info )
392 ainvnm = slansp( '1', uplo, n, ainv, rwork )
393*
394* Compute the 1-norm condition number of A.
395*
396 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
397 rcondc = one
398 ELSE
399 rcondc = ( one / anorm ) / ainvnm
400 END IF
401 END IF
402*
403* Form an exact solution and set the right hand side.
404*
405 srnamt = 'SLARHS'
406 CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
407 $ nrhs, a, lda, xact, lda, b, lda, iseed,
408 $ info )
409 xtype = 'C'
410*
411* --- Test SSPSV ---
412*
413 IF( ifact.EQ.2 ) THEN
414 CALL scopy( npp, a, 1, afac, 1 )
415 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
416*
417* Factor the matrix and solve the system using SSPSV.
418*
419 srnamt = 'SSPSV '
420 CALL sspsv( uplo, n, nrhs, afac, iwork, x, lda,
421 $ info )
422*
423* Adjust the expected value of INFO to account for
424* pivoting.
425*
426 k = izero
427 IF( k.GT.0 ) THEN
428 100 CONTINUE
429 IF( iwork( k ).LT.0 ) THEN
430 IF( iwork( k ).NE.-k ) THEN
431 k = -iwork( k )
432 GO TO 100
433 END IF
434 ELSE IF( iwork( k ).NE.k ) THEN
435 k = iwork( k )
436 GO TO 100
437 END IF
438 END IF
439*
440* Check error code from SSPSV .
441*
442 IF( info.NE.k ) THEN
443 CALL alaerh( path, 'SSPSV ', info, k, uplo, n,
444 $ n, -1, -1, nrhs, imat, nfail,
445 $ nerrs, nout )
446 GO TO 120
447 ELSE IF( info.NE.0 ) THEN
448 GO TO 120
449 END IF
450*
451* Reconstruct matrix from factors and compute
452* residual.
453*
454 CALL sspt01( uplo, n, a, afac, iwork, ainv, lda,
455 $ rwork, result( 1 ) )
456*
457* Compute residual of the computed solution.
458*
459 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
460 CALL sppt02( uplo, n, nrhs, a, x, lda, work, lda,
461 $ rwork, result( 2 ) )
462*
463* Check solution from generated exact solution.
464*
465 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
466 $ result( 3 ) )
467 nt = 3
468*
469* Print information about the tests that did not pass
470* the threshold.
471*
472 DO 110 k = 1, nt
473 IF( result( k ).GE.thresh ) THEN
474 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
475 $ CALL aladhd( nout, path )
476 WRITE( nout, fmt = 9999 )'SSPSV ', uplo, n,
477 $ imat, k, result( k )
478 nfail = nfail + 1
479 END IF
480 110 CONTINUE
481 nrun = nrun + nt
482 120 CONTINUE
483 END IF
484*
485* --- Test SSPSVX ---
486*
487 IF( ifact.EQ.2 .AND. npp.GT.0 )
488 $ CALL slaset( 'Full', npp, 1, zero, zero, afac,
489 $ npp )
490 CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
491*
492* Solve the system and compute the condition number and
493* error bounds using SSPSVX.
494*
495 srnamt = 'SSPSVX'
496 CALL sspsvx( fact, uplo, n, nrhs, a, afac, iwork, b,
497 $ lda, x, lda, rcond, rwork,
498 $ rwork( nrhs+1 ), work, iwork( n+1 ),
499 $ info )
500*
501* Adjust the expected value of INFO to account for
502* pivoting.
503*
504 k = izero
505 IF( k.GT.0 ) THEN
506 130 CONTINUE
507 IF( iwork( k ).LT.0 ) THEN
508 IF( iwork( k ).NE.-k ) THEN
509 k = -iwork( k )
510 GO TO 130
511 END IF
512 ELSE IF( iwork( k ).NE.k ) THEN
513 k = iwork( k )
514 GO TO 130
515 END IF
516 END IF
517*
518* Check the error code from SSPSVX.
519*
520 IF( info.NE.k ) THEN
521 CALL alaerh( path, 'SSPSVX', info, k, fact // uplo,
522 $ n, n, -1, -1, nrhs, imat, nfail,
523 $ nerrs, nout )
524 GO TO 150
525 END IF
526*
527 IF( info.EQ.0 ) THEN
528 IF( ifact.GE.2 ) THEN
529*
530* Reconstruct matrix from factors and compute
531* residual.
532*
533 CALL sspt01( uplo, n, a, afac, iwork, ainv, lda,
534 $ rwork( 2*nrhs+1 ), result( 1 ) )
535 k1 = 1
536 ELSE
537 k1 = 2
538 END IF
539*
540* Compute residual of the computed solution.
541*
542 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
543 CALL sppt02( uplo, n, nrhs, a, x, lda, work, lda,
544 $ rwork( 2*nrhs+1 ), result( 2 ) )
545*
546* Check solution from generated exact solution.
547*
548 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
549 $ result( 3 ) )
550*
551* Check the error bounds from iterative refinement.
552*
553 CALL sppt05( uplo, n, nrhs, a, b, lda, x, lda,
554 $ xact, lda, rwork, rwork( nrhs+1 ),
555 $ result( 4 ) )
556 ELSE
557 k1 = 6
558 END IF
559*
560* Compare RCOND from SSPSVX with the computed value
561* in RCONDC.
562*
563 result( 6 ) = sget06( rcond, rcondc )
564*
565* Print information about the tests that did not pass
566* the threshold.
567*
568 DO 140 k = k1, 6
569 IF( result( k ).GE.thresh ) THEN
570 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
571 $ CALL aladhd( nout, path )
572 WRITE( nout, fmt = 9998 )'SSPSVX', fact, uplo,
573 $ n, imat, k, result( k )
574 nfail = nfail + 1
575 END IF
576 140 CONTINUE
577 nrun = nrun + 7 - k1
578*
579 150 CONTINUE
580*
581 160 CONTINUE
582 170 CONTINUE
583 180 CONTINUE
584*
585* Print a summary of the results.
586*
587 CALL alasvm( path, nout, nfail, nrun, nerrs )
588*
589 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
590 $ ', test ', i2, ', ratio =', g12.5 )
591 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N =', i5,
592 $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
593 RETURN
594*
595* End of SDRVSP
596*
597 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 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 scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sspsv(uplo, n, nrhs, ap, ipiv, b, ldb, info)
SSPSV computes the solution to system of linear equations A * X = B for OTHER matrices
Definition sspsv.f:162
subroutine sspsvx(fact, uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
SSPSVX computes the solution to system of linear equations A * X = B for OTHER matrices
Definition sspsvx.f:276
subroutine ssptrf(uplo, n, ap, ipiv, info)
SSPTRF
Definition ssptrf.f:157
subroutine ssptri(uplo, n, ap, ipiv, work, info)
SSPTRI
Definition ssptri.f:109
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 sdrvsp(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
SDRVSP
Definition sdrvsp.f:156
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 sppt02(uplo, n, nrhs, a, x, ldx, b, ldb, rwork, resid)
SPPT02
Definition sppt02.f:122
subroutine sppt05(uplo, n, nrhs, ap, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
SPPT05
Definition sppt05.f:156
subroutine sspt01(uplo, n, a, afac, ipiv, c, ldc, rwork, resid)
SSPT01
Definition sspt01.f:110