File: pchrddriver.f

package info (click to toggle)
scalapack 2.2.2-1
  • links: PTS, VCS
  • area: main
  • in suites: sid, trixie
  • size: 37,012 kB
  • sloc: fortran: 339,113; ansic: 74,517; makefile: 1,494; sh: 34
file content (506 lines) | stat: -rw-r--r-- 19,229 bytes parent folder | download | duplicates (8)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
      PROGRAM PCHRDDRIVER
*
*  -- ScaLAPACK testing driver (version 1.7) --
*     University of Tennessee, Knoxville, Oak Ridge National Laboratory,
*     and University of California, Berkeley.
*     March 13, 2000
*
*  Purpose
*  =======
*
*  PCHRDDRIVER is the main test program for the COMPLEX
*  ScaLAPACK HRD (Hessenberg Reduction) routines.
*
*  The program must be driven by a short data file.  An annotated
*  example of a data file can be obtained by deleting the first 3
*  characters from the following 14 lines:
*  'ScaLAPACK HRD input file'
*  'PVM machine'
*  'HRD.out'            output file name (if any)
*  6                    device out
*  2                    number of problems sizes
*  100 101              values of N
*  2   1                values of ILO
*  99  101              values of IHI
*  3                    number of NB's
*  2 3 5                values of NB
*  7                    number of process grids (ordered pairs of P & Q)
*  1 2 1 4 2 3 8        values of P
*  1 2 4 1 3 2 1        values of Q
*  3.0                  threshold
*
*  Internal Parameters
*  ===================
*
*  TOTMEM   INTEGER, default = 2000000
*           TOTMEM is a machine-specific parameter indicating the
*           maximum amount of available memory in bytes.
*           The user should customize TOTMEM to his platform.  Remember
*           to leave room in memory for the operating system, the BLACS
*           buffer, etc.  For example, on a system with 8 MB of memory
*           per process (e.g., one processor on an Intel iPSC/860), the
*           parameters we use are TOTMEM=6200000 (leaving 1.8 MB for OS,
*           code, BLACS buffer, etc).  However, for PVM, we usually set
*           TOTMEM = 2000000.  Some experimenting with the maximum value
*           of TOTMEM may be required.
*
*  INTGSZ   INTEGER, default = 4 bytes.
*  CPLXSZ   INTEGER, default = 8 bytes.
*           INTGSZ and CPLXSZ indicate the length in bytes on the
*           given platform for an integer and a single precision
*           complex.
*  MEM      COMPLEX array, dimension ( TOTMEM / CPLXSZ )
*
*           All arrays used by SCALAPACK routines are allocated from
*           this array and referenced by pointers.  The integer IPA,
*           for example, is a pointer to the starting element of MEM for
*           the matrix A.
*
*  =====================================================================
*
*     .. Parameters ..
      INTEGER            BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
     $                   LLD_, MB_, M_, NB_, N_, RSRC_
      PARAMETER          ( BLOCK_CYCLIC_2D = 1, DLEN_ = 9, DTYPE_ = 1,
     $                     CTXT_ = 2, M_ = 3, N_ = 4, MB_ = 5, NB_ = 6,
     $                     RSRC_ = 7, CSRC_ = 8, LLD_ = 9 )
      INTEGER            CPLXSZ, MEMSIZ, NTESTS, TOTMEM
      COMPLEX            PADVAL
      PARAMETER          ( CPLXSZ = 8, TOTMEM = 2000000,
     $                     MEMSIZ = TOTMEM / CPLXSZ, NTESTS = 20,
     $                     PADVAL = ( -9923.0E+0, -9923.0E+0 ) )
