LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
ddrvgex.f
Go to the documentation of this file.
1 *> \brief \b DDRVGEX
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 DDRVGE( 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 *> DDRVGE tests the driver routines DGESV, -SVX, and -SVXX.
35 *>
36 *> Note that this file is used only when the XBLAS are available,
37 *> otherwise ddrvge.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 column 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 (2*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 (2*NRHS+NMAX)
140 *> \endverbatim
141 *>
142 *> \param[out] IWORK
143 *> \verbatim
144 *> IWORK is INTEGER array, dimension (2*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 April 2012
162 *
163 *> \ingroup double_lin
164 *
165 * =====================================================================
166  SUBROUTINE ddrvge( 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.4.1) --
171 * -- LAPACK is a software package provided by Univ. of Tennessee, --
172 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
173 * April 2012
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 = 11 )
195  INTEGER ntests
196  parameter ( ntests = 7 )
197  INTEGER ntran
198  parameter ( ntran = 3 )
199 * ..
200 * .. Local Scalars ..
201  LOGICAL equil, nofact, prefac, trfcon, zerot
202  CHARACTER dist, equed, fact, trans, TYPE, xtype
203  CHARACTER*3 path
204  INTEGER i, iequed, ifact, imat, in, info, ioff, itran,
205  $ izero, k, k1, kl, ku, lda, lwork, mode, n, nb,
206  $ nbmin, nerrs, nfact, nfail, nimat, nrun, nt,
207  $ n_err_bnds
208  DOUBLE PRECISION ainvnm, amax, anorm, anormi, anormo, cndnum,
209  $ colcnd, rcond, rcondc, rcondi, rcondo, roldc,
210  $ roldi, roldo, rowcnd, rpvgrw, rpvgrw_svxx
211 * ..
212 * .. Local Arrays ..
213  CHARACTER equeds( 4 ), facts( 3 ), transs( ntran )
214  INTEGER iseed( 4 ), iseedy( 4 )
215  DOUBLE PRECISION result( ntests ), berr( nrhs ),
216  $ errbnds_n( nrhs, 3 ), errbnds_c( nrhs, 3 )
217 * ..
218 * .. External Functions ..
219  LOGICAL lsame
220  DOUBLE PRECISION dget06, dlamch, dlange, dlantr, dla_gerpvgrw
221  EXTERNAL lsame, dget06, dlamch, dlange, dlantr,
222  $ dla_gerpvgrw
223 * ..
224 * .. External Subroutines ..
225  EXTERNAL aladhd, alaerh, alasvm, derrvx, dgeequ, dgesv,
228  $ dlatms, xlaenv, dgesvxx
229 * ..
230 * .. Intrinsic Functions ..
231  INTRINSIC abs, max
232 * ..
233 * .. Scalars in Common ..
234  LOGICAL lerr, ok
235  CHARACTER*32 srnamt
236  INTEGER infot, nunit
237 * ..
238 * .. Common blocks ..
239  COMMON / infoc / infot, nunit, ok, lerr
240  COMMON / srnamc / srnamt
241 * ..
242 * .. Data statements ..
243  DATA iseedy / 1988, 1989, 1990, 1991 /
244  DATA transs / 'N', 'T', 'C' /
245  DATA facts / 'F', 'N', 'E' /
246  DATA equeds / 'N', 'R', 'C', 'B' /
247 * ..
248 * .. Executable Statements ..
249 *
250 * Initialize constants and the random number seed.
251 *
252  path( 1: 1 ) = 'Double precision'
253  path( 2: 3 ) = 'GE'
254  nrun = 0
255  nfail = 0
256  nerrs = 0
257  DO 10 i = 1, 4
258  iseed( i ) = iseedy( i )
259  10 CONTINUE
260 *
261 * Test the error exits
262 *
263  IF( tsterr )
264  $ CALL derrvx( path, nout )
265  infot = 0
266 *
267 * Set the block size and minimum block size for testing.
268 *
269  nb = 1
270  nbmin = 2
271  CALL xlaenv( 1, nb )
272  CALL xlaenv( 2, nbmin )
273 *
274 * Do for each value of N in NVAL
275 *
276  DO 90 in = 1, nn
277  n = nval( in )
278  lda = max( n, 1 )
279  xtype = 'N'
280  nimat = ntypes
281  IF( n.LE.0 )
282  $ nimat = 1
283 *
284  DO 80 imat = 1, nimat
285 *
286 * Do the tests only if DOTYPE( IMAT ) is true.
287 *
288  IF( .NOT.dotype( imat ) )
289  $ GO TO 80
290 *
291 * Skip types 5, 6, or 7 if the matrix size is too small.
292 *
293  zerot = imat.GE.5 .AND. imat.LE.7
294  IF( zerot .AND. n.LT.imat-4 )
295  $ GO TO 80
296 *
297 * Set up parameters with DLATB4 and generate a test matrix
298 * with DLATMS.
299 *
300  CALL dlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
301  $ cndnum, dist )
302  rcondc = one / cndnum
303 *
304  srnamt = 'DLATMS'
305  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode, cndnum,
306  $ anorm, kl, ku, 'No packing', 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, ' ', n, n, -1, -1,
313  $ -1, imat, nfail, nerrs, nout )
314  GO TO 80
315  END IF
316 *
317 * For types 5-7, zero one or more columns of the matrix to
318 * test that INFO is returned correctly.
319 *
320  IF( zerot ) THEN
321  IF( imat.EQ.5 ) THEN
322  izero = 1
323  ELSE IF( imat.EQ.6 ) THEN
324  izero = n
325  ELSE
326  izero = n / 2 + 1
327  END IF
328  ioff = ( izero-1 )*lda
329  IF( imat.LT.7 ) THEN
330  DO 20 i = 1, n
331  a( ioff+i ) = zero
332  20 CONTINUE
333  ELSE
334  CALL dlaset( 'Full', n, n-izero+1, zero, zero,
335  $ a( ioff+1 ), lda )
336  END IF
337  ELSE
338  izero = 0
339  END IF
340 *
341 * Save a copy of the matrix A in ASAV.
342 *
343  CALL dlacpy( 'Full', n, n, a, lda, asav, lda )
344 *
345  DO 70 iequed = 1, 4
346  equed = equeds( iequed )
347  IF( iequed.EQ.1 ) THEN
348  nfact = 3
349  ELSE
350  nfact = 1
351  END IF
352 *
353  DO 60 ifact = 1, nfact
354  fact = facts( ifact )
355  prefac = lsame( fact, 'F' )
356  nofact = lsame( fact, 'N' )
357  equil = lsame( fact, 'E' )
358 *
359  IF( zerot ) THEN
360  IF( prefac )
361  $ GO TO 60
362  rcondo = zero
363  rcondi = zero
364 *
365  ELSE IF( .NOT.nofact ) THEN
366 *
367 * Compute the condition number for comparison with
368 * the value returned by DGESVX (FACT = 'N' reuses
369 * the condition number from the previous iteration
370 * with FACT = 'F').
371 *
372  CALL dlacpy( 'Full', n, n, asav, lda, afac, lda )
373  IF( equil .OR. iequed.GT.1 ) THEN
374 *
375 * Compute row and column scale factors to
376 * equilibrate the matrix A.
377 *
378  CALL dgeequ( n, n, afac, lda, s, s( n+1 ),
379  $ rowcnd, colcnd, amax, info )
380  IF( info.EQ.0 .AND. n.GT.0 ) THEN
381  IF( lsame( equed, 'R' ) ) THEN
382  rowcnd = zero
383  colcnd = one
384  ELSE IF( lsame( equed, 'C' ) ) THEN
385  rowcnd = one
386  colcnd = zero
387  ELSE IF( lsame( equed, 'B' ) ) THEN
388  rowcnd = zero
389  colcnd = zero
390  END IF
391 *
392 * Equilibrate the matrix.
393 *
394  CALL dlaqge( n, n, afac, lda, s, s( n+1 ),
395  $ rowcnd, colcnd, amax, equed )
396  END IF
397  END IF
398 *
399 * Save the condition number of the non-equilibrated
400 * system for use in DGET04.
401 *
402  IF( equil ) THEN
403  roldo = rcondo
404  roldi = rcondi
405  END IF
406 *
407 * Compute the 1-norm and infinity-norm of A.
408 *
409  anormo = dlange( '1', n, n, afac, lda, rwork )
410  anormi = dlange( 'I', n, n, afac, lda, rwork )
411 *
412 * Factor the matrix A.
413 *
414  CALL dgetrf( n, n, afac, lda, iwork, info )
415 *
416 * Form the inverse of A.
417 *
418  CALL dlacpy( 'Full', n, n, afac, lda, a, lda )
419  lwork = nmax*max( 3, nrhs )
420  CALL dgetri( n, a, lda, iwork, work, lwork, info )
421 *
422 * Compute the 1-norm condition number of A.
423 *
424  ainvnm = dlange( '1', n, n, a, lda, rwork )
425  IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
426  rcondo = one
427  ELSE
428  rcondo = ( one / anormo ) / ainvnm
429  END IF
430 *
431 * Compute the infinity-norm condition number of A.
432 *
433  ainvnm = dlange( 'I', n, n, a, lda, rwork )
434  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
435  rcondi = one
436  ELSE
437  rcondi = ( one / anormi ) / ainvnm
438  END IF
439  END IF
440 *
441  DO 50 itran = 1, ntran
442 *
443 * Do for each value of TRANS.
444 *
445  trans = transs( itran )
446  IF( itran.EQ.1 ) THEN
447  rcondc = rcondo
448  ELSE
449  rcondc = rcondi
450  END IF
451 *
452 * Restore the matrix A.
453 *
454  CALL dlacpy( 'Full', n, n, asav, lda, a, lda )
455 *
456 * Form an exact solution and set the right hand side.
457 *
458  srnamt = 'DLARHS'
459  CALL dlarhs( path, xtype, 'Full', trans, n, n, kl,
460  $ ku, nrhs, a, lda, xact, lda, b, lda,
461  $ iseed, info )
462  xtype = 'C'
463  CALL dlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
464 *
465  IF( nofact .AND. itran.EQ.1 ) THEN
466 *
467 * --- Test DGESV ---
468 *
469 * Compute the LU factorization of the matrix and
470 * solve the system.
471 *
472  CALL dlacpy( 'Full', n, n, a, lda, afac, lda )
473  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
474 *
475  srnamt = 'DGESV '
476  CALL dgesv( n, nrhs, afac, lda, iwork, x, lda,
477  $ info )
478 *
479 * Check error code from DGESV .
480 *
481  IF( info.NE.izero )
482  $ CALL alaerh( path, 'DGESV ', info, izero,
483  $ ' ', n, n, -1, -1, nrhs, imat,
484  $ nfail, nerrs, nout )
485 *
486 * Reconstruct matrix from factors and compute
487 * residual.
488 *
489  CALL dget01( n, n, a, lda, afac, lda, iwork,
490  $ rwork, result( 1 ) )
491  nt = 1
492  IF( izero.EQ.0 ) THEN
493 *
494 * Compute residual of the computed solution.
495 *
496  CALL dlacpy( 'Full', n, nrhs, b, lda, work,
497  $ lda )
498  CALL dget02( 'No transpose', n, n, nrhs, a,
499  $ lda, x, lda, work, lda, rwork,
500  $ result( 2 ) )
501 *
502 * Check solution from generated exact solution.
503 *
504  CALL dget04( n, nrhs, x, lda, xact, lda,
505  $ rcondc, result( 3 ) )
506  nt = 3
507  END IF
508 *
509 * Print information about the tests that did not
510 * pass the threshold.
511 *
512  DO 30 k = 1, nt
513  IF( result( k ).GE.thresh ) THEN
514  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
515  $ CALL aladhd( nout, path )
516  WRITE( nout, fmt = 9999 )'DGESV ', n,
517  $ imat, k, result( k )
518  nfail = nfail + 1
519  END IF
520  30 CONTINUE
521  nrun = nrun + nt
522  END IF
523 *
524 * --- Test DGESVX ---
525 *
526  IF( .NOT.prefac )
527  $ CALL dlaset( 'Full', n, n, zero, zero, afac,
528  $ lda )
529  CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
530  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
531 *
532 * Equilibrate the matrix if FACT = 'F' and
533 * EQUED = 'R', 'C', or 'B'.
534 *
535  CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
536  $ colcnd, amax, equed )
537  END IF
538 *
539 * Solve the system and compute the condition number
540 * and error bounds using DGESVX.
541 *
542  srnamt = 'DGESVX'
543  CALL dgesvx( fact, trans, n, nrhs, a, lda, afac,
544  $ lda, iwork, equed, s, s( n+1 ), b,
545  $ lda, x, lda, rcond, rwork,
546  $ rwork( nrhs+1 ), work, iwork( n+1 ),
547  $ info )
548 *
549 * Check the error code from DGESVX.
550 *
551  IF( info.NE.izero )
552  $ CALL alaerh( path, 'DGESVX', info, izero,
553  $ fact // trans, n, n, -1, -1, nrhs,
554  $ imat, nfail, nerrs, nout )
555 *
556 * Compare WORK(1) from DGESVX with the computed
557 * reciprocal pivot growth factor RPVGRW
558 *
559  IF( info.NE.0 ) THEN
560  rpvgrw = dlantr( 'M', 'U', 'N', info, info,
561  $ afac, lda, work )
562  IF( rpvgrw.EQ.zero ) THEN
563  rpvgrw = one
564  ELSE
565  rpvgrw = dlange( 'M', n, info, a, lda,
566  $ work ) / rpvgrw
567  END IF
568  ELSE
569  rpvgrw = dlantr( 'M', 'U', 'N', n, n, afac, lda,
570  $ work )
571  IF( rpvgrw.EQ.zero ) THEN
572  rpvgrw = one
573  ELSE
574  rpvgrw = dlange( 'M', n, n, a, lda, work ) /
575  $ rpvgrw
576  END IF
577  END IF
578  result( 7 ) = abs( rpvgrw-work( 1 ) ) /
579  $ max( work( 1 ), rpvgrw ) /
580  $ dlamch( 'E' )
581 *
582  IF( .NOT.prefac ) THEN
583 *
584 * Reconstruct matrix from factors and compute
585 * residual.
586 *
587  CALL dget01( n, n, a, lda, afac, lda, iwork,
588  $ rwork( 2*nrhs+1 ), result( 1 ) )
589  k1 = 1
590  ELSE
591  k1 = 2
592  END IF
593 *
594  IF( info.EQ.0 ) THEN
595  trfcon = .false.
596 *
597 * Compute residual of the computed solution.
598 *
599  CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
600  $ lda )
601  CALL dget02( trans, n, n, nrhs, asav, lda, x,
602  $ lda, work, lda, rwork( 2*nrhs+1 ),
603  $ result( 2 ) )
604 *
605 * Check solution from generated exact solution.
606 *
607  IF( nofact .OR. ( prefac .AND. lsame( equed,
608  $ 'N' ) ) ) THEN
609  CALL dget04( n, nrhs, x, lda, xact, lda,
610  $ rcondc, result( 3 ) )
611  ELSE
612  IF( itran.EQ.1 ) THEN
613  roldc = roldo
614  ELSE
615  roldc = roldi
616  END IF
617  CALL dget04( n, nrhs, x, lda, xact, lda,
618  $ roldc, result( 3 ) )
619  END IF
620 *
621 * Check the error bounds from iterative
622 * refinement.
623 *
624  CALL dget07( trans, n, nrhs, asav, lda, b, lda,
625  $ x, lda, xact, lda, rwork, .true.,
626  $ rwork( nrhs+1 ), result( 4 ) )
627  ELSE
628  trfcon = .true.
629  END IF
630 *
631 * Compare RCOND from DGESVX with the computed value
632 * in RCONDC.
633 *
634  result( 6 ) = dget06( rcond, rcondc )
635 *
636 * Print information about the tests that did not pass
637 * the threshold.
638 *
639  IF( .NOT.trfcon ) THEN
640  DO 40 k = k1, ntests
641  IF( result( k ).GE.thresh ) THEN
642  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
643  $ CALL aladhd( nout, path )
644  IF( prefac ) THEN
645  WRITE( nout, fmt = 9997 )'DGESVX',
646  $ fact, trans, n, equed, imat, k,
647  $ result( k )
648  ELSE
649  WRITE( nout, fmt = 9998 )'DGESVX',
650  $ fact, trans, n, imat, k, result( k )
651  END IF
652  nfail = nfail + 1
653  END IF
654  40 CONTINUE
655  nrun = nrun + 7 - k1
656  ELSE
657  IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
658  $ THEN
659  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
660  $ CALL aladhd( nout, path )
661  IF( prefac ) THEN
662  WRITE( nout, fmt = 9997 )'DGESVX', fact,
663  $ trans, n, equed, imat, 1, result( 1 )
664  ELSE
665  WRITE( nout, fmt = 9998 )'DGESVX', fact,
666  $ trans, n, imat, 1, result( 1 )
667  END IF
668  nfail = nfail + 1
669  nrun = nrun + 1
670  END IF
671  IF( result( 6 ).GE.thresh ) THEN
672  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
673  $ CALL aladhd( nout, path )
674  IF( prefac ) THEN
675  WRITE( nout, fmt = 9997 )'DGESVX', fact,
676  $ trans, n, equed, imat, 6, result( 6 )
677  ELSE
678  WRITE( nout, fmt = 9998 )'DGESVX', fact,
679  $ trans, n, imat, 6, result( 6 )
680  END IF
681  nfail = nfail + 1
682  nrun = nrun + 1
683  END IF
684  IF( result( 7 ).GE.thresh ) THEN
685  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
686  $ CALL aladhd( nout, path )
687  IF( prefac ) THEN
688  WRITE( nout, fmt = 9997 )'DGESVX', fact,
689  $ trans, n, equed, imat, 7, result( 7 )
690  ELSE
691  WRITE( nout, fmt = 9998 )'DGESVX', fact,
692  $ trans, n, imat, 7, result( 7 )
693  END IF
694  nfail = nfail + 1
695  nrun = nrun + 1
696  END IF
697 *
698  END IF
699 *
700 * --- Test DGESVXX ---
701 *
702 * Restore the matrices A and B.
703 *
704  CALL dlacpy( 'Full', n, n, asav, lda, a, lda )
705  CALL dlacpy( 'Full', n, nrhs, bsav, lda, b, lda )
706 
707  IF( .NOT.prefac )
708  $ CALL dlaset( 'Full', n, n, zero, zero, afac,
709  $ lda )
710  CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
711  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
712 *
713 * Equilibrate the matrix if FACT = 'F' and
714 * EQUED = 'R', 'C', or 'B'.
715 *
716  CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
717  $ colcnd, amax, equed )
718  END IF
719 *
720 * Solve the system and compute the condition number
721 * and error bounds using DGESVXX.
722 *
723  srnamt = 'DGESVXX'
724  n_err_bnds = 3
725  CALL dgesvxx( fact, trans, n, nrhs, a, lda, afac,
726  $ lda, iwork, equed, s, s( n+1 ), b, lda, x,
727  $ lda, rcond, rpvgrw_svxx, berr, n_err_bnds,
728  $ errbnds_n, errbnds_c, 0, zero, work,
729  $ iwork( n+1 ), info )
730 *
731 * Check the error code from DGESVXX.
732 *
733  IF( info.EQ.n+1 ) GOTO 50
734  IF( info.NE.izero ) THEN
735  CALL alaerh( path, 'DGESVXX', info, izero,
736  $ fact // trans, n, n, -1, -1, nrhs,
737  $ imat, nfail, nerrs, nout )
738  GOTO 50
739  END IF
740 *
741 * Compare rpvgrw_svxx from DGESVXX with the computed
742 * reciprocal pivot growth factor RPVGRW
743 *
744 
745  IF ( info .GT. 0 .AND. info .LT. n+1 ) THEN
746  rpvgrw = dla_gerpvgrw
747  $ (n, info, a, lda, afac, lda)
748  ELSE
749  rpvgrw = dla_gerpvgrw
750  $ (n, n, a, lda, afac, lda)
751  ENDIF
752 
753  result( 7 ) = abs( rpvgrw-rpvgrw_svxx ) /
754  $ max( rpvgrw_svxx, rpvgrw ) /
755  $ dlamch( 'E' )
756 *
757  IF( .NOT.prefac ) THEN
758 *
759 * Reconstruct matrix from factors and compute
760 * residual.
761 *
762  CALL dget01( n, n, a, lda, afac, lda, iwork,
763  $ rwork( 2*nrhs+1 ), result( 1 ) )
764  k1 = 1
765  ELSE
766  k1 = 2
767  END IF
768 *
769  IF( info.EQ.0 ) THEN
770  trfcon = .false.
771 *
772 * Compute residual of the computed solution.
773 *
774  CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
775  $ lda )
776  CALL dget02( trans, n, n, nrhs, asav, lda, x,
777  $ lda, work, lda, rwork( 2*nrhs+1 ),
778  $ result( 2 ) )
779 *
780 * Check solution from generated exact solution.
781 *
782  IF( nofact .OR. ( prefac .AND. lsame( equed,
783  $ 'N' ) ) ) THEN
784  CALL dget04( n, nrhs, x, lda, xact, lda,
785  $ rcondc, result( 3 ) )
786  ELSE
787  IF( itran.EQ.1 ) THEN
788  roldc = roldo
789  ELSE
790  roldc = roldi
791  END IF
792  CALL dget04( n, nrhs, x, lda, xact, lda,
793  $ roldc, result( 3 ) )
794  END IF
795  ELSE
796  trfcon = .true.
797  END IF
798 *
799 * Compare RCOND from DGESVXX with the computed value
800 * in RCONDC.
801 *
802  result( 6 ) = dget06( rcond, rcondc )
803 *
804 * Print information about the tests that did not pass
805 * the threshold.
806 *
807  IF( .NOT.trfcon ) THEN
808  DO 45 k = k1, ntests
809  IF( result( k ).GE.thresh ) THEN
810  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
811  $ CALL aladhd( nout, path )
812  IF( prefac ) THEN
813  WRITE( nout, fmt = 9997 )'DGESVXX',
814  $ fact, trans, n, equed, imat, k,
815  $ result( k )
816  ELSE
817  WRITE( nout, fmt = 9998 )'DGESVXX',
818  $ fact, trans, n, imat, k, result( k )
819  END IF
820  nfail = nfail + 1
821  END IF
822  45 CONTINUE
823  nrun = nrun + 7 - k1
824  ELSE
825  IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
826  $ THEN
827  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
828  $ CALL aladhd( nout, path )
829  IF( prefac ) THEN
830  WRITE( nout, fmt = 9997 )'DGESVXX', fact,
831  $ trans, n, equed, imat, 1, result( 1 )
832  ELSE
833  WRITE( nout, fmt = 9998 )'DGESVXX', fact,
834  $ trans, n, imat, 1, result( 1 )
835  END IF
836  nfail = nfail + 1
837  nrun = nrun + 1
838  END IF
839  IF( result( 6 ).GE.thresh ) THEN
840  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
841  $ CALL aladhd( nout, path )
842  IF( prefac ) THEN
843  WRITE( nout, fmt = 9997 )'DGESVXX', fact,
844  $ trans, n, equed, imat, 6, result( 6 )
845  ELSE
846  WRITE( nout, fmt = 9998 )'DGESVXX', fact,
847  $ trans, n, imat, 6, result( 6 )
848  END IF
849  nfail = nfail + 1
850  nrun = nrun + 1
851  END IF
852  IF( result( 7 ).GE.thresh ) THEN
853  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
854  $ CALL aladhd( nout, path )
855  IF( prefac ) THEN
856  WRITE( nout, fmt = 9997 )'DGESVXX', fact,
857  $ trans, n, equed, imat, 7, result( 7 )
858  ELSE
859  WRITE( nout, fmt = 9998 )'DGESVXX', fact,
860  $ trans, n, imat, 7, result( 7 )
861  END IF
862  nfail = nfail + 1
863  nrun = nrun + 1
864  END IF
865 *
866  END IF
867 *
868  50 CONTINUE
869  60 CONTINUE
870  70 CONTINUE
871  80 CONTINUE
872  90 CONTINUE
873 *
874 * Print a summary of the results.
875 *
876  CALL alasvm( path, nout, nfail, nrun, nerrs )
877 *
878 
879 * Test Error Bounds from DGESVXX
880 
881  CALL debchvxx( thresh, path )
882 
883  9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test(', i2, ') =',
884  $ g12.5 )
885  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
886  $ ', type ', i2, ', test(', i1, ')=', g12.5 )
887  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
888  $ ', EQUED=''', a1, ''', type ', i2, ', test(', i1, ')=',
889  $ g12.5 )
890  RETURN
891 *
892 * End of DDRVGE
893 *
894  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
subroutine dgetrf(M, N, A, LDA, IPIV, INFO)
DGETRF
Definition: dgetrf.f:110
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 ddrvge(DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX, A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK, RWORK, IWORK, NOUT)
DDRVGE
Definition: ddrvge.f:166
subroutine dlaqge(M, N, A, LDA, R, C, ROWCND, COLCND, AMAX, EQUED)
DLAQGE scales a general rectangular matrix, using row and column scaling factors computed by sgeequ...
Definition: dlaqge.f:144
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
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 dgetri(N, A, LDA, IPIV, WORK, LWORK, INFO)
DGETRI
Definition: dgetri.f:116
subroutine dgeequ(M, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFO)
DGEEQU
Definition: dgeequ.f:141
subroutine dget01(M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK, RESID)
DGET01
Definition: dget01.f:109
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
double precision function dlange(NORM, M, N, A, LDA, WORK)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlange.f:116
subroutine dget07(TRANS, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, CHKFERR, BERR, RESLTS)
DGET07
Definition: dget07.f:167
subroutine dgesvxx(FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, EQUED, R, C, B, LDB, X, LDX, RCOND, RPVGRW, BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO)
DGESVXX computes the solution to system of linear equations A * X = B for GE matrices ...
Definition: dgesvxx.f:542
subroutine dgesvx(FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
DGESVX computes the solution to system of linear equations A * X = B for GE matrices ...
Definition: dgesvx.f:351
subroutine dgesv(N, NRHS, A, LDA, IPIV, B, LDB, INFO)
DGESV computes the solution to system of linear equations A * X = B for GE matrices ...
Definition: dgesv.f:124
subroutine dget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DGET02
Definition: dget02.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 dlantr(NORM, UPLO, DIAG, M, N, A, LDA, WORK)
DLANTR returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a trapezoidal or triangular matrix.
Definition: dlantr.f:143
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
double precision function dla_gerpvgrw(N, NCOLS, A, LDA, AF, LDAF)
DLA_GERPVGRW
Definition: dla_gerpvgrw.f:102
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55