1*> \brief \b DCHKAA
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*  Definition:
9*  ===========
10*
11*       PROGRAM DCHKAA
12*
13*
14*> \par Purpose:
15*  =============
16*>
17*> \verbatim
18*>
19*> DCHKAA is the main test program for the DOUBLE PRECISION LAPACK
20*> linear equation routines
21*>
22*> The program must be driven by a short data file. The first 15 records
23*> (not including the first comment  line) specify problem dimensions
24*> and program options using list-directed input. The remaining lines
25*> specify the LAPACK test paths and the number of matrix types to use
26*> in testing.  An annotated example of a data file can be obtained by
27*> deleting the first 3 characters from the following 40 lines:
28*> Data file for testing DOUBLE PRECISION LAPACK linear eqn. routines
29*> 7                      Number of values of M
30*> 0 1 2 3 5 10 16        Values of M (row dimension)
31*> 7                      Number of values of N
32*> 0 1 2 3 5 10 16        Values of N (column dimension)
33*> 1                      Number of values of NRHS
34*> 2                      Values of NRHS (number of right hand sides)
35*> 5                      Number of values of NB
36*> 1 3 3 3 20             Values of NB (the blocksize)
37*> 1 0 5 9 1              Values of NX (crossover point)
38*> 3                      Number of values of RANK
39*> 30 50 90               Values of rank (as a % of N)
40*> 20.0                   Threshold value of test ratio
41*> T                      Put T to test the LAPACK routines
42*> T                      Put T to test the driver routines
43*> T                      Put T to test the error exits
44*> DGE   11               List types on next line if 0 < NTYPES < 11
45*> DGB    8               List types on next line if 0 < NTYPES <  8
46*> DGT   12               List types on next line if 0 < NTYPES < 12
47*> DPO    9               List types on next line if 0 < NTYPES <  9
48*> DPS    9               List types on next line if 0 < NTYPES <  9
49*> DPP    9               List types on next line if 0 < NTYPES <  9
50*> DPB    8               List types on next line if 0 < NTYPES <  8
51*> DPT   12               List types on next line if 0 < NTYPES < 12
52*> DSY   10               List types on next line if 0 < NTYPES < 10
53*> DSR   10               List types on next line if 0 < NTYPES < 10
54*> DSP   10               List types on next line if 0 < NTYPES < 10
55*> DTR   18               List types on next line if 0 < NTYPES < 18
56*> DTP   18               List types on next line if 0 < NTYPES < 18
57*> DTB   17               List types on next line if 0 < NTYPES < 17
58*> DQR    8               List types on next line if 0 < NTYPES <  8
59*> DRQ    8               List types on next line if 0 < NTYPES <  8
60*> DLQ    8               List types on next line if 0 < NTYPES <  8
61*> DQL    8               List types on next line if 0 < NTYPES <  8
62*> DQP    6               List types on next line if 0 < NTYPES <  6
63*> DTZ    3               List types on next line if 0 < NTYPES <  3
64*> DLS    6               List types on next line if 0 < NTYPES <  6
65*> DEQ
66*> DQT
67*> DQX
68*> \endverbatim
69*
70*  Parameters:
71*  ==========
72*
73*> \verbatim
74*>  NMAX    INTEGER
75*>          The maximum allowable value for M and N.
76*>
77*>  MAXIN   INTEGER
78*>          The number of different values that can be used for each of
79*>          M, N, NRHS, NB, NX and RANK
80*>
81*>  MAXRHS  INTEGER
82*>          The maximum number of right hand sides
83*>
84*>  MATMAX  INTEGER
85*>          The maximum number of matrix types to use for testing
86*>
87*>  NIN     INTEGER
88*>          The unit number for input
89*>
90*>  NOUT    INTEGER
91*>          The unit number for output
92*> \endverbatim
93*
94*  Authors:
95*  ========
96*
97*> \author Univ. of Tennessee
98*> \author Univ. of California Berkeley
99*> \author Univ. of Colorado Denver
100*> \author NAG Ltd.
101*
102*> \date April 2012
103*
104*> \ingroup double_lin
105*
106*  =====================================================================
107      PROGRAM DCHKAA
108*
109*  -- LAPACK test routine (version 3.4.1) --
110*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
111*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
112*     April 2012
113*
114*  =====================================================================
115*
116*     .. Parameters ..
117      INTEGER            NMAX
118      PARAMETER          ( NMAX = 132 )
119      INTEGER            MAXIN
120      PARAMETER          ( MAXIN = 12 )
121      INTEGER            MAXRHS
122      PARAMETER          ( MAXRHS = 16 )
123      INTEGER            MATMAX
124      PARAMETER          ( MATMAX = 30 )
125      INTEGER            NIN, NOUT
126      PARAMETER          ( NIN = 5, NOUT = 6 )
127      INTEGER            KDMAX
128      PARAMETER          ( KDMAX = NMAX+( NMAX+1 ) / 4 )
129*     ..
130*     .. Local Scalars ..
131      LOGICAL            FATAL, TSTCHK, TSTDRV, TSTERR
132      CHARACTER          C1
133      CHARACTER*2        C2
134      CHARACTER*3        PATH
135      CHARACTER*10       INTSTR
136      CHARACTER*72       ALINE
137      INTEGER            I, IC, J, K, LA, LAFAC, LDA, NB, NM, NMATS, NN,
138     $                   NNB, NNB2, NNS, NRHS, NTYPES, NRANK,
139     $                   VERS_MAJOR, VERS_MINOR, VERS_PATCH
140      DOUBLE PRECISION   EPS, S1, S2, THREQ, THRESH
141*     ..
142*     .. Local Arrays ..
143      LOGICAL            DOTYPE( MATMAX )
144      INTEGER            IWORK( 25*NMAX ), MVAL( MAXIN ),
145     $                   NBVAL( MAXIN ), NBVAL2( MAXIN ),
146     $                   NSVAL( MAXIN ), NVAL( MAXIN ), NXVAL( MAXIN ),
147     $                   RANKVAL( MAXIN ), PIV( NMAX )
148      DOUBLE PRECISION   A( ( KDMAX+1 )*NMAX, 7 ), B( NMAX*MAXRHS, 4 ),
149     $                   RWORK( 5*NMAX+2*MAXRHS ), S( 2*NMAX ),
150     $                   WORK( NMAX, NMAX+MAXRHS+30 )
151*     ..
152*     .. External Functions ..
153      LOGICAL            LSAME, LSAMEN
154      DOUBLE PRECISION   DLAMCH, DSECND
155      EXTERNAL           LSAME, LSAMEN, DLAMCH, DSECND
156*     ..
157*     .. External Subroutines ..
158      EXTERNAL           ALAREQ, DCHKEQ, DCHKGB, DCHKGE, DCHKGT, DCHKLQ,
159     $                   DCHKPB, DCHKPO, DCHKPS, DCHKPP, DCHKPT, DCHKQ3,
160     $                   DCHKQL, DCHKQP, DCHKQR, DCHKRQ, DCHKSP, DCHKSY,
161     $                   DCHKSY_ROOK, DCHKTB, DCHKTP, DCHKTR, DCHKTZ,
162     $                   DDRVGB, DDRVGE, DDRVGT, DDRVLS, DDRVPB, DDRVPO,
163     $                   DDRVPP, DDRVPT, DDRVSP, DDRVSY, DDRVSY_ROOK,
164     $                   ILAVER, DCHKQRT, DCHKQRTP
165*     ..
166*     .. Scalars in Common ..
167      LOGICAL            LERR, OK
168      CHARACTER*32       SRNAMT
169      INTEGER            INFOT, NUNIT
170*     ..
171*     .. Arrays in Common ..
172      INTEGER            IPARMS( 100 )
173*     ..
174*     .. Common blocks ..
175      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
176      COMMON             / SRNAMC / SRNAMT
177      COMMON             / CLAENV / IPARMS
178*     ..
179*     .. Data statements ..
180      DATA               THREQ / 2.0D0 / , INTSTR / '0123456789' /
181*     ..
182*     .. Executable Statements ..
183*
184#ifdef __FLAME__
185      CALL FLA_INIT
186#endif
187      S1 = DSECND( )
188      LDA = NMAX
189      FATAL = .FALSE.
190*
191*     Read a dummy line.
192*
193      READ( NIN, FMT = * )
194*
195*     Report values of parameters.
196*
197      CALL ILAVER( VERS_MAJOR, VERS_MINOR, VERS_PATCH )
198      WRITE( NOUT, FMT = 9994 ) VERS_MAJOR, VERS_MINOR, VERS_PATCH
199*
200*     Read the values of M
201*
202      READ( NIN, FMT = * )NM
203      IF( NM.LT.1 ) THEN
204         WRITE( NOUT, FMT = 9996 )' NM ', NM, 1
205         NM = 0
206         FATAL = .TRUE.
207      ELSE IF( NM.GT.MAXIN ) THEN
208         WRITE( NOUT, FMT = 9995 )' NM ', NM, MAXIN
209         NM = 0
210         FATAL = .TRUE.
211      END IF
212      READ( NIN, FMT = * )( MVAL( I ), I = 1, NM )
213      DO 10 I = 1, NM
214         IF( MVAL( I ).LT.0 ) THEN
215            WRITE( NOUT, FMT = 9996 )' M  ', MVAL( I ), 0
216            FATAL = .TRUE.
217         ELSE IF( MVAL( I ).GT.NMAX ) THEN
218            WRITE( NOUT, FMT = 9995 )' M  ', MVAL( I ), NMAX
219            FATAL = .TRUE.
220         END IF
221   10 CONTINUE
222      IF( NM.GT.0 )
223     $   WRITE( NOUT, FMT = 9993 )'M   ', ( MVAL( I ), I = 1, NM )
224*
225*     Read the values of N
226*
227      READ( NIN, FMT = * )NN
228      IF( NN.LT.1 ) THEN
229         WRITE( NOUT, FMT = 9996 )' NN ', NN, 1
230         NN = 0
231         FATAL = .TRUE.
232      ELSE IF( NN.GT.MAXIN ) THEN
233         WRITE( NOUT, FMT = 9995 )' NN ', NN, MAXIN
234         NN = 0
235         FATAL = .TRUE.
236      END IF
237      READ( NIN, FMT = * )( NVAL( I ), I = 1, NN )
238      DO 20 I = 1, NN
239         IF( NVAL( I ).LT.0 ) THEN
240            WRITE( NOUT, FMT = 9996 )' N  ', NVAL( I ), 0
241            FATAL = .TRUE.
242         ELSE IF( NVAL( I ).GT.NMAX ) THEN
243            WRITE( NOUT, FMT = 9995 )' N  ', NVAL( I ), NMAX
244            FATAL = .TRUE.
245         END IF
246   20 CONTINUE
247      IF( NN.GT.0 )
248     $   WRITE( NOUT, FMT = 9993 )'N   ', ( NVAL( I ), I = 1, NN )
249*
250*     Read the values of NRHS
251*
252      READ( NIN, FMT = * )NNS
253      IF( NNS.LT.1 ) THEN
254         WRITE( NOUT, FMT = 9996 )' NNS', NNS, 1
255         NNS = 0
256         FATAL = .TRUE.
257      ELSE IF( NNS.GT.MAXIN ) THEN
258         WRITE( NOUT, FMT = 9995 )' NNS', NNS, MAXIN
259         NNS = 0
260         FATAL = .TRUE.
261      END IF
262      READ( NIN, FMT = * )( NSVAL( I ), I = 1, NNS )
263      DO 30 I = 1, NNS
264         IF( NSVAL( I ).LT.0 ) THEN
265            WRITE( NOUT, FMT = 9996 )'NRHS', NSVAL( I ), 0
266            FATAL = .TRUE.
267         ELSE IF( NSVAL( I ).GT.MAXRHS ) THEN
268            WRITE( NOUT, FMT = 9995 )'NRHS', NSVAL( I ), MAXRHS
269            FATAL = .TRUE.
270         END IF
271   30 CONTINUE
272      IF( NNS.GT.0 )
273     $   WRITE( NOUT, FMT = 9993 )'NRHS', ( NSVAL( I ), I = 1, NNS )
274*
275*     Read the values of NB
276*
277      READ( NIN, FMT = * )NNB
278      IF( NNB.LT.1 ) THEN
279         WRITE( NOUT, FMT = 9996 )'NNB ', NNB, 1
280         NNB = 0
281         FATAL = .TRUE.
282      ELSE IF( NNB.GT.MAXIN ) THEN
283         WRITE( NOUT, FMT = 9995 )'NNB ', NNB, MAXIN
284         NNB = 0
285         FATAL = .TRUE.
286      END IF
287      READ( NIN, FMT = * )( NBVAL( I ), I = 1, NNB )
288      DO 40 I = 1, NNB
289         IF( NBVAL( I ).LT.0 ) THEN
290            WRITE( NOUT, FMT = 9996 )' NB ', NBVAL( I ), 0
291            FATAL = .TRUE.
292         END IF
293   40 CONTINUE
294      IF( NNB.GT.0 )
295     $   WRITE( NOUT, FMT = 9993 )'NB  ', ( NBVAL( I ), I = 1, NNB )
296*
297*     Set NBVAL2 to be the set of unique values of NB
298*
299      NNB2 = 0
300      DO 60 I = 1, NNB
301         NB = NBVAL( I )
302         DO 50 J = 1, NNB2
303            IF( NB.EQ.NBVAL2( J ) )
304     $         GO TO 60
305   50    CONTINUE
306         NNB2 = NNB2 + 1
307         NBVAL2( NNB2 ) = NB
308   60 CONTINUE
309*
310*     Read the values of NX
311*
312      READ( NIN, FMT = * )( NXVAL( I ), I = 1, NNB )
313      DO 70 I = 1, NNB
314         IF( NXVAL( I ).LT.0 ) THEN
315            WRITE( NOUT, FMT = 9996 )' NX ', NXVAL( I ), 0
316            FATAL = .TRUE.
317         END IF
318   70 CONTINUE
319      IF( NNB.GT.0 )
320     $   WRITE( NOUT, FMT = 9993 )'NX  ', ( NXVAL( I ), I = 1, NNB )
321*
322*     Read the values of RANKVAL
323*
324      READ( NIN, FMT = * )NRANK
325      IF( NN.LT.1 ) THEN
326         WRITE( NOUT, FMT = 9996 )' NRANK ', NRANK, 1
327         NRANK = 0
328         FATAL = .TRUE.
329      ELSE IF( NN.GT.MAXIN ) THEN
330         WRITE( NOUT, FMT = 9995 )' NRANK ', NRANK, MAXIN
331         NRANK = 0
332         FATAL = .TRUE.
333      END IF
334      READ( NIN, FMT = * )( RANKVAL( I ), I = 1, NRANK )
335      DO I = 1, NRANK
336         IF( RANKVAL( I ).LT.0 ) THEN
337            WRITE( NOUT, FMT = 9996 )' RANK  ', RANKVAL( I ), 0
338            FATAL = .TRUE.
339         ELSE IF( RANKVAL( I ).GT.100 ) THEN
340            WRITE( NOUT, FMT = 9995 )' RANK  ', RANKVAL( I ), 100
341            FATAL = .TRUE.
342         END IF
343      END DO
344      IF( NRANK.GT.0 )
345     $   WRITE( NOUT, FMT = 9993 )'RANK % OF N',
346     $   ( RANKVAL( I ), I = 1, NRANK )
347*
348*     Read the threshold value for the test ratios.
349*
350      READ( NIN, FMT = * )THRESH
351      WRITE( NOUT, FMT = 9992 )THRESH
352*
353*     Read the flag that indicates whether to test the LAPACK routines.
354*
355      READ( NIN, FMT = * )TSTCHK
356*
357*     Read the flag that indicates whether to test the driver routines.
358*
359      READ( NIN, FMT = * )TSTDRV
360*
361*     Read the flag that indicates whether to test the error exits.
362*
363      READ( NIN, FMT = * )TSTERR
364*
365      IF( FATAL ) THEN
366         WRITE( NOUT, FMT = 9999 )
367         STOP
368      END IF
369*
370*     Calculate and print the machine dependent constants.
371*
372      EPS = DLAMCH( 'Underflow threshold' )
373      WRITE( NOUT, FMT = 9991 )'underflow', EPS
374      EPS = DLAMCH( 'Overflow threshold' )
375      WRITE( NOUT, FMT = 9991 )'overflow ', EPS
376      EPS = DLAMCH( 'Epsilon' )
377      WRITE( NOUT, FMT = 9991 )'precision', EPS
378      WRITE( NOUT, FMT = * )
379*
380   80 CONTINUE
381*
382*     Read a test path and the number of matrix types to use.
383*
384      READ( NIN, FMT = '(A72)', END = 140 )ALINE
385      PATH = ALINE( 1: 3 )
386      NMATS = MATMAX
387      I = 3
388   90 CONTINUE
389      I = I + 1
390      IF( I.GT.72 ) THEN
391         NMATS = MATMAX
392         GO TO 130
393      END IF
394      IF( ALINE( I: I ).EQ.' ' )
395     $   GO TO 90
396      NMATS = 0
397  100 CONTINUE
398      C1 = ALINE( I: I )
399      DO 110 K = 1, 10
400         IF( C1.EQ.INTSTR( K: K ) ) THEN
401            IC = K - 1
402            GO TO 120
403         END IF
404  110 CONTINUE
405      GO TO 130
406  120 CONTINUE
407      NMATS = NMATS*10 + IC
408      I = I + 1
409      IF( I.GT.72 )
410     $   GO TO 130
411      GO TO 100
412  130 CONTINUE
413      C1 = PATH( 1: 1 )
414      C2 = PATH( 2: 3 )
415      NRHS = NSVAL( 1 )
416*
417*     Check first character for correct precision.
418*
419      IF( .NOT.LSAME( C1, 'Double precision' ) ) THEN
420         WRITE( NOUT, FMT = 9990 )PATH
421*
422      ELSE IF( NMATS.LE.0 ) THEN
423*
424*        Check for a positive number of tests requested.
425*
426         WRITE( NOUT, FMT = 9989 )PATH
427*
428      ELSE IF( LSAMEN( 2, C2, 'GE' ) ) THEN
429*
430*        GE:  general matrices
431*
432         NTYPES = 11
433         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
434*
435         IF( TSTCHK ) THEN
436            CALL DCHKGE( DOTYPE, NM, MVAL, NN, NVAL, NNB2, NBVAL2, NNS,
437     $                   NSVAL, THRESH, TSTERR, LDA, A( 1, 1 ),
438     $                   A( 1, 2 ), A( 1, 3 ), B( 1, 1 ), B( 1, 2 ),
439     $                   B( 1, 3 ), WORK, RWORK, IWORK, NOUT )
440         ELSE
441            WRITE( NOUT, FMT = 9989 )PATH
442         END IF
443*
444         IF( TSTDRV ) THEN
445            CALL DDRVGE( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, LDA,
446     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
447     $                   B( 1, 2 ), B( 1, 3 ), B( 1, 4 ), S, WORK,
448     $                   RWORK, IWORK, NOUT )
449         ELSE
450            WRITE( NOUT, FMT = 9988 )PATH
451         END IF
452*
453      ELSE IF( LSAMEN( 2, C2, 'GB' ) ) THEN
454*
455*        GB:  general banded matrices
456*
457         LA = ( 2*KDMAX+1 )*NMAX
458         LAFAC = ( 3*KDMAX+1 )*NMAX
459         NTYPES = 8
460         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
461*
462         IF( TSTCHK ) THEN
463            CALL DCHKGB( DOTYPE, NM, MVAL, NN, NVAL, NNB2, NBVAL2, NNS,
464     $                   NSVAL, THRESH, TSTERR, A( 1, 1 ), LA,
465     $                   A( 1, 3 ), LAFAC, B( 1, 1 ), B( 1, 2 ),
466     $                   B( 1, 3 ), WORK, RWORK, IWORK, NOUT )
467         ELSE
468            WRITE( NOUT, FMT = 9989 )PATH
469         END IF
470*
471         IF( TSTDRV ) THEN
472            CALL DDRVGB( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
473     $                   A( 1, 1 ), LA, A( 1, 3 ), LAFAC, A( 1, 6 ),
474     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), B( 1, 4 ), S,
475     $                   WORK, RWORK, IWORK, NOUT )
476         ELSE
477            WRITE( NOUT, FMT = 9988 )PATH
478         END IF
479*
480      ELSE IF( LSAMEN( 2, C2, 'GT' ) ) THEN
481*
482*        GT:  general tridiagonal matrices
483*
484         NTYPES = 12
485         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
486*
487         IF( TSTCHK ) THEN
488            CALL DCHKGT( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
489     $                   A( 1, 1 ), A( 1, 2 ), B( 1, 1 ), B( 1, 2 ),
490     $                   B( 1, 3 ), WORK, RWORK, IWORK, NOUT )
491         ELSE
492            WRITE( NOUT, FMT = 9989 )PATH
493         END IF
494*
495         IF( TSTDRV ) THEN
496            CALL DDRVGT( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
497     $                   A( 1, 1 ), A( 1, 2 ), B( 1, 1 ), B( 1, 2 ),
498     $                   B( 1, 3 ), WORK, RWORK, IWORK, NOUT )
499         ELSE
500            WRITE( NOUT, FMT = 9988 )PATH
501         END IF
502*
503      ELSE IF( LSAMEN( 2, C2, 'PO' ) ) THEN
504*
505*        PO:  positive definite matrices
506*
507         NTYPES = 9
508         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
509*
510         IF( TSTCHK ) THEN
511            CALL DCHKPO( DOTYPE, NN, NVAL, NNB2, NBVAL2, NNS, NSVAL,
512     $                   THRESH, TSTERR, LDA, A( 1, 1 ), A( 1, 2 ),
513     $                   A( 1, 3 ), B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
514     $                   WORK, RWORK, IWORK, NOUT )
515         ELSE
516            WRITE( NOUT, FMT = 9989 )PATH
517         END IF
518*
519         IF( TSTDRV ) THEN
520            CALL DDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, LDA,
521     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
522     $                   B( 1, 2 ), B( 1, 3 ), B( 1, 4 ), S, WORK,
523     $                   RWORK, IWORK, NOUT )
524         ELSE
525            WRITE( NOUT, FMT = 9988 )PATH
526         END IF
527*
528      ELSE IF( LSAMEN( 2, C2, 'PS' ) ) THEN
529*
530*        PS:  positive semi-definite matrices
531*
532         NTYPES = 9
533*
534         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
535*
536         IF( TSTCHK ) THEN
537            CALL DCHKPS( DOTYPE, NN, NVAL, NNB2, NBVAL2, NRANK,
538     $                   RANKVAL, THRESH, TSTERR, LDA, A( 1, 1 ),
539     $                   A( 1, 2 ), A( 1, 3 ), PIV, WORK, RWORK,
540     $                   NOUT )
541         ELSE
542            WRITE( NOUT, FMT = 9989 )PATH
543         END IF
544*
545      ELSE IF( LSAMEN( 2, C2, 'PP' ) ) THEN
546*
547*        PP:  positive definite packed matrices
548*
549         NTYPES = 9
550         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
551*
552         IF( TSTCHK ) THEN
553            CALL DCHKPP( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
554     $                   LDA, A( 1, 1 ), A( 1, 2 ), A( 1, 3 ),
555     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), WORK, RWORK,
556     $                   IWORK, NOUT )
557         ELSE
558            WRITE( NOUT, FMT = 9989 )PATH
559         END IF
560*
561         IF( TSTDRV ) THEN
562            CALL DDRVPP( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, LDA,
563     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
564     $                   B( 1, 2 ), B( 1, 3 ), B( 1, 4 ), S, WORK,
565     $                   RWORK, IWORK, NOUT )
566         ELSE
567            WRITE( NOUT, FMT = 9988 )PATH
568         END IF
569*
570      ELSE IF( LSAMEN( 2, C2, 'PB' ) ) THEN
571*
572*        PB:  positive definite banded matrices
573*
574         NTYPES = 8
575         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
576*
577         IF( TSTCHK ) THEN
578            CALL DCHKPB( DOTYPE, NN, NVAL, NNB2, NBVAL2, NNS, NSVAL,
579     $                   THRESH, TSTERR, LDA, A( 1, 1 ), A( 1, 2 ),
580     $                   A( 1, 3 ), B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
581     $                   WORK, RWORK, IWORK, NOUT )
582         ELSE
583            WRITE( NOUT, FMT = 9989 )PATH
584         END IF
585*
586         IF( TSTDRV ) THEN
587            CALL DDRVPB( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, LDA,
588     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
589     $                   B( 1, 2 ), B( 1, 3 ), B( 1, 4 ), S, WORK,
590     $                   RWORK, IWORK, NOUT )
591         ELSE
592            WRITE( NOUT, FMT = 9988 )PATH
593         END IF
594*
595      ELSE IF( LSAMEN( 2, C2, 'PT' ) ) THEN
596*
597*        PT:  positive definite tridiagonal matrices
598*
599         NTYPES = 12
600         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
601*
602         IF( TSTCHK ) THEN
603            CALL DCHKPT( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
604     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
605     $                   B( 1, 2 ), B( 1, 3 ), WORK, RWORK, NOUT )
606         ELSE
607            WRITE( NOUT, FMT = 9989 )PATH
608         END IF
609*
610         IF( TSTDRV ) THEN
611            CALL DDRVPT( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
612     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
613     $                   B( 1, 2 ), B( 1, 3 ), WORK, RWORK, NOUT )
614         ELSE
615            WRITE( NOUT, FMT = 9988 )PATH
616         END IF
617*
618      ELSE IF( LSAMEN( 2, C2, 'SY' ) ) THEN
619*
620*        SY:  symmetric indefinite matrices,
621*             with partial (Bunch-Kaufman) pivoting algorithm
622*
623         NTYPES = 10
624         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
625*
626         IF( TSTCHK ) THEN
627            CALL DCHKSY( DOTYPE, NN, NVAL, NNB2, NBVAL2, NNS, NSVAL,
628     $                   THRESH, TSTERR, LDA, A( 1, 1 ), A( 1, 2 ),
629     $                   A( 1, 3 ), B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
630     $                   WORK, RWORK, IWORK, NOUT )
631         ELSE
632            WRITE( NOUT, FMT = 9989 )PATH
633         END IF
634*
635         IF( TSTDRV ) THEN
636            CALL DDRVSY( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, LDA,
637     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
638     $                   B( 1, 2 ), B( 1, 3 ), WORK, RWORK, IWORK,
639     $                   NOUT )
640         ELSE
641            WRITE( NOUT, FMT = 9988 )PATH
642         END IF
643*
644      ELSE IF( LSAMEN( 2, C2, 'SR' ) ) THEN
645*
646*        SR:  symmetric indefinite matrices with Rook pivoting,
647*             with rook (bounded Bunch-Kaufman) pivoting algorithm
648*
649         NTYPES = 10
650         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
651*
652         IF( TSTCHK ) THEN
653            CALL DCHKSY_ROOK(DOTYPE, NN, NVAL, NNB2, NBVAL2, NNS, NSVAL,
654     $                       THRESH, TSTERR, LDA, A( 1, 1 ), A( 1, 2 ),
655     $                       A( 1, 3 ), B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
656     $                       WORK, RWORK, IWORK, NOUT )
657         ELSE
658            WRITE( NOUT, FMT = 9989 )PATH
659         END IF
660*
661         IF( TSTDRV ) THEN
662            CALL DDRVSY_ROOK( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
663     $                        LDA, A( 1, 1 ), A( 1, 2 ), A( 1, 3 ),
664     $                        B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
665     $                        WORK, RWORK, IWORK, NOUT )
666         ELSE
667            WRITE( NOUT, FMT = 9988 )PATH
668         END IF
669*
670      ELSE IF( LSAMEN( 2, C2, 'SP' ) ) THEN
671*
672*        SP:  symmetric indefinite packed matrices,
673*             with partial (Bunch-Kaufman) pivoting algorithm
674*
675         NTYPES = 10
676         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
677*
678         IF( TSTCHK ) THEN
679            CALL DCHKSP( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
680     $                   LDA, A( 1, 1 ), A( 1, 2 ), A( 1, 3 ),
681     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), WORK, RWORK,
682     $                   IWORK, NOUT )
683         ELSE
684            WRITE( NOUT, FMT = 9989 )PATH
685         END IF
686*
687         IF( TSTDRV ) THEN
688            CALL DDRVSP( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, LDA,
689     $                   A( 1, 1 ), A( 1, 2 ), A( 1, 3 ), B( 1, 1 ),
690     $                   B( 1, 2 ), B( 1, 3 ), WORK, RWORK, IWORK,
691     $                   NOUT )
692         ELSE
693            WRITE( NOUT, FMT = 9988 )PATH
694         END IF
695*
696      ELSE IF( LSAMEN( 2, C2, 'TR' ) ) THEN
697*
698*        TR:  triangular matrices
699*
700         NTYPES = 18
701         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
702*
703         IF( TSTCHK ) THEN
704            CALL DCHKTR( DOTYPE, NN, NVAL, NNB2, NBVAL2, NNS, NSVAL,
705     $                   THRESH, TSTERR, LDA, A( 1, 1 ), A( 1, 2 ),
706     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), WORK, RWORK,
707     $                   IWORK, NOUT )
708         ELSE
709            WRITE( NOUT, FMT = 9989 )PATH
710         END IF
711*
712      ELSE IF( LSAMEN( 2, C2, 'TP' ) ) THEN
713*
714*        TP:  triangular packed matrices
715*
716         NTYPES = 18
717         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
718*
719         IF( TSTCHK ) THEN
720            CALL DCHKTP( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
721     $                   LDA, A( 1, 1 ), A( 1, 2 ), B( 1, 1 ),
722     $                   B( 1, 2 ), B( 1, 3 ), WORK, RWORK, IWORK,
723     $                   NOUT )
724         ELSE
725            WRITE( NOUT, FMT = 9989 )PATH
726         END IF
727*
728      ELSE IF( LSAMEN( 2, C2, 'TB' ) ) THEN
729*
730*        TB:  triangular banded matrices
731*
732         NTYPES = 17
733         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
734*
735         IF( TSTCHK ) THEN
736            CALL DCHKTB( DOTYPE, NN, NVAL, NNS, NSVAL, THRESH, TSTERR,
737     $                   LDA, A( 1, 1 ), A( 1, 2 ), B( 1, 1 ),
738     $                   B( 1, 2 ), B( 1, 3 ), WORK, RWORK, IWORK,
739     $                   NOUT )
740         ELSE
741            WRITE( NOUT, FMT = 9989 )PATH
742         END IF
743*
744      ELSE IF( LSAMEN( 2, C2, 'QR' ) ) THEN
745*
746*        QR:  QR factorization
747*
748         NTYPES = 8
749         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
750*
751         IF( TSTCHK ) THEN
752            CALL DCHKQR( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
753     $                   NRHS, THRESH, TSTERR, NMAX, A( 1, 1 ),
754     $                   A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), A( 1, 5 ),
755     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), B( 1, 4 ),
756     $                   WORK, RWORK, IWORK, NOUT )
757         ELSE
758            WRITE( NOUT, FMT = 9989 )PATH
759         END IF
760*
761      ELSE IF( LSAMEN( 2, C2, 'LQ' ) ) THEN
762*
763*        LQ:  LQ factorization
764*
765         NTYPES = 8
766         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
767*
768         IF( TSTCHK ) THEN
769            CALL DCHKLQ( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
770     $                   NRHS, THRESH, TSTERR, NMAX, A( 1, 1 ),
771     $                   A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), A( 1, 5 ),
772     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), B( 1, 4 ),
773     $                   WORK, RWORK, NOUT )
774         ELSE
775            WRITE( NOUT, FMT = 9989 )PATH
776         END IF
777*
778      ELSE IF( LSAMEN( 2, C2, 'QL' ) ) THEN
779*
780*        QL:  QL factorization
781*
782         NTYPES = 8
783         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
784*
785         IF( TSTCHK ) THEN
786            CALL DCHKQL( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
787     $                   NRHS, THRESH, TSTERR, NMAX, A( 1, 1 ),
788     $                   A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), A( 1, 5 ),
789     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), B( 1, 4 ),
790     $                   WORK, RWORK, IWORK, NOUT )
791         ELSE
792            WRITE( NOUT, FMT = 9989 )PATH
793         END IF
794*
795      ELSE IF( LSAMEN( 2, C2, 'RQ' ) ) THEN
796*
797*        RQ:  RQ factorization
798*
799         NTYPES = 8
800         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
801*
802         IF( TSTCHK ) THEN
803            CALL DCHKRQ( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
804     $                   NRHS, THRESH, TSTERR, NMAX, A( 1, 1 ),
805     $                   A( 1, 2 ), A( 1, 3 ), A( 1, 4 ), A( 1, 5 ),
806     $                   B( 1, 1 ), B( 1, 2 ), B( 1, 3 ), B( 1, 4 ),
807     $                   WORK, RWORK, IWORK, NOUT )
808         ELSE
809            WRITE( NOUT, FMT = 9989 )PATH
810         END IF
811*
812      ELSE IF( LSAMEN( 2, C2, 'QP' ) ) THEN
813*
814*        QP:  QR factorization with pivoting
815*
816         NTYPES = 6
817         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
818*
819         IF( TSTCHK ) THEN
820            CALL DCHKQP( DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR,
821     $                   A( 1, 1 ), A( 1, 2 ), B( 1, 1 ),
822     $                   B( 1, 3 ), WORK, IWORK, NOUT )
823            CALL DCHKQ3( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
824     $                   THRESH, A( 1, 1 ), A( 1, 2 ), B( 1, 1 ),
825     $                   B( 1, 3 ), WORK, IWORK, NOUT )
826         ELSE
827            WRITE( NOUT, FMT = 9989 )PATH
828         END IF
829*
830      ELSE IF( LSAMEN( 2, C2, 'TZ' ) ) THEN
831*
832*        TZ:  Trapezoidal matrix
833*
834         NTYPES = 3
835         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
836*
837         IF( TSTCHK ) THEN
838            CALL DCHKTZ( DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR,
839     $                   A( 1, 1 ), A( 1, 2 ), B( 1, 1 ),
840     $                   B( 1, 3 ), WORK, NOUT )
841         ELSE
842            WRITE( NOUT, FMT = 9989 )PATH
843         END IF
844*
845      ELSE IF( LSAMEN( 2, C2, 'LS' ) ) THEN
846*
847*        LS:  Least squares drivers
848*
849         NTYPES = 6
850         CALL ALAREQ( PATH, NMATS, DOTYPE, NTYPES, NIN, NOUT )
851*
852         IF( TSTDRV ) THEN
853            CALL DDRVLS( DOTYPE, NM, MVAL, NN, NVAL, NNS, NSVAL, NNB,
854     $                   NBVAL, NXVAL, THRESH, TSTERR, A( 1, 1 ),
855     $                   A( 1, 2 ), B( 1, 1 ), B( 1, 2 ), B( 1, 3 ),
856     $                   RWORK, RWORK( NMAX+1 ), WORK, IWORK, NOUT )
857         ELSE
858            WRITE( NOUT, FMT = 9988 )PATH
859         END IF
860*
861      ELSE IF( LSAMEN( 2, C2, 'EQ' ) ) THEN
862*
863*        EQ:  Equilibration routines for general and positive definite
864*             matrices (THREQ should be between 2 and 10)
865*
866         IF( TSTCHK ) THEN
867            CALL DCHKEQ( THREQ, NOUT )
868         ELSE
869            WRITE( NOUT, FMT = 9989 )PATH
870         END IF
871*
872      ELSE IF( LSAMEN( 2, C2, 'QT' ) ) THEN
873*
874*        QT:  QRT routines for general matrices
875*
876         IF( TSTCHK ) THEN
877            CALL DCHKQRT( THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB,
878     $                    NBVAL, NOUT )
879         ELSE
880            WRITE( NOUT, FMT = 9989 )PATH
881         END IF
882*
883      ELSE IF( LSAMEN( 2, C2, 'QX' ) ) THEN
884*
885*        QX:  QRT routines for triangular-pentagonal matrices
886*
887         IF( TSTCHK ) THEN
888            CALL DCHKQRTP( THRESH, TSTERR, NM, MVAL, NN, NVAL, NNB,
889     $                     NBVAL, NOUT )
890         ELSE
891            WRITE( NOUT, FMT = 9989 )PATH
892         END IF
893*
894      ELSE
895*
896         WRITE( NOUT, FMT = 9990 )PATH
897      END IF
898*
899*     Go back to get another input line.
900*
901      GO TO 80
902*
903*     Branch to this line when the last record is read.
904*
905  140 CONTINUE
906      CLOSE ( NIN )
907      S2 = DSECND( )
908      WRITE( NOUT, FMT = 9998 )
909      WRITE( NOUT, FMT = 9997 )S2 - S1
910*
911 9999 FORMAT( / ' Execution not attempted due to input errors' )
912 9998 FORMAT( / ' End of tests' )
913 9997 FORMAT( ' Total time used = ', F12.2, ' seconds', / )
914 9996 FORMAT( ' Invalid input value: ', A4, '=', I6, '; must be >=',
915     $      I6 )
916 9995 FORMAT( ' Invalid input value: ', A4, '=', I6, '; must be <=',
917     $      I6 )
918 9994 FORMAT( ' Tests of the DOUBLE PRECISION LAPACK routines ',
919     $      / ' LAPACK VERSION ', I1, '.', I1, '.', I1,
920     $      / / ' The following parameter values will be used:' )
921 9993 FORMAT( 4X, A4, ':  ', 10I6, / 11X, 10I6 )
922 9992 FORMAT( / ' Routines pass computational tests if test ratio is ',
923     $      'less than', F8.2, / )
924 9991 FORMAT( ' Relative machine ', A, ' is taken to be', D16.6 )
925 9990 FORMAT( / 1X, A3, ':  Unrecognized path name' )
926 9989 FORMAT( / 1X, A3, ' routines were not tested' )
927 9988 FORMAT( / 1X, A3, ' driver routines were not tested' )
928*
929#ifdef __FLAME__
930      CALL FLA_FINALIZE
931#endif
932*
933*     End of DCHKAA
934*
935      END
936