LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 *> \date November 2013
162 *
163 *> \ingroup double_lin
164 *
165 * =====================================================================
166  SUBROUTINE ddrvpo( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
167  $ a, afac, asav, b, bsav, x, xact, s, work,
168  $ rwork, iwork, nout )
169 *
170 * -- LAPACK test routine (version 3.5.0) --
171 * -- LAPACK is a software package provided by Univ. of Tennessee, --
172 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
173 * November 2013
174 *
175 * .. Scalar Arguments ..
176  LOGICAL tsterr
177  INTEGER nmax, nn, nout, nrhs
178  DOUBLE PRECISION thresh
179 * ..
180 * .. Array Arguments ..
181  LOGICAL dotype( * )
182  INTEGER iwork( * ), nval( * )
183  DOUBLE PRECISION a( * ), afac( * ), asav( * ), b( * ),
184  $ bsav( * ), rwork( * ), s( * ), work( * ),
185  $ x( * ), xact( * )
186 * ..
187 *
188 * =====================================================================
189 *
190 * .. Parameters ..
191  DOUBLE PRECISION one, zero
192  parameter ( one = 1.0d+0, zero = 0.0d+0 )
193  INTEGER ntypes
194  parameter ( ntypes = 9 )
195  INTEGER ntests
196  parameter ( ntests = 6 )
197 * ..
198 * .. Local Scalars ..
199  LOGICAL equil, nofact, prefac, zerot
200  CHARACTER dist, equed, fact, TYPE, uplo, xtype
201  CHARACTER*3 path
202  INTEGER i, iequed, ifact, imat, in, info, ioff, iuplo,
203  $ izero, k, k1, kl, ku, lda, mode, n, nb, nbmin,
204  $ nerrs, nfact, nfail, nimat, nrun, nt,
205  $ n_err_bnds
206  DOUBLE PRECISION ainvnm, amax, anorm, cndnum, rcond, rcondc,
207  $ roldc, scond, rpvgrw_svxx
208 * ..
209 * .. Local Arrays ..
210  CHARACTER equeds( 2 ), facts( 3 ), uplos( 2 )
211  INTEGER iseed( 4 ), iseedy( 4 )
212  DOUBLE PRECISION result( ntests ), berr( nrhs ),
213  $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
214 * ..
215 * .. External Functions ..
216  LOGICAL lsame
217  DOUBLE PRECISION dget06, dlansy
218  EXTERNAL lsame, dget06, dlansy
219 * ..
220 * .. External Subroutines ..
221  EXTERNAL aladhd, alaerh, alasvm, derrvx, dget04, dlacpy,
224  $ dpotri, xlaenv
225 * ..
226 * .. Intrinsic Functions ..
227  INTRINSIC max
228 * ..
229 * .. Scalars in Common ..
230  LOGICAL lerr, ok
231  CHARACTER*32 srnamt
232  INTEGER infot, nunit
233 * ..
234 * .. Common blocks ..
235  COMMON / infoc / infot, nunit, ok, lerr
236  COMMON / srnamc / srnamt
237 * ..
238 * .. Data statements ..
239  DATA iseedy / 1988, 1989, 1990, 1991 /
240  DATA uplos / 'U', 'L' /
241  DATA facts / 'F', 'N', 'E' /
242  DATA equeds / 'N', 'Y' /
243 * ..
244 * .. Executable Statements ..
245 *
246 * Initialize constants and the random number seed.
247 *
248  path( 1: 1 ) = 'Double precision'
249  path( 2: 3 ) = 'PO'
250  nrun = 0
251  nfail = 0
252  nerrs = 0
253  DO 10 i = 1, 4
254  iseed( i ) = iseedy( i )
255  10 CONTINUE
256 *
257 * Test the error exits
258 *
259  IF( tsterr )
260  $ CALL derrvx( path, nout )
261  infot = 0
262 *
263 * Set the block size and minimum block size for testing.
264 *
265  nb = 1
266  nbmin = 2
267  CALL xlaenv( 1, nb )
268  CALL xlaenv( 2, nbmin )
269 *
270 * Do for each value of N in NVAL
271 *
272  DO 130 in = 1, nn
273  n = nval( in )
274  lda = max( n, 1 )
275  xtype = 'N'
276  nimat = ntypes
277  IF( n.LE.0 )
278  $ nimat = 1
279 *
280  DO 120 imat = 1, nimat
281 *
282 * Do the tests only if DOTYPE( IMAT ) is true.
283 *
284  IF( .NOT.dotype( imat ) )
285  $ GO TO 120
286 *
287 * Skip types 3, 4, or 5 if the matrix size is too small.
288 *
289  zerot = imat.GE.3 .AND. imat.LE.5
290  IF( zerot .AND. n.LT.imat-2 )
291  $ GO TO 120
292 *
293 * Do first for UPLO = 'U', then for UPLO = 'L'
294 *
295  DO 110 iuplo = 1, 2
296  uplo = uplos( iuplo )
297 *
298 * Set up parameters with DLATB4 and generate a test matrix
299 * with DLATMS.
300 *
301  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
302  $ cndnum, dist )
303 *
304  srnamt = 'DLATMS'
305  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
306  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
307  $ info )
308 *
309 * Check error code from DLATMS.
310 *
311  IF( info.NE.0 ) THEN
312  CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
313  $ -1, -1, imat, nfail, nerrs, nout )
314  GO TO 110
315  END IF
316 *
317 * For types 3-5, zero one row and column of the matrix to
318 * test that INFO is returned correctly.
319 *
320  IF( zerot ) THEN
321  IF( imat.EQ.3 ) THEN
322  izero = 1
323  ELSE IF( imat.EQ.4 ) THEN
324  izero = n
325  ELSE
326  izero = n / 2 + 1
327  END IF
328  ioff = ( izero-1 )*lda
329 *
330 * Set row and column IZERO of A to 0.
331 *
332  IF( iuplo.EQ.1 ) THEN
333  DO 20 i = 1, izero - 1
334  a( ioff+i ) = zero
335  20 CONTINUE
336  ioff = ioff + izero
337  DO 30 i = izero, n
338  a( ioff ) = zero
339  ioff = ioff + lda
340  30 CONTINUE
341  ELSE
342  ioff = izero
343  DO 40 i = 1, izero - 1
344  a( ioff ) = zero
345  ioff = ioff + lda
346  40 CONTINUE
347  ioff = ioff - izero
348  DO 50 i = izero, n
349  a( ioff+i ) = zero
350  50 CONTINUE
351  END IF
352  ELSE
353  izero = 0
354  END IF
355 *
356 * Save a copy of the matrix A in ASAV.
357 *
358  CALL dlacpy( uplo, n, n, a, lda, asav, lda )
359 *
360  DO 100 iequed = 1, 2
361  equed = equeds( iequed )
362  IF( iequed.EQ.1 ) THEN
363  nfact = 3
364  ELSE
365  nfact = 1
366  END IF
367 *
368  DO 90 ifact = 1, nfact
369  fact = facts( ifact )
370  prefac = lsame( fact, 'F' )
371  nofact = lsame( fact, 'N' )
372  equil = lsame( fact, 'E' )
373 *
374  IF( zerot ) THEN
375  IF( prefac )
376  $ GO TO 90
377  rcondc = zero
378 *
379  ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
380 *
381 * Compute the condition number for comparison with
382 * the value returned by DPOSVX (FACT = 'N' reuses
383 * the condition number from the previous iteration
384 * with FACT = 'F').
385 *
386  CALL dlacpy( uplo, n, n, asav, lda, afac, lda )
387  IF( equil .OR. iequed.GT.1 ) THEN
388 *
389 * Compute row and column scale factors to
390 * equilibrate the matrix A.
391 *
392  CALL dpoequ( n, afac, lda, s, scond, amax,
393  $ info )
394  IF( info.EQ.0 .AND. n.GT.0 ) THEN
395  IF( iequed.GT.1 )
396  $ scond = zero
397 *
398 * Equilibrate the matrix.
399 *
400  CALL dlaqsy( uplo, n, afac, lda, s, scond,
401  $ amax, equed )
402  END IF
403  END IF
404 *
405 * Save the condition number of the
406 * non-equilibrated system for use in DGET04.
407 *
408  IF( equil )
409  $ roldc = rcondc
410 *
411 * Compute the 1-norm of A.
412 *
413  anorm = dlansy( '1', uplo, n, afac, lda, rwork )
414 *
415 * Factor the matrix A.
416 *
417  CALL dpotrf( uplo, n, afac, lda, info )
418 *
419 * Form the inverse of A.
420 *
421  CALL dlacpy( uplo, n, n, afac, lda, a, lda )
422  CALL dpotri( uplo, n, a, lda, info )
423 *
424 * Compute the 1-norm condition number of A.
425 *
426  ainvnm = dlansy( '1', uplo, n, a, lda, rwork )
427  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
428  rcondc = one
429  ELSE
430  rcondc = ( one / anorm ) / ainvnm
431  END IF
432  END IF
433 *
434 * Restore the matrix A.
435 *
436  CALL dlacpy( uplo, n, n, asav, lda, a, lda )
437 *
438 * Form an exact solution and set the right hand side.
439 *
440  srnamt = 'DLARHS'
441  CALL dlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
442  $ nrhs, a, lda, xact, lda, b, lda,
443  $ iseed, info )
444  xtype = 'C'
445  CALL dlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
446 *
447  IF( nofact ) THEN
448 *
449 * --- Test DPOSV ---
450 *
451 * Compute the L*L' or U'*U factorization of the
452 * matrix and solve the system.
453 *
454  CALL dlacpy( uplo, n, n, a, lda, afac, lda )
455  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
456 *
457  srnamt = 'DPOSV '
458  CALL dposv( uplo, n, nrhs, afac, lda, x, lda,
459  $ info )
460 *
461 * Check error code from DPOSV .
462 *
463  IF( info.NE.izero ) THEN
464  CALL alaerh( path, 'DPOSV ', info, izero,
465  $ uplo, n, n, -1, -1, nrhs, imat,
466  $ nfail, nerrs, nout )
467  GO TO 70
468  ELSE IF( info.NE.0 ) THEN
469  GO TO 70
470  END IF
471 *
472 * Reconstruct matrix from factors and compute
473 * residual.
474 *
475  CALL dpot01( uplo, n, a, lda, afac, lda, rwork,
476  $ result( 1 ) )
477 *
478 * Compute residual of the computed solution.
479 *
480  CALL dlacpy( 'Full', n, nrhs, b, lda, work,
481  $ lda )
482  CALL dpot02( uplo, n, nrhs, a, lda, x, lda,
483  $ work, lda, rwork, result( 2 ) )
484 *
485 * Check solution from generated exact solution.
486 *
487  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
488  $ result( 3 ) )
489  nt = 3
490 *
491 * Print information about the tests that did not
492 * pass the threshold.
493 *
494  DO 60 k = 1, nt
495  IF( result( k ).GE.thresh ) THEN
496  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
497  $ CALL aladhd( nout, path )
498  WRITE( nout, fmt = 9999 )'DPOSV ', uplo,
499  $ n, imat, k, result( k )
500  nfail = nfail + 1
501  END IF
502  60 CONTINUE
503  nrun = nrun + nt
504  70 CONTINUE
505  END IF
506 *
507 * --- Test DPOSVX ---
508 *
509  IF( .NOT.prefac )
510  $ CALL dlaset( uplo, n, n, zero, zero, afac, lda )
511  CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
512  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
513 *
514 * Equilibrate the matrix if FACT='F' and
515 * EQUED='Y'.
516 *
517  CALL dlaqsy( uplo, n, a, lda, s, scond, amax,
518  $ equed )
519  END IF
520 *
521 * Solve the system and compute the condition number
522 * and error bounds using DPOSVX.
523 *
524  srnamt = 'DPOSVX'
525  CALL dposvx( fact, uplo, n, nrhs, a, lda, afac,
526  $ lda, equed, s, b, lda, x, lda, rcond,
527  $ rwork, rwork( nrhs+1 ), work, iwork,
528  $ info )
529 *
530 * Check the error code from DPOSVX.
531 *
532  IF( info.NE.izero ) THEN
533  CALL alaerh( path, 'DPOSVX', info, izero,
534  $ fact // uplo, n, n, -1, -1, nrhs,
535  $ imat, nfail, nerrs, nout )
536  GO TO 90
537  END IF
538 *
539  IF( info.EQ.0 ) THEN
540  IF( .NOT.prefac ) THEN
541 *
542 * Reconstruct matrix from factors and compute
543 * residual.
544 *
545  CALL dpot01( uplo, n, a, lda, afac, lda,
546  $ rwork( 2*nrhs+1 ), result( 1 ) )
547  k1 = 1
548  ELSE
549  k1 = 2
550  END IF
551 *
552 * Compute residual of the computed solution.
553 *
554  CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
555  $ lda )
556  CALL dpot02( uplo, n, nrhs, asav, lda, x, lda,
557  $ work, lda, rwork( 2*nrhs+1 ),
558  $ result( 2 ) )
559 *
560 * Check solution from generated exact solution.
561 *
562  IF( nofact .OR. ( prefac .AND. lsame( equed,
563  $ 'N' ) ) ) THEN
564  CALL dget04( n, nrhs, x, lda, xact, lda,
565  $ rcondc, result( 3 ) )
566  ELSE
567  CALL dget04( n, nrhs, x, lda, xact, lda,
568  $ roldc, result( 3 ) )
569  END IF
570 *
571 * Check the error bounds from iterative
572 * refinement.
573 *
574  CALL dpot05( uplo, n, nrhs, asav, lda, b, lda,
575  $ x, lda, xact, lda, rwork,
576  $ rwork( nrhs+1 ), result( 4 ) )
577  ELSE
578  k1 = 6
579  END IF
580 *
581 * Compare RCOND from DPOSVX with the computed value
582 * in RCONDC.
583 *
584  result( 6 ) = dget06( rcond, rcondc )
585 *
586 * Print information about the tests that did not pass
587 * the threshold.
588 *
589  DO 80 k = k1, 6
590  IF( result( k ).GE.thresh ) THEN
591  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
592  $ CALL aladhd( nout, path )
593  IF( prefac ) THEN
594  WRITE( nout, fmt = 9997 )'DPOSVX', fact,
595  $ uplo, n, equed, imat, k, result( k )
596  ELSE
597  WRITE( nout, fmt = 9998 )'DPOSVX', fact,
598  $ uplo, n, imat, k, result( k )
599  END IF
600  nfail = nfail + 1
601  END IF
602  80 CONTINUE
603  nrun = nrun + 7 - k1
604 *
605 * --- Test DPOSVXX ---
606 *
607 * Restore the matrices A and B.
608 *
609  CALL dlacpy( 'Full', n, n, asav, lda, a, lda )
610  CALL dlacpy( 'Full', n, nrhs, bsav, lda, b, lda )
611 
612  IF( .NOT.prefac )
613  $ CALL dlaset( uplo, n, n, zero, zero, afac, lda )
614  CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
615  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
616 *
617 * Equilibrate the matrix if FACT='F' and
618 * EQUED='Y'.
619 *
620  CALL dlaqsy( uplo, n, a, lda, s, scond, amax,
621  $ equed )
622  END IF
623 *
624 * Solve the system and compute the condition number
625 * and error bounds using DPOSVXX.
626 *
627  srnamt = 'DPOSVXX'
628  n_err_bnds = 3
629  CALL dposvxx( fact, uplo, n, nrhs, a, lda, afac,
630  $ lda, equed, s, b, lda, x,
631  $ lda, rcond, rpvgrw_svxx, berr, n_err_bnds,
632  $ errbnds_n, errbnds_c, 0, zero, work,
633  $ iwork, info )
634 *
635 * Check the error code from DPOSVXX.
636 *
637  IF( info.EQ.n+1 ) GOTO 90
638  IF( info.NE.izero ) THEN
639  CALL alaerh( path, 'DPOSVXX', info, izero,
640  $ fact // uplo, n, n, -1, -1, nrhs,
641  $ imat, nfail, nerrs, nout )
642  GO TO 90
643  END IF
644 *
645  IF( info.EQ.0 ) THEN
646  IF( .NOT.prefac ) THEN
647 *
648 * Reconstruct matrix from factors and compute
649 * residual.
650 *
651  CALL dpot01( uplo, n, a, lda, afac, lda,
652  $ rwork( 2*nrhs+1 ), result( 1 ) )
653  k1 = 1
654  ELSE
655  k1 = 2
656  END IF
657 *
658 * Compute residual of the computed solution.
659 *
660  CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
661  $ lda )
662  CALL dpot02( uplo, n, nrhs, asav, lda, x, lda,
663  $ work, lda, rwork( 2*nrhs+1 ),
664  $ result( 2 ) )
665 *
666 * Check solution from generated exact solution.
667 *
668  IF( nofact .OR. ( prefac .AND. lsame( equed,
669  $ 'N' ) ) ) THEN
670  CALL dget04( n, nrhs, x, lda, xact, lda,
671  $ rcondc, result( 3 ) )
672  ELSE
673  CALL dget04( n, nrhs, x, lda, xact, lda,
674  $ roldc, result( 3 ) )
675  END IF
676 *
677 * Check the error bounds from iterative
678 * refinement.
679 *
680  CALL dpot05( uplo, n, nrhs, asav, lda, b, lda,
681  $ x, lda, xact, lda, rwork,
682  $ rwork( nrhs+1 ), result( 4 ) )
683  ELSE
684  k1 = 6
685  END IF
686 *
687 * Compare RCOND from DPOSVXX with the computed value
688 * in RCONDC.
689 *
690  result( 6 ) = dget06( rcond, rcondc )
691 *
692 * Print information about the tests that did not pass
693 * the threshold.
694 *
695  DO 85 k = k1, 6
696  IF( result( k ).GE.thresh ) THEN
697  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
698  $ CALL aladhd( nout, path )
699  IF( prefac ) THEN
700  WRITE( nout, fmt = 9997 )'DPOSVXX', fact,
701  $ uplo, n, equed, imat, k, result( k )
702  ELSE
703  WRITE( nout, fmt = 9998 )'DPOSVXX', fact,
704  $ uplo, n, imat, k, result( k )
705  END IF
706  nfail = nfail + 1
707  END IF
708  85 CONTINUE
709  nrun = nrun + 7 - k1
710  90 CONTINUE
711  100 CONTINUE
712  110 CONTINUE
713  120 CONTINUE
714  130 CONTINUE
715 *
716 * Print a summary of the results.
717 *
718  CALL alasvm( path, nout, nfail, nrun, nerrs )
719 *
720 
721 * Test Error Bounds from DPOSVXX
722 
723  CALL debchvxx( thresh, path )
724 
725  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
726  $ ', test(', i1, ')=', g12.5 )
727  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
728  $ ', type ', i1, ', test(', i1, ')=', g12.5 )
729  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
730  $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ') =',
731  $ g12.5 )
732  RETURN
733 *
734 * End of DDRVPO
735 *
736  END
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
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:112
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric matrix.
Definition: dlansy.f:124
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:109
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:166
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine debchvxx(THRESH, PATH)
DEBCHVXX
Definition: debchvxx.f:98
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:132
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
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:135
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine derrvx(PATH, NUNIT)
DERRVX
Definition: derrvx.f:57
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dpot01(UPLO, N, A, LDA, AFAC, LDAFAC, RWORK, RESID)
DPOT01
Definition: dpot01.f:106
subroutine dpot05(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
DPOT05
Definition: dpot05.f:166
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:309
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dpoequ(N, A, LDA, S, SCOND, AMAX, INFO)
DPOEQU
Definition: dpoequ.f:114
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine dpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DPOT02
Definition: dpot02.f:129
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:496
subroutine dpotri(UPLO, N, A, LDA, INFO)
DPOTRI
Definition: dpotri.f:97