LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
ddrvge.f
Go to the documentation of this file.
1 *> \brief \b DDRVGE
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 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 column 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 DOUBLE PRECISION
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 DOUBLE PRECISION array, dimension (NMAX*NMAX)
91 *> \endverbatim
92 *>
93 *> \param[out] AFAC
94 *> \verbatim
95 *> AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
96 *> \endverbatim
97 *>
98 *> \param[out] ASAV
99 *> \verbatim
100 *> ASAV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
101 *> \endverbatim
102 *>
103 *> \param[out] B
104 *> \verbatim
105 *> B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
106 *> \endverbatim
107 *>
108 *> \param[out] BSAV
109 *> \verbatim
110 *> BSAV is DOUBLE PRECISION array, dimension (NMAX*NRHS)
111 *> \endverbatim
112 *>
113 *> \param[out] X
114 *> \verbatim
115 *> X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
116 *> \endverbatim
117 *>
118 *> \param[out] XACT
119 *> \verbatim
120 *> XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
121 *> \endverbatim
122 *>
123 *> \param[out] S
124 *> \verbatim
125 *> S is DOUBLE PRECISION array, dimension (2*NMAX)
126 *> \endverbatim
127 *>
128 *> \param[out] WORK
129 *> \verbatim
130 *> WORK is DOUBLE PRECISION array, dimension
131 *> (NMAX*max(3,NRHS))
132 *> \endverbatim
133 *>
134 *> \param[out] RWORK
135 *> \verbatim
136 *> RWORK is DOUBLE PRECISION array, dimension (2*NRHS+NMAX)
137 *> \endverbatim
138 *>
139 *> \param[out] IWORK
140 *> \verbatim
141 *> IWORK is INTEGER array, dimension (2*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 *> \date November 2015
159 *
160 *> \ingroup double_lin
161 *
162 * =====================================================================
163  SUBROUTINE ddrvge( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
164  $ a, afac, asav, b, bsav, x, xact, s, work,
165  $ rwork, iwork, nout )
166 *
167 * -- LAPACK test routine (version 3.6.0) --
168 * -- LAPACK is a software package provided by Univ. of Tennessee, --
169 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
170 * November 2015
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 = 11 )
192  INTEGER NTESTS
193  parameter ( ntests = 7 )
194  INTEGER NTRAN
195  parameter ( ntran = 3 )
196 * ..
197 * .. Local Scalars ..
198  LOGICAL EQUIL, NOFACT, PREFAC, TRFCON, ZEROT
199  CHARACTER DIST, EQUED, FACT, TRANS, TYPE, XTYPE
200  CHARACTER*3 PATH
201  INTEGER I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, ITRAN,
202  $ izero, k, k1, kl, ku, lda, lwork, mode, n, nb,
203  $ nbmin, nerrs, nfact, nfail, nimat, nrun, nt
204  DOUBLE PRECISION AINVNM, AMAX, ANORM, ANORMI, ANORMO, CNDNUM,
205  $ colcnd, rcond, rcondc, rcondi, rcondo, roldc,
206  $ roldi, roldo, rowcnd, rpvgrw
207 * ..
208 * .. Local Arrays ..
209  CHARACTER EQUEDS( 4 ), FACTS( 3 ), TRANSS( ntran )
210  INTEGER ISEED( 4 ), ISEEDY( 4 )
211  DOUBLE PRECISION RESULT( ntests )
212 * ..
213 * .. External Functions ..
214  LOGICAL LSAME
215  DOUBLE PRECISION DGET06, DLAMCH, DLANGE, DLANTR
216  EXTERNAL lsame, dget06, dlamch, dlange, dlantr
217 * ..
218 * .. External Subroutines ..
219  EXTERNAL aladhd, alaerh, alasvm, derrvx, dgeequ, dgesv,
222  $ dlatms, xlaenv
223 * ..
224 * .. Intrinsic Functions ..
225  INTRINSIC abs, max
226 * ..
227 * .. Scalars in Common ..
228  LOGICAL LERR, OK
229  CHARACTER*32 SRNAMT
230  INTEGER INFOT, NUNIT
231 * ..
232 * .. Common blocks ..
233  COMMON / infoc / infot, nunit, ok, lerr
234  COMMON / srnamc / srnamt
235 * ..
236 * .. Data statements ..
237  DATA iseedy / 1988, 1989, 1990, 1991 /
238  DATA transs / 'N', 'T', 'C' /
239  DATA facts / 'F', 'N', 'E' /
240  DATA equeds / 'N', 'R', 'C', 'B' /
241 * ..
242 * .. Executable Statements ..
243 *
244 * Initialize constants and the random number seed.
245 *
246  path( 1: 1 ) = 'Double precision'
247  path( 2: 3 ) = 'GE'
248  nrun = 0
249  nfail = 0
250  nerrs = 0
251  DO 10 i = 1, 4
252  iseed( i ) = iseedy( i )
253  10 CONTINUE
254 *
255 * Test the error exits
256 *
257  IF( tsterr )
258  $ CALL derrvx( path, nout )
259  infot = 0
260 *
261 * Set the block size and minimum block size for testing.
262 *
263  nb = 1
264  nbmin = 2
265  CALL xlaenv( 1, nb )
266  CALL xlaenv( 2, nbmin )
267 *
268 * Do for each value of N in NVAL
269 *
270  DO 90 in = 1, nn
271  n = nval( in )
272  lda = max( n, 1 )
273  xtype = 'N'
274  nimat = ntypes
275  IF( n.LE.0 )
276  $ nimat = 1
277 *
278  DO 80 imat = 1, nimat
279 *
280 * Do the tests only if DOTYPE( IMAT ) is true.
281 *
282  IF( .NOT.dotype( imat ) )
283  $ GO TO 80
284 *
285 * Skip types 5, 6, or 7 if the matrix size is too small.
286 *
287  zerot = imat.GE.5 .AND. imat.LE.7
288  IF( zerot .AND. n.LT.imat-4 )
289  $ GO TO 80
290 *
291 * Set up parameters with DLATB4 and generate a test matrix
292 * with DLATMS.
293 *
294  CALL dlatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
295  $ cndnum, dist )
296  rcondc = one / cndnum
297 *
298  srnamt = 'DLATMS'
299  CALL dlatms( n, n, dist, iseed, TYPE, RWORK, MODE, CNDNUM,
300  $ anorm, kl, ku, 'No packing', a, lda, work,
301  $ info )
302 *
303 * Check error code from DLATMS.
304 *
305  IF( info.NE.0 ) THEN
306  CALL alaerh( path, 'DLATMS', info, 0, ' ', n, n, -1, -1,
307  $ -1, imat, nfail, nerrs, nout )
308  GO TO 80
309  END IF
310 *
311 * For types 5-7, zero one or more columns of the matrix to
312 * test that INFO is returned correctly.
313 *
314  IF( zerot ) THEN
315  IF( imat.EQ.5 ) THEN
316  izero = 1
317  ELSE IF( imat.EQ.6 ) THEN
318  izero = n
319  ELSE
320  izero = n / 2 + 1
321  END IF
322  ioff = ( izero-1 )*lda
323  IF( imat.LT.7 ) THEN
324  DO 20 i = 1, n
325  a( ioff+i ) = zero
326  20 CONTINUE
327  ELSE
328  CALL dlaset( 'Full', n, n-izero+1, zero, zero,
329  $ a( ioff+1 ), lda )
330  END IF
331  ELSE
332  izero = 0
333  END IF
334 *
335 * Save a copy of the matrix A in ASAV.
336 *
337  CALL dlacpy( 'Full', n, n, a, lda, asav, lda )
338 *
339  DO 70 iequed = 1, 4
340  equed = equeds( iequed )
341  IF( iequed.EQ.1 ) THEN
342  nfact = 3
343  ELSE
344  nfact = 1
345  END IF
346 *
347  DO 60 ifact = 1, nfact
348  fact = facts( ifact )
349  prefac = lsame( fact, 'F' )
350  nofact = lsame( fact, 'N' )
351  equil = lsame( fact, 'E' )
352 *
353  IF( zerot ) THEN
354  IF( prefac )
355  $ GO TO 60
356  rcondo = zero
357  rcondi = zero
358 *
359  ELSE IF( .NOT.nofact ) THEN
360 *
361 * Compute the condition number for comparison with
362 * the value returned by DGESVX (FACT = 'N' reuses
363 * the condition number from the previous iteration
364 * with FACT = 'F').
365 *
366  CALL dlacpy( 'Full', n, n, asav, lda, afac, lda )
367  IF( equil .OR. iequed.GT.1 ) THEN
368 *
369 * Compute row and column scale factors to
370 * equilibrate the matrix A.
371 *
372  CALL dgeequ( n, n, afac, lda, s, s( n+1 ),
373  $ rowcnd, colcnd, amax, info )
374  IF( info.EQ.0 .AND. n.GT.0 ) THEN
375  IF( lsame( equed, 'R' ) ) THEN
376  rowcnd = zero
377  colcnd = one
378  ELSE IF( lsame( equed, 'C' ) ) THEN
379  rowcnd = one
380  colcnd = zero
381  ELSE IF( lsame( equed, 'B' ) ) THEN
382  rowcnd = zero
383  colcnd = zero
384  END IF
385 *
386 * Equilibrate the matrix.
387 *
388  CALL dlaqge( n, n, afac, lda, s, s( n+1 ),
389  $ rowcnd, colcnd, amax, equed )
390  END IF
391  END IF
392 *
393 * Save the condition number of the non-equilibrated
394 * system for use in DGET04.
395 *
396  IF( equil ) THEN
397  roldo = rcondo
398  roldi = rcondi
399  END IF
400 *
401 * Compute the 1-norm and infinity-norm of A.
402 *
403  anormo = dlange( '1', n, n, afac, lda, rwork )
404  anormi = dlange( 'I', n, n, afac, lda, rwork )
405 *
406 * Factor the matrix A.
407 *
408  srnamt = 'DGETRF'
409  CALL dgetrf( n, n, afac, lda, iwork, info )
410 *
411 * Form the inverse of A.
412 *
413  CALL dlacpy( 'Full', n, n, afac, lda, a, lda )
414  lwork = nmax*max( 3, nrhs )
415  srnamt = 'DGETRI'
416  CALL dgetri( n, a, lda, iwork, work, lwork, info )
417 *
418 * Compute the 1-norm condition number of A.
419 *
420  ainvnm = dlange( '1', n, n, a, lda, rwork )
421  IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
422  rcondo = one
423  ELSE
424  rcondo = ( one / anormo ) / ainvnm
425  END IF
426 *
427 * Compute the infinity-norm condition number of A.
428 *
429  ainvnm = dlange( 'I', n, n, a, lda, rwork )
430  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
431  rcondi = one
432  ELSE
433  rcondi = ( one / anormi ) / ainvnm
434  END IF
435  END IF
436 *
437  DO 50 itran = 1, ntran
438 *
439 * Do for each value of TRANS.
440 *
441  trans = transs( itran )
442  IF( itran.EQ.1 ) THEN
443  rcondc = rcondo
444  ELSE
445  rcondc = rcondi
446  END IF
447 *
448 * Restore the matrix A.
449 *
450  CALL dlacpy( 'Full', n, n, asav, lda, a, lda )
451 *
452 * Form an exact solution and set the right hand side.
453 *
454  srnamt = 'DLARHS'
455  CALL dlarhs( path, xtype, 'Full', trans, n, n, kl,
456  $ ku, nrhs, a, lda, xact, lda, b, lda,
457  $ iseed, info )
458  xtype = 'C'
459  CALL dlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
460 *
461  IF( nofact .AND. itran.EQ.1 ) THEN
462 *
463 * --- Test DGESV ---
464 *
465 * Compute the LU factorization of the matrix and
466 * solve the system.
467 *
468  CALL dlacpy( 'Full', n, n, a, lda, afac, lda )
469  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
470 *
471  srnamt = 'DGESV '
472  CALL dgesv( n, nrhs, afac, lda, iwork, x, lda,
473  $ info )
474 *
475 * Check error code from DGESV .
476 *
477  IF( info.NE.izero )
478  $ CALL alaerh( path, 'DGESV ', info, izero,
479  $ ' ', n, n, -1, -1, nrhs, imat,
480  $ nfail, nerrs, nout )
481 *
482 * Reconstruct matrix from factors and compute
483 * residual.
484 *
485  CALL dget01( n, n, a, lda, afac, lda, iwork,
486  $ rwork, result( 1 ) )
487  nt = 1
488  IF( izero.EQ.0 ) THEN
489 *
490 * Compute residual of the computed solution.
491 *
492  CALL dlacpy( 'Full', n, nrhs, b, lda, work,
493  $ lda )
494  CALL dget02( 'No transpose', n, n, nrhs, a,
495  $ lda, x, lda, work, lda, rwork,
496  $ result( 2 ) )
497 *
498 * Check solution from generated exact solution.
499 *
500  CALL dget04( n, nrhs, x, lda, xact, lda,
501  $ rcondc, result( 3 ) )
502  nt = 3
503  END IF
504 *
505 * Print information about the tests that did not
506 * pass the threshold.
507 *
508  DO 30 k = 1, nt
509  IF( result( k ).GE.thresh ) THEN
510  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
511  $ CALL aladhd( nout, path )
512  WRITE( nout, fmt = 9999 )'DGESV ', n,
513  $ imat, k, result( k )
514  nfail = nfail + 1
515  END IF
516  30 CONTINUE
517  nrun = nrun + nt
518  END IF
519 *
520 * --- Test DGESVX ---
521 *
522  IF( .NOT.prefac )
523  $ CALL dlaset( 'Full', n, n, zero, zero, afac,
524  $ lda )
525  CALL dlaset( 'Full', n, nrhs, zero, zero, x, lda )
526  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
527 *
528 * Equilibrate the matrix if FACT = 'F' and
529 * EQUED = 'R', 'C', or 'B'.
530 *
531  CALL dlaqge( n, n, a, lda, s, s( n+1 ), rowcnd,
532  $ colcnd, amax, equed )
533  END IF
534 *
535 * Solve the system and compute the condition number
536 * and error bounds using DGESVX.
537 *
538  srnamt = 'DGESVX'
539  CALL dgesvx( fact, trans, n, nrhs, a, lda, afac,
540  $ lda, iwork, equed, s, s( n+1 ), b,
541  $ lda, x, lda, rcond, rwork,
542  $ rwork( nrhs+1 ), work, iwork( n+1 ),
543  $ info )
544 *
545 * Check the error code from DGESVX.
546 *
547  IF( info.NE.izero )
548  $ CALL alaerh( path, 'DGESVX', info, izero,
549  $ fact // trans, n, n, -1, -1, nrhs,
550  $ imat, nfail, nerrs, nout )
551 *
552 * Compare WORK(1) from DGESVX with the computed
553 * reciprocal pivot growth factor RPVGRW
554 *
555  IF( info.NE.0 .AND. info.LE.n) THEN
556  rpvgrw = dlantr( 'M', 'U', 'N', info, info,
557  $ afac, lda, work )
558  IF( rpvgrw.EQ.zero ) THEN
559  rpvgrw = one
560  ELSE
561  rpvgrw = dlange( 'M', n, info, a, lda,
562  $ work ) / rpvgrw
563  END IF
564  ELSE
565  rpvgrw = dlantr( 'M', 'U', 'N', n, n, afac, lda,
566  $ work )
567  IF( rpvgrw.EQ.zero ) THEN
568  rpvgrw = one
569  ELSE
570  rpvgrw = dlange( 'M', n, n, a, lda, work ) /
571  $ rpvgrw
572  END IF
573  END IF
574  result( 7 ) = abs( rpvgrw-work( 1 ) ) /
575  $ max( work( 1 ), rpvgrw ) /
576  $ dlamch( 'E' )
577 *
578  IF( .NOT.prefac ) THEN
579 *
580 * Reconstruct matrix from factors and compute
581 * residual.
582 *
583  CALL dget01( n, n, a, lda, afac, lda, iwork,
584  $ rwork( 2*nrhs+1 ), result( 1 ) )
585  k1 = 1
586  ELSE
587  k1 = 2
588  END IF
589 *
590  IF( info.EQ.0 ) THEN
591  trfcon = .false.
592 *
593 * Compute residual of the computed solution.
594 *
595  CALL dlacpy( 'Full', n, nrhs, bsav, lda, work,
596  $ lda )
597  CALL dget02( trans, n, n, nrhs, asav, lda, x,
598  $ lda, work, lda, rwork( 2*nrhs+1 ),
599  $ result( 2 ) )
600 *
601 * Check solution from generated exact solution.
602 *
603  IF( nofact .OR. ( prefac .AND. lsame( equed,
604  $ 'N' ) ) ) THEN
605  CALL dget04( n, nrhs, x, lda, xact, lda,
606  $ rcondc, result( 3 ) )
607  ELSE
608  IF( itran.EQ.1 ) THEN
609  roldc = roldo
610  ELSE
611  roldc = roldi
612  END IF
613  CALL dget04( n, nrhs, x, lda, xact, lda,
614  $ roldc, result( 3 ) )
615  END IF
616 *
617 * Check the error bounds from iterative
618 * refinement.
619 *
620  CALL dget07( trans, n, nrhs, asav, lda, b, lda,
621  $ x, lda, xact, lda, rwork, .true.,
622  $ rwork( nrhs+1 ), result( 4 ) )
623  ELSE
624  trfcon = .true.
625  END IF
626 *
627 * Compare RCOND from DGESVX with the computed value
628 * in RCONDC.
629 *
630  result( 6 ) = dget06( rcond, rcondc )
631 *
632 * Print information about the tests that did not pass
633 * the threshold.
634 *
635  IF( .NOT.trfcon ) THEN
636  DO 40 k = k1, ntests
637  IF( result( k ).GE.thresh ) THEN
638  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
639  $ CALL aladhd( nout, path )
640  IF( prefac ) THEN
641  WRITE( nout, fmt = 9997 )'DGESVX',
642  $ fact, trans, n, equed, imat, k,
643  $ result( k )
644  ELSE
645  WRITE( nout, fmt = 9998 )'DGESVX',
646  $ fact, trans, n, imat, k, result( k )
647  END IF
648  nfail = nfail + 1
649  END IF
650  40 CONTINUE
651  nrun = nrun + ntests - k1 + 1
652  ELSE
653  IF( result( 1 ).GE.thresh .AND. .NOT.prefac )
654  $ THEN
655  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
656  $ CALL aladhd( nout, path )
657  IF( prefac ) THEN
658  WRITE( nout, fmt = 9997 )'DGESVX', fact,
659  $ trans, n, equed, imat, 1, result( 1 )
660  ELSE
661  WRITE( nout, fmt = 9998 )'DGESVX', fact,
662  $ trans, n, imat, 1, result( 1 )
663  END IF
664  nfail = nfail + 1
665  nrun = nrun + 1
666  END IF
667  IF( result( 6 ).GE.thresh ) THEN
668  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
669  $ CALL aladhd( nout, path )
670  IF( prefac ) THEN
671  WRITE( nout, fmt = 9997 )'DGESVX', fact,
672  $ trans, n, equed, imat, 6, result( 6 )
673  ELSE
674  WRITE( nout, fmt = 9998 )'DGESVX', fact,
675  $ trans, n, imat, 6, result( 6 )
676  END IF
677  nfail = nfail + 1
678  nrun = nrun + 1
679  END IF
680  IF( result( 7 ).GE.thresh ) THEN
681  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
682  $ CALL aladhd( nout, path )
683  IF( prefac ) THEN
684  WRITE( nout, fmt = 9997 )'DGESVX', fact,
685  $ trans, n, equed, imat, 7, result( 7 )
686  ELSE
687  WRITE( nout, fmt = 9998 )'DGESVX', fact,
688  $ trans, n, imat, 7, result( 7 )
689  END IF
690  nfail = nfail + 1
691  nrun = nrun + 1
692  END IF
693 *
694  END IF
695 *
696  50 CONTINUE
697  60 CONTINUE
698  70 CONTINUE
699  80 CONTINUE
700  90 CONTINUE
701 *
702 * Print a summary of the results.
703 *
704  CALL alasvm( path, nout, nfail, nrun, nerrs )
705 *
706  9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test(', i2, ') =',
707  $ g12.5 )
708  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
709  $ ', type ', i2, ', test(', i1, ')=', g12.5 )
710  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N=', i5,
711  $ ', EQUED=''', a1, ''', type ', i2, ', test(', i1, ')=',
712  $ g12.5 )
713  RETURN
714 *
715 * End of DDRVGE
716 *
717  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
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 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
subroutine dget07(TRANS, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, CHKFERR, BERR, RESLTS)
DGET07
Definition: dget07.f:167
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
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323