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