*     ..
*     .. Local Scalars ..
      LOGICAL            CHECK
      CHARACTER*6        PASSED
      CHARACTER*80       OUTFILE
      INTEGER            I, IAM, IASEED, ICTXT, IHI, IHIP, IHLP, IHLQ,
     $                   ILCOL, ILO, ILROW, INFO, INLQ, IMIDPAD, IPA,
     $                   IPT, IPW, IPOSTPAD, IPREPAD, ITEMP, J, K,
     $                   KFAIL, KPASS, KSKIP, KTESTS, LCM, LCMQ, LOFF,
     $                   LWORK, MYCOL, MYROW, N, NB, NGRIDS, NMAT, NNB,
     $                   NPROCS, NOUT, NP, NPCOL, NPROW, NQ, WORKHRD,
     $                   WORKSIZ
      REAL               ANORM, FRESID, THRESH
      DOUBLE PRECISION   NOPS, TMFLOPS
*     ..
*     .. Local Arrays ..
      INTEGER            DESCA( DLEN_ ), IERR( 1 ), NBVAL( NTESTS ),
     $                   NVAL( NTESTS ), NVHI( NTESTS ), NVLO( NTESTS ),
     $                   PVAL( NTESTS ), QVAL( NTESTS )
      DOUBLE PRECISION   CTIME( 1 ), WTIME( 1 )
      COMPLEX            MEM( MEMSIZ )
*     ..
*     .. External Subroutines ..
      EXTERNAL           BLACS_BARRIER, BLACS_EXIT, BLACS_GET,
     $                   BLACS_GRIDEXIT, BLACS_GRIDINIT, BLACS_GRIDINFO,
     $                   DESCINIT, IGSUM2D, BLACS_PINFO, PCFILLPAD,
     $                   PCLAFCHK, PCGEHDRV, PCGEHRD,
     $                   PCHRDINFO, PCMATGEN, SLBOOT,
     $                   SLCOMBINE, SLTIMER
*     ..
*     .. External Functions ..
      INTEGER            ILCM, INDXG2P, NUMROC
      REAL               PCLANGE
      EXTERNAL           ILCM, INDXG2P, NUMROC, PCLANGE
*     ..
*     .. Intrinsic Functions ..
      INTRINSIC          DBLE, MAX
*     ..
*     .. Data statements ..
      DATA               KTESTS, KPASS, KFAIL, KSKIP / 4*0 /
*     ..
*     .. Executable Statements ..
*
*     Get starting information
*
      CALL BLACS_PINFO( IAM, NPROCS )
      IASEED = 100
      CALL PCHRDINFO( OUTFILE, NOUT, NMAT, NVAL, NVLO, NVHI, NTESTS,
     $                NNB, NBVAL, NTESTS, NGRIDS, PVAL, NTESTS, QVAL,
     $                NTESTS, THRESH, MEM, IAM, NPROCS )
      CHECK = ( THRESH.GE.0.0E+0 )
*
*     Print headings
*
      IF( IAM.EQ.0 ) THEN
         WRITE( NOUT, FMT = * )
         WRITE( NOUT, FMT = 9995 )
         WRITE( NOUT, FMT = 9994 )
         WRITE( NOUT, FMT = * )
      END IF
*
*     Loop over different process grids
*
      DO 30 I = 1, NGRIDS
*
         NPROW = PVAL( I )
         NPCOL = QVAL( I )
*
*        Make sure grid information is correct
*
         IERR( 1 ) = 0
         IF( NPROW.LT.1 ) THEN
            IF( IAM.EQ.0 )
     $         WRITE( NOUT, FMT = 9999 ) 'GRID', 'nprow', NPROW
            IERR( 1 ) = 1
         ELSE IF( NPCOL.LT.1 ) THEN
            IF( IAM.EQ.0 )
     $         WRITE( NOUT, FMT = 9999 ) 'GRID', 'npcol', NPCOL
            IERR( 1 ) = 1
         ELSE IF( NPROW*NPCOL.GT.NPROCS ) THEN
            IF( IAM.EQ.0 )
     $         WRITE( NOUT, FMT = 9998 )NPROW*NPCOL, NPROCS
            IERR( 1 ) = 1
         END IF
*
         IF( IERR( 1 ).GT.0 ) THEN
            IF( IAM.EQ.0 )
     $         WRITE( NOUT, FMT = 9997 ) 'grid'
            KSKIP = KSKIP + 1
            GO TO 30
         END IF
