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*> \ingroup double_lin
159*
160*  =====================================================================
161      SUBROUTINE DDRVGE( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
162     $                   A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
163     $                   RWORK, IWORK, NOUT )
164*
165*  -- LAPACK test routine --
166*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
167*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
168*
169*     .. Scalar Arguments ..
170      LOGICAL            TSTERR
171      INTEGER            NMAX, NN, NOUT, NRHS
172      DOUBLE PRECISION   THRESH
173*     ..
174*     .. Array Arguments ..
175      LOGICAL            DOTYPE( * )
176      INTEGER            IWORK( * ), NVAL( * )
177      DOUBLE PRECISION   A( * ), AFAC( * ), ASAV( * ), B( * ),
178     $                   BSAV( * ), RWORK( * ), S( * ), WORK( * ),
179     $                   X( * ), XACT( * )
180*     ..
181*
182*  =====================================================================
183*
184*     .. Parameters ..
185      DOUBLE PRECISION   ONE, ZERO
186      PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
187      INTEGER            NTYPES
188      PARAMETER          ( NTYPES = 11 )
189      INTEGER            NTESTS
190      PARAMETER          ( NTESTS = 7 )
191      INTEGER            NTRAN
192      PARAMETER          ( NTRAN = 3 )
193*     ..
194*     .. Local Scalars ..
195      LOGICAL            EQUIL, NOFACT, PREFAC, TRFCON, ZEROT
196      CHARACTER          DIST, EQUED, FACT, TRANS, TYPE, XTYPE
197      CHARACTER*3        PATH
198      INTEGER            I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, ITRAN,
199     $                   IZERO, K, K1, KL, KU, LDA, LWORK, MODE, N, NB,
200     $                   NBMIN, NERRS, NFACT, NFAIL, NIMAT, NRUN, NT
201      DOUBLE PRECISION   AINVNM, AMAX, ANORM, ANORMI, ANORMO, CNDNUM,
202     $                   COLCND, RCOND, RCONDC, RCONDI, RCONDO, ROLDC,
203     $                   ROLDI, ROLDO, ROWCND, RPVGRW
204*     ..
205*     .. Local Arrays ..
206      CHARACTER          EQUEDS( 4 ), FACTS( 3 ), TRANSS( NTRAN )
207      INTEGER            ISEED( 4 ), ISEEDY( 4 )
208      DOUBLE PRECISION   RESULT( NTESTS )
209*     ..
210*     .. External Functions ..
211      LOGICAL            LSAME
212      DOUBLE PRECISION   DGET06, DLAMCH, DLANGE, DLANTR
213      EXTERNAL           LSAME, DGET06, DLAMCH, DLANGE, DLANTR
214*     ..
215*     .. External Subroutines ..
216      EXTERNAL           ALADHD, ALAERH, ALASVM, DERRVX, DGEEQU, DGESV,
217     $                   DGESVX, DGET01, DGET02, DGET04, DGET07, DGETRF,
218     $                   DGETRI, DLACPY, DLAQGE, DLARHS, DLASET, DLATB4,
219     $                   DLATMS, XLAENV
220*     ..
221*     .. Intrinsic Functions ..
222      INTRINSIC          ABS, MAX
223*     ..
224*     .. Scalars in Common ..
225      LOGICAL            LERR, OK
226      CHARACTER*32       SRNAMT
227      INTEGER            INFOT, NUNIT
228*     ..
229*     .. Common blocks ..
230      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
231      COMMON             / SRNAMC / SRNAMT
232*     ..
233*     .. Data statements ..
234      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
235      DATA               TRANSS / 'N', 'T', 'C' /
236      DATA               FACTS / 'F', 'N', 'E' /
237      DATA               EQUEDS / 'N', 'R', 'C', 'B' /
238*     ..
239*     .. Executable Statements ..
240*
241*     Initialize constants and the random number seed.
242*
243      PATH( 1: 1 ) = 'Double precision'
244      PATH( 2: 3 ) = 'GE'
245      NRUN = 0
246      NFAIL = 0
247      NERRS = 0
248      DO 10 I = 1, 4
249         ISEED( I ) = ISEEDY( I )
250   10 CONTINUE
251*
252*     Test the error exits
253*
254      IF( TSTERR )
255     $   CALL DERRVX( PATH, NOUT )
256      INFOT = 0
257*
258*     Set the block size and minimum block size for testing.
259*
260      NB = 1
261      NBMIN = 2
262      CALL XLAENV( 1, NB )
263      CALL XLAENV( 2, NBMIN )
264*
265*     Do for each value of N in NVAL
266*
267      DO 90 IN = 1, NN
268         N = NVAL( IN )
269         LDA = MAX( N, 1 )
270         XTYPE = 'N'
271         NIMAT = NTYPES
272         IF( N.LE.0 )
273     $      NIMAT = 1
274*
275         DO 80 IMAT = 1, NIMAT
276*
277*           Do the tests only if DOTYPE( IMAT ) is true.
278*
279            IF( .NOT.DOTYPE( IMAT ) )
280     $         GO TO 80
281*
282*           Skip types 5, 6, or 7 if the matrix size is too small.
283*
284            ZEROT = IMAT.GE.5 .AND. IMAT.LE.7
285            IF( ZEROT .AND. N.LT.IMAT-4 )
286     $         GO TO 80
287*
288*           Set up parameters with DLATB4 and generate a test matrix
289*           with DLATMS.
290*
291            CALL DLATB4( PATH, IMAT, N, N, TYPE, KL, KU, ANORM, MODE,
292     $                   CNDNUM, DIST )
293            RCONDC = ONE / CNDNUM
294*
295            SRNAMT = 'DLATMS'
296            CALL DLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE, CNDNUM,
297     $                   ANORM, KL, KU, 'No packing', A, LDA, WORK,
298     $                   INFO )
299*
300*           Check error code from DLATMS.
301*
302            IF( INFO.NE.0 ) THEN
303               CALL ALAERH( PATH, 'DLATMS', INFO, 0, ' ', N, N, -1, -1,
304     $                      -1, IMAT, NFAIL, NERRS, NOUT )
305               GO TO 80
306            END IF
307*
308*           For types 5-7, zero one or more columns of the matrix to
309*           test that INFO is returned correctly.
310*
311            IF( ZEROT ) THEN
312               IF( IMAT.EQ.5 ) THEN
313                  IZERO = 1
314               ELSE IF( IMAT.EQ.6 ) THEN
315                  IZERO = N
316               ELSE
317                  IZERO = N / 2 + 1
318               END IF
319               IOFF = ( IZERO-1 )*LDA
320               IF( IMAT.LT.7 ) THEN
321                  DO 20 I = 1, N
322                     A( IOFF+I ) = ZERO
323   20             CONTINUE
324               ELSE
325                  CALL DLASET( 'Full', N, N-IZERO+1, ZERO, ZERO,
326     $                         A( IOFF+1 ), LDA )
327               END IF
328            ELSE
329               IZERO = 0
330            END IF
331*
332*           Save a copy of the matrix A in ASAV.
333*
334            CALL DLACPY( 'Full', N, N, A, LDA, ASAV, LDA )
335*
336            DO 70 IEQUED = 1, 4
337               EQUED = EQUEDS( IEQUED )
338               IF( IEQUED.EQ.1 ) THEN
339                  NFACT = 3
340               ELSE
341                  NFACT = 1
342               END IF
343*
344               DO 60 IFACT = 1, NFACT
345                  FACT = FACTS( IFACT )
346                  PREFAC = LSAME( FACT, 'F' )
347                  NOFACT = LSAME( FACT, 'N' )
348                  EQUIL = LSAME( FACT, 'E' )
349*
350                  IF( ZEROT ) THEN
351                     IF( PREFAC )
352     $                  GO TO 60
353                     RCONDO = ZERO
354                     RCONDI = ZERO
355*
356                  ELSE IF( .NOT.NOFACT ) THEN
357*
358*                    Compute the condition number for comparison with
359*                    the value returned by DGESVX (FACT = 'N' reuses
360*                    the condition number from the previous iteration
361*                    with FACT = 'F').
362*
363                     CALL DLACPY( 'Full', N, N, ASAV, LDA, AFAC, LDA )
364                     IF( EQUIL .OR. IEQUED.GT.1 ) THEN
365*
366*                       Compute row and column scale factors to
367*                       equilibrate the matrix A.
368*
369                        CALL DGEEQU( N, N, AFAC, LDA, S, S( N+1 ),
370     $                               ROWCND, COLCND, AMAX, INFO )
371                        IF( INFO.EQ.0 .AND. N.GT.0 ) THEN
372                           IF( LSAME( EQUED, 'R' ) ) THEN
373                              ROWCND = ZERO
374                              COLCND = ONE
375                           ELSE IF( LSAME( EQUED, 'C' ) ) THEN
376                              ROWCND = ONE
377                              COLCND = ZERO
378                           ELSE IF( LSAME( EQUED, 'B' ) ) THEN
379                              ROWCND = ZERO
380                              COLCND = ZERO
381                           END IF
382*
383*                          Equilibrate the matrix.
384*
385                           CALL DLAQGE( N, N, AFAC, LDA, S, S( N+1 ),
386     $                                  ROWCND, COLCND, AMAX, EQUED )
387                        END IF
388                     END IF
389*
390*                    Save the condition number of the non-equilibrated
391*                    system for use in DGET04.
392*
393                     IF( EQUIL ) THEN
394                        ROLDO = RCONDO
395                        ROLDI = RCONDI
396                     END IF
397*
398*                    Compute the 1-norm and infinity-norm of A.
399*
400                     ANORMO = DLANGE( '1', N, N, AFAC, LDA, RWORK )
401                     ANORMI = DLANGE( 'I', N, N, AFAC, LDA, RWORK )
402*
403*                    Factor the matrix A.
404*
405                     SRNAMT = 'DGETRF'
406                     CALL DGETRF( N, N, AFAC, LDA, IWORK, INFO )
407*
408*                    Form the inverse of A.
409*
410                     CALL DLACPY( 'Full', N, N, AFAC, LDA, A, LDA )
411                     LWORK = NMAX*MAX( 3, NRHS )
412                     SRNAMT = 'DGETRI'
413                     CALL DGETRI( N, A, LDA, IWORK, WORK, LWORK, INFO )
414*
415*                    Compute the 1-norm condition number of A.
416*
417                     AINVNM = DLANGE( '1', N, N, A, LDA, RWORK )
418                     IF( ANORMO.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
419                        RCONDO = ONE
420                     ELSE
421                        RCONDO = ( ONE / ANORMO ) / AINVNM
422                     END IF
423*
424*                    Compute the infinity-norm condition number of A.
425*
426                     AINVNM = DLANGE( 'I', N, N, A, LDA, RWORK )
427                     IF( ANORMI.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
428                        RCONDI = ONE
429                     ELSE
430                        RCONDI = ( ONE / ANORMI ) / AINVNM
431                     END IF
432                  END IF
433*
434                  DO 50 ITRAN = 1, NTRAN
435*
436*                    Do for each value of TRANS.
437*
438                     TRANS = TRANSS( ITRAN )
439                     IF( ITRAN.EQ.1 ) THEN
440                        RCONDC = RCONDO
441                     ELSE
442                        RCONDC = RCONDI
443                     END IF
444*
445*                    Restore the matrix A.
446*
447                     CALL DLACPY( 'Full', N, N, ASAV, LDA, A, LDA )
448*
449*                    Form an exact solution and set the right hand side.
450*
451                     SRNAMT = 'DLARHS'
452                     CALL DLARHS( PATH, XTYPE, 'Full', TRANS, N, N, KL,
453     $                            KU, NRHS, A, LDA, XACT, LDA, B, LDA,
454     $                            ISEED, INFO )
455                     XTYPE = 'C'
456                     CALL DLACPY( 'Full', N, NRHS, B, LDA, BSAV, LDA )
457*
458                     IF( NOFACT .AND. ITRAN.EQ.1 ) THEN
459*
460*                       --- Test DGESV  ---
461*
462*                       Compute the LU factorization of the matrix and
463*                       solve the system.
464*
465                        CALL DLACPY( 'Full', N, N, A, LDA, AFAC, LDA )
466                        CALL DLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
467*
468                        SRNAMT = 'DGESV '
469                        CALL DGESV( N, NRHS, AFAC, LDA, IWORK, X, LDA,
470     $                              INFO )
471*
472*                       Check error code from DGESV .
473*
474                        IF( INFO.NE.IZERO )
475     $                     CALL ALAERH( PATH, 'DGESV ', INFO, IZERO,
476     $                                  ' ', N, N, -1, -1, NRHS, IMAT,
477     $                                  NFAIL, NERRS, NOUT )
478*
479*                       Reconstruct matrix from factors and compute
480*                       residual.
481*
482                        CALL DGET01( N, N, A, LDA, AFAC, LDA, IWORK,
483     $                               RWORK, RESULT( 1 ) )
484                        NT = 1
485                        IF( IZERO.EQ.0 ) THEN
486*
487*                          Compute residual of the computed solution.
488*
489                           CALL DLACPY( 'Full', N, NRHS, B, LDA, WORK,
490     $                                  LDA )
491                           CALL DGET02( 'No transpose', N, N, NRHS, A,
492     $                                  LDA, X, LDA, WORK, LDA, RWORK,
493     $                                  RESULT( 2 ) )
494*
495*                          Check solution from generated exact solution.
496*
497                           CALL DGET04( N, NRHS, X, LDA, XACT, LDA,
498     $                                  RCONDC, RESULT( 3 ) )
499                           NT = 3
500                        END IF
501*
502*                       Print information about the tests that did not
503*                       pass the threshold.
504*
505                        DO 30 K = 1, NT
506                           IF( RESULT( K ).GE.THRESH ) THEN
507                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
508     $                           CALL ALADHD( NOUT, PATH )
509                              WRITE( NOUT, FMT = 9999 )'DGESV ', N,
510     $                           IMAT, K, RESULT( K )
511                              NFAIL = NFAIL + 1
512                           END IF
513   30                   CONTINUE
514                        NRUN = NRUN + NT
515                     END IF
516*
517*                    --- Test DGESVX ---
518*
519                     IF( .NOT.PREFAC )
520     $                  CALL DLASET( 'Full', N, N, ZERO, ZERO, AFAC,
521     $                               LDA )
522                     CALL DLASET( 'Full', N, NRHS, ZERO, ZERO, X, LDA )
523                     IF( IEQUED.GT.1 .AND. N.GT.0 ) THEN
524*
525*                       Equilibrate the matrix if FACT = 'F' and
526*                       EQUED = 'R', 'C', or 'B'.
527*
528                        CALL DLAQGE( N, N, A, LDA, S, S( N+1 ), ROWCND,
529     $                               COLCND, AMAX, EQUED )
530                     END IF
531*
532*                    Solve the system and compute the condition number
533*                    and error bounds using DGESVX.
534*
535                     SRNAMT = 'DGESVX'
536                     CALL DGESVX( FACT, TRANS, N, NRHS, A, LDA, AFAC,
537     $                            LDA, IWORK, EQUED, S, S( N+1 ), B,
538     $                            LDA, X, LDA, RCOND, RWORK,
539     $                            RWORK( NRHS+1 ), WORK, IWORK( N+1 ),
540     $                            INFO )
541*
542*                    Check the error code from DGESVX.
543*
544                     IF( INFO.NE.IZERO )
545     $                  CALL ALAERH( PATH, 'DGESVX', INFO, IZERO,
546     $                               FACT // TRANS, N, N, -1, -1, NRHS,
547     $                               IMAT, NFAIL, NERRS, NOUT )
548*
549*                    Compare WORK(1) from DGESVX with the computed
550*                    reciprocal pivot growth factor RPVGRW
551*
552                     IF( INFO.NE.0 .AND. INFO.LE.N) THEN
553                        RPVGRW = DLANTR( 'M', 'U', 'N', INFO, INFO,
554     $                           AFAC, LDA, WORK )
555                        IF( RPVGRW.EQ.ZERO ) THEN
556                           RPVGRW = ONE
557                        ELSE
558                           RPVGRW = DLANGE( 'M', N, INFO, A, LDA,
559     $                              WORK ) / RPVGRW
560                        END IF
561                     ELSE
562                        RPVGRW = DLANTR( 'M', 'U', 'N', N, N, AFAC, LDA,
563     $                           WORK )
564                        IF( RPVGRW.EQ.ZERO ) THEN
565                           RPVGRW = ONE
566                        ELSE
567                           RPVGRW = DLANGE( 'M', N, N, A, LDA, WORK ) /
568     $                              RPVGRW
569                        END IF
570                     END IF
571                     RESULT( 7 ) = ABS( RPVGRW-WORK( 1 ) ) /
572     $                             MAX( WORK( 1 ), RPVGRW ) /
573     $                             DLAMCH( 'E' )
574*
575                     IF( .NOT.PREFAC ) THEN
576*
577*                       Reconstruct matrix from factors and compute
578*                       residual.
579*
580                        CALL DGET01( N, N, A, LDA, AFAC, LDA, IWORK,
581     $                               RWORK( 2*NRHS+1 ), RESULT( 1 ) )
582                        K1 = 1
583                     ELSE
584                        K1 = 2
585                     END IF
586*
587                     IF( INFO.EQ.0 ) THEN
588                        TRFCON = .FALSE.
589*
590*                       Compute residual of the computed solution.
591*
592                        CALL DLACPY( 'Full', N, NRHS, BSAV, LDA, WORK,
593     $                               LDA )
594                        CALL DGET02( TRANS, N, N, NRHS, ASAV, LDA, X,
595     $                               LDA, WORK, LDA, RWORK( 2*NRHS+1 ),
596     $                               RESULT( 2 ) )
597*
598*                       Check solution from generated exact solution.
599*
600                        IF( NOFACT .OR. ( PREFAC .AND. LSAME( EQUED,
601     $                      'N' ) ) ) THEN
602                           CALL DGET04( N, NRHS, X, LDA, XACT, LDA,
603     $                                  RCONDC, RESULT( 3 ) )
604                        ELSE
605                           IF( ITRAN.EQ.1 ) THEN
606                              ROLDC = ROLDO
607                           ELSE
608                              ROLDC = ROLDI
609                           END IF
610                           CALL DGET04( N, NRHS, X, LDA, XACT, LDA,
611     $                                  ROLDC, RESULT( 3 ) )
612                        END IF
613*
614*                       Check the error bounds from iterative
615*                       refinement.
616*
617                        CALL DGET07( TRANS, N, NRHS, ASAV, LDA, B, LDA,
618     $                               X, LDA, XACT, LDA, RWORK, .TRUE.,
619     $                               RWORK( NRHS+1 ), RESULT( 4 ) )
620                     ELSE
621                        TRFCON = .TRUE.
622                     END IF
623*
624*                    Compare RCOND from DGESVX with the computed value
625*                    in RCONDC.
626*
627                     RESULT( 6 ) = DGET06( RCOND, RCONDC )
628*
629*                    Print information about the tests that did not pass
630*                    the threshold.
631*
632                     IF( .NOT.TRFCON ) THEN
633                        DO 40 K = K1, NTESTS
634                           IF( RESULT( K ).GE.THRESH ) THEN
635                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
636     $                           CALL ALADHD( NOUT, PATH )
637                              IF( PREFAC ) THEN
638                                 WRITE( NOUT, FMT = 9997 )'DGESVX',
639     $                              FACT, TRANS, N, EQUED, IMAT, K,
640     $                              RESULT( K )
641                              ELSE
642                                 WRITE( NOUT, FMT = 9998 )'DGESVX',
643     $                              FACT, TRANS, N, IMAT, K, RESULT( K )
644                              END IF
645                              NFAIL = NFAIL + 1
646                           END IF
647   40                   CONTINUE
648                        NRUN = NRUN + NTESTS - K1 + 1
649                     ELSE
650                        IF( RESULT( 1 ).GE.THRESH .AND. .NOT.PREFAC )
651     $                       THEN
652                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
653     $                        CALL ALADHD( NOUT, PATH )
654                           IF( PREFAC ) THEN
655                              WRITE( NOUT, FMT = 9997 )'DGESVX', FACT,
656     $                           TRANS, N, EQUED, IMAT, 1, RESULT( 1 )
657                           ELSE
658                              WRITE( NOUT, FMT = 9998 )'DGESVX', FACT,
659     $                           TRANS, N, IMAT, 1, RESULT( 1 )
660                           END IF
661                           NFAIL = NFAIL + 1
662                           NRUN = NRUN + 1
663                        END IF
664                        IF( RESULT( 6 ).GE.THRESH ) THEN
665                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
666     $                        CALL ALADHD( NOUT, PATH )
667                           IF( PREFAC ) THEN
668                              WRITE( NOUT, FMT = 9997 )'DGESVX', FACT,
669     $                           TRANS, N, EQUED, IMAT, 6, RESULT( 6 )
670                           ELSE
671                              WRITE( NOUT, FMT = 9998 )'DGESVX', FACT,
672     $                           TRANS, N, IMAT, 6, RESULT( 6 )
673                           END IF
674                           NFAIL = NFAIL + 1
675                           NRUN = NRUN + 1
676                        END IF
677                        IF( RESULT( 7 ).GE.THRESH ) THEN
678                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
679     $                        CALL ALADHD( NOUT, PATH )
680                           IF( PREFAC ) THEN
681                              WRITE( NOUT, FMT = 9997 )'DGESVX', FACT,
682     $                           TRANS, N, EQUED, IMAT, 7, RESULT( 7 )
683                           ELSE
684                              WRITE( NOUT, FMT = 9998 )'DGESVX', FACT,
685     $                           TRANS, N, IMAT, 7, RESULT( 7 )
686                           END IF
687                           NFAIL = NFAIL + 1
688                           NRUN = NRUN + 1
689                        END IF
690*
691                     END IF
692*
693   50             CONTINUE
694   60          CONTINUE
695   70       CONTINUE
696   80    CONTINUE
697   90 CONTINUE
698*
699*     Print a summary of the results.
700*
701      CALL ALASVM( PATH, NOUT, NFAIL, NRUN, NERRS )
702*
703 9999 FORMAT( 1X, A, ', N =', I5, ', type ', I2, ', test(', I2, ') =',
704     $      G12.5 )
705 9998 FORMAT( 1X, A, ', FACT=''', A1, ''', TRANS=''', A1, ''', N=', I5,
706     $      ', type ', I2, ', test(', I1, ')=', G12.5 )
707 9997 FORMAT( 1X, A, ', FACT=''', A1, ''', TRANS=''', A1, ''', N=', I5,
708     $      ', EQUED=''', A1, ''', type ', I2, ', test(', I1, ')=',
709     $      G12.5 )
710      RETURN
711*
712*     End of DDRVGE
713*
714      END
715