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,
226     $                   DGESVX, DGET01, DGET02, DGET04, DGET07, DGETRF,
227     $                   DGETRI, DLACPY, DLAQGE, DLARHS, DLASET, DLATB4,
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
895