*
*        Define process grid
*
         CALL BLACS_GET( -1, 0, ICTXT )
         CALL BLACS_GRIDINIT( ICTXT, 'Row-major', NPROW, NPCOL )
         CALL BLACS_GRIDINFO( ICTXT, NPROW, NPCOL, MYROW, MYCOL )
*
*        Go to bottom of loop if this case doesn't use my process
*
         IF( MYROW.GE.NPROW .OR. MYCOL.GE.NPCOL )
     $      GOTO 30
*
         DO 20 J = 1, NMAT
*
            N = NVAL( J )
            ILO = NVLO( J )
            IHI = NVHI( J )
*
*           Make sure matrix information is correct
*
            IERR( 1 ) = 0
            IF( N.LT.1 ) THEN
               IF( IAM.EQ.0 )
     $            WRITE( NOUT, FMT = 9999 ) 'MATRIX', 'N', N
               IERR( 1 ) = 1
            END IF
*
*           Check all processes for an error
*
            CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
            IF( IERR( 1 ).GT.0 ) THEN
               IF( IAM.EQ.0 )
     $            WRITE( NOUT, FMT = 9997 ) 'matrix'
               KSKIP = KSKIP + 1
               GO TO 20
            END IF
*
            DO 10 K = 1, NNB
               NB = NBVAL( K )
*
*              Make sure nb is legal
*
               IERR( 1 ) = 0
               IF( NB.LT.1 ) THEN
                  IERR( 1 ) = 1
                  IF( IAM.EQ.0 )
     $               WRITE( NOUT, FMT = 9999 ) 'NB', 'NB', NB
               END IF
*
*              Check all processes for an error
*
               CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
               IF( IERR( 1 ).GT.0 ) THEN
                  IF( IAM.EQ.0 )
     $               WRITE( NOUT, FMT = 9997 ) 'NB'
                  KSKIP = KSKIP + 1
                  GO TO 10
               END IF
*
               NP = NUMROC( N, NB, MYROW, 0, NPROW )
               NQ = NUMROC( N, NB, MYCOL, 0, NPCOL )
               IF( CHECK ) THEN
                  IPREPAD  = MAX( NB, NP )
                  IMIDPAD  = NB
                  IPOSTPAD = MAX( NB, NQ )
               ELSE
                  IPREPAD  = 0
                  IMIDPAD  = 0
                  IPOSTPAD = 0
               END IF
*
*              Initialize the array descriptor for the matrix A
*
               CALL DESCINIT( DESCA, N, N, NB, NB, 0, 0, ICTXT,
     $                        MAX( 1, NP ) + IMIDPAD, INFO )
*
*              Check all processes for an error
*
               CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
               IF( IERR( 1 ).LT.0 ) THEN
                  IF( IAM.EQ.0 )
     $               WRITE( NOUT, FMT = 9997 ) 'descriptor'
                  KSKIP = KSKIP + 1
                  GO TO 10
               END IF
*
*              Assign pointers into MEM for SCALAPACK arrays, A is
*              allocated starting at position MEM( IPREPAD+1 )
*
               IPA = IPREPAD + 1
               IPT = IPA + DESCA( LLD_ )*NQ + IPOSTPAD + IPREPAD
               IPW = IPT + NQ + IPOSTPAD + IPREPAD
*
*              Calculate the amount of workspace required for the
*              reduction
*
               IHIP = NUMROC( IHI, NB, MYROW, DESCA( RSRC_ ), NPROW )
               LOFF = MOD( ILO-1, NB )
               ILROW = INDXG2P( ILO, NB, MYROW, DESCA( RSRC_ ), NPROW )
               ILCOL = INDXG2P( ILO, NB, MYCOL, DESCA( CSRC_ ), NPCOL )
               IHLP = NUMROC( IHI-ILO+LOFF+1, NB, MYROW, ILROW, NPROW )
               INLQ = NUMROC( N-ILO+LOFF+1, NB, MYCOL, ILCOL, NPCOL )
               LWORK = NB*( NB + MAX( IHIP+1, IHLP+INLQ ) )
               WORKHRD = LWORK + IPOSTPAD
               WORKSIZ = WORKHRD
