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