*
*              Figure the amount of workspace required by the check
*
               IF( CHECK ) THEN
                  LCM = ILCM( NPROW, NPCOL )
                  LCMQ = LCM / NPCOL
                  IHLQ = NUMROC( IHI-ILO+LOFF+1, NB, MYCOL, ILCOL,
     $                           NPCOL )
                  ITEMP = NB*MAX( IHLP+INLQ, IHLQ+MAX( IHIP,
     $                    IHLP+NUMROC( NUMROC( IHI-ILO+LOFF+1, NB, 0, 0,
     $                    NPCOL ), NB, 0, 0, LCMQ ) ) )
                  WORKSIZ = MAX( NB*NB + NB*IHLP + ITEMP, NB * NP ) +
     $                           IPOSTPAD
               END IF
*
*              Check for adequate memory for problem size
*
               IERR( 1 ) = 0
               IF( IPW+WORKSIZ.GT.MEMSIZ ) THEN
                  IF( IAM.EQ.0 )
     $               WRITE( NOUT, FMT = 9996 ) 'Hessenberg reduction',
     $                      ( IPW+WORKSIZ )*CPLXSZ
                  IERR( 1 ) = 1
               END IF
*
*              Check all processes for an error
*
               CALL IGSUM2D( ICTXT, 'All', ' ', 1, 1, IERR, 1, -1, 0 )
*
               IF( IERR( 1 ).GT.0 ) THEN
                  IF( IAM.EQ.0 )
     $               WRITE( NOUT, FMT = 9997 ) 'MEMORY'
                  KSKIP = KSKIP + 1
                  GO TO 10
               END IF
*
*              Generate A
*
               CALL PCMATGEN( ICTXT, 'No', 'No', DESCA( M_ ),
     $                        DESCA( N_ ), DESCA( MB_ ), DESCA( NB_ ),
     $                        MEM( IPA ), DESCA( LLD_ ), DESCA( RSRC_ ),
     $                        DESCA( CSRC_ ),
     $                        IASEED, 0, NP, 0, NQ, MYROW, MYCOL,
     $                        NPROW, NPCOL )
*
*              Need Infinity-norm of A for checking
*
               IF( CHECK ) THEN
                  CALL PCFILLPAD( ICTXT, NP, NQ, MEM( IPA-IPREPAD ),
     $                            DESCA( LLD_ ), IPREPAD, IPOSTPAD,
     $                            PADVAL )
                  CALL PCFILLPAD( ICTXT, NQ, 1, MEM( IPT-IPREPAD ),
     $                            NQ, IPREPAD, IPOSTPAD, PADVAL )
                  CALL PCFILLPAD( ICTXT, WORKSIZ-IPOSTPAD, 1,
     $                            MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
     $                            IPREPAD, IPOSTPAD, PADVAL )
                  ANORM = PCLANGE( 'I', N, N, MEM( IPA ), 1, 1, DESCA,
     $                             MEM( IPW ) )
                  CALL PCCHEKPAD( ICTXT, 'PCLANGE', NP, NQ,
     $                            MEM( IPA-IPREPAD ), DESCA( LLD_ ),
     $                            IPREPAD, IPOSTPAD, PADVAL )
                  CALL PCCHEKPAD( ICTXT, 'PCLANGE',
     $                            WORKSIZ-IPOSTPAD, 1,
     $                            MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
     $                            IPREPAD, IPOSTPAD, PADVAL )
                  CALL PCFILLPAD( ICTXT, WORKHRD-IPOSTPAD, 1,
     $                            MEM( IPW-IPREPAD ), WORKHRD-IPOSTPAD,
     $                            IPREPAD, IPOSTPAD, PADVAL )
               END IF
*
               CALL SLBOOT()
               CALL BLACS_BARRIER( ICTXT, 'All' )
               CALL SLTIMER( 1 )
*
*              Reduce Hessenberg form
*
               CALL PCGEHRD( N, ILO, IHI, MEM( IPA ), 1, 1, DESCA,
     $                       MEM( IPT ), MEM( IPW ), LWORK, INFO )
               CALL SLTIMER( 1 )
*
               IF( CHECK ) THEN
*
*                 Check for memory overwrite
*
                  CALL PCCHEKPAD( ICTXT, 'PCGEHRD', NP, NQ,
     $                            MEM( IPA-IPREPAD ), DESCA( LLD_ ),
     $                            IPREPAD, IPOSTPAD, PADVAL )
                  CALL PCCHEKPAD( ICTXT, 'PCGEHRD', NQ, 1,
     $                            MEM( IPT-IPREPAD ), NQ, IPREPAD,
     $                            IPOSTPAD, PADVAL )
                  CALL PCCHEKPAD( ICTXT, 'PCGEHRD', WORKHRD-IPOSTPAD,
     $                            1, MEM( IPW-IPREPAD ),
     $                            WORKHRD-IPOSTPAD, IPREPAD,
     $                            IPOSTPAD, PADVAL )
                  CALL PCFILLPAD( ICTXT, WORKSIZ-IPOSTPAD, 1,
     $                            MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
     $                            IPREPAD, IPOSTPAD, PADVAL )
*
*                 Compute fctres = ||A - Q H Q'|| / (||A||*N*eps)
*
                  CALL PCGEHDRV( N, ILO, IHI, MEM( IPA ), 1, 1, DESCA,
     $                           MEM( IPT ), MEM( IPW ) )
                  CALL PCLAFCHK( 'No', 'No', N, N, MEM( IPA ), 1, 1,
     $                           DESCA, IASEED, ANORM, FRESID,
     $                           MEM( IPW ) )
*
*                 Check for memory overwrite
*
                  CALL PCCHEKPAD( ICTXT, 'PCGEHDRV', NP, NQ,
     $                            MEM( IPA-IPREPAD ), DESCA( LLD_ ),
     $                            IPREPAD, IPOSTPAD, PADVAL )
                  CALL PCCHEKPAD( ICTXT, 'PCGEHDRV', NQ, 1,
     $                            MEM( IPT-IPREPAD ), NQ, IPREPAD,
     $                            IPOSTPAD, PADVAL )
                  CALL PCCHEKPAD( ICTXT, 'PCGEHDRV',
     $                            WORKSIZ-IPOSTPAD, 1,
     $                            MEM( IPW-IPREPAD ), WORKSIZ-IPOSTPAD,
     $                            IPREPAD, IPOSTPAD, PADVAL )
*
*                 Test residual and detect NaN result
*
                  IF( FRESID.LE.THRESH .AND. FRESID-FRESID.EQ.0.0E+0 )
     $                THEN
                     KPASS = KPASS + 1
                     PASSED = 'PASSED'
                  ELSE
                     IF( MYROW.EQ.0 .AND. MYCOL.EQ.0 )
     $                  WRITE( NOUT, FMT = 9986 ) FRESID
                     KFAIL = KFAIL + 1
                     PASSED = 'FAILED'
                  END IF
               ELSE
*
*                 Don't perform the checking, only the timing operation
*
                  KPASS = KPASS + 1
                  FRESID = FRESID - FRESID
                  PASSED = 'BYPASS'
               END IF
*
*              Gather max. of all CPU and WALL clock timings
*
               CALL SLCOMBINE( ICTXT, 'All', '>', 'W', 1, 1, WTIME )
               CALL SLCOMBINE( ICTXT, 'All', '>', 'C', 1, 1, CTIME )
*
*              Print results
*
               IF( MYROW.EQ.0 .AND. MYCOL.EQ.0 ) THEN
*
*                 HRD requires 40/3 * N^3 floating point ops. (flops)
*                 more precisely,
*                 HRD requires 16/3*(IHI-ILO)^3+8*IHI*(IHI-ILO)^2 flops
*
                  NOPS = DBLE( IHI-ILO )
                  NOPS = NOPS * NOPS *
     $                   ( 8.0D0*DBLE( IHI ) + (16.0D0/3.0D0)*NOPS )
                  NOPS = NOPS / 1.0D+6
*
*                 Print WALL time
*
                  IF( WTIME( 1 ).GT.0.0D+0 ) THEN
                     TMFLOPS = NOPS / WTIME( 1 )
                  ELSE
                     TMFLOPS = 0.0D+0
                  END IF
                  IF( WTIME( 1 ).GE.0.0D+0 )
     $               WRITE( NOUT, FMT = 9993 ) 'WALL',  N, ILO, IHI, NB,
     $                      NPROW, NPCOL, WTIME( 1 ), TMFLOPS, FRESID,
     $                      PASSED
*
*                 Print CPU time
*
                  IF( CTIME( 1 ).GT.0.0D+0 ) THEN
                     TMFLOPS = NOPS / CTIME( 1 )
                  ELSE
                     TMFLOPS = 0.0D+0
                  END IF
                  IF( CTIME( 1 ).GE.0.0D+0 )
     $               WRITE( NOUT, FMT = 9993 ) 'CPU ', N, ILO, IHI, NB,
     $                      NPROW, NPCOL, CTIME( 1 ), TMFLOPS, FRESID,
     $                      PASSED
               END IF
   10       CONTINUE
   20    CONTINUE
*
         CALL BLACS_GRIDEXIT( ICTXT )
   30 CONTINUE
*
*     Print ending messages and close output file
*
      IF( IAM.EQ.0 ) THEN
         KTESTS = KPASS + KFAIL + KSKIP
         WRITE( NOUT, FMT = * )
         WRITE( NOUT, FMT = 9992 ) KTESTS
         IF( CHECK ) THEN
            WRITE( NOUT, FMT = 9991 ) KPASS
            WRITE( NOUT, FMT = 9989 ) KFAIL
         ELSE
            WRITE( NOUT, FMT = 9990 ) KPASS
         END IF
         WRITE( NOUT, FMT = 9988 ) KSKIP
         WRITE( NOUT, FMT = * )
         WRITE( NOUT, FMT = * )
         WRITE( NOUT, FMT = 9987 )
         IF( NOUT.NE.6 .AND. NOUT.NE.0 )
     $      CLOSE( NOUT )
      END IF
*
      CALL BLACS_EXIT( 0 )
*
 9999 FORMAT( 'ILLEGAL ', A6, ': ', A5, ' = ', I3,
     $      '; It should be at least 1' )
 9998 FORMAT( 'ILLEGAL GRID: nprow*npcol = ', I4, '. It can be at most',
     $      I4 )
 9997 FORMAT( 'Bad ', A6, ' parameters: going on to next test case.' )
 9996 FORMAT( 'Unable to perform ', A, ': need TOTMEM of at least',
     $      I11 )
 9995 FORMAT( 'TIME      N    ILO    IHI  NB     P     Q  HRD Time ',
     $        '     MFLOPS Residual  CHECK' )
 9994 FORMAT( '---- ------ ------ ------ --- ----- ----- --------- ',
     $        '----------- -------- ------' )
 9993 FORMAT( A4, 1X, I6, 1X, I6, 1X, I6, 1X, I3, 1X, I5, 1X, I5, 1X,
     $        F9.2, 1X, F11.2, 1X, F8.2, 1X, A6 )
 9992 FORMAT( 'Finished', I4, ' tests, with the following results:' )
 9991 FORMAT( I5, ' tests completed and passed residual checks.' )
 9990 FORMAT( I5, ' tests completed without checking.' )
 9989 FORMAT( I5, ' tests completed and failed residual checks.' )
 9988 FORMAT( I5, ' tests skipped because of illegal input values.' )
 9987 FORMAT( 'END OF TESTS.' )
 9986 FORMAT( '||A - Q*H*Q''|| / (||A|| * N * eps) = ', G25.7 )
*
      STOP
*
*     End of PCHRDDRIVER
*
      END