1*> \brief \b ZCHKHE
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 ZCHKHE( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
12*                          THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
13*                          XACT, WORK, RWORK, IWORK, NOUT )
14*
15*       .. Scalar Arguments ..
16*       LOGICAL            TSTERR
17*       INTEGER            NMAX, NN, NNB, NNS, NOUT
18*       DOUBLE PRECISION   THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
23*       DOUBLE PRECISION   RWORK( * )
24*       COMPLEX*16         A( * ), AFAC( * ), AINV( * ), B( * ),
25*      $                   WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> ZCHKHE tests ZHETRF, -TRI2, -TRS, -TRS2, -RFS, and -CON.
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 dimension N.
58*> \endverbatim
59*>
60*> \param[in] NNB
61*> \verbatim
62*>          NNB is INTEGER
63*>          The number of values of NB contained in the vector NBVAL.
64*> \endverbatim
65*>
66*> \param[in] NBVAL
67*> \verbatim
68*>          NBVAL is INTEGER array, dimension (NNB)
69*>          The values of the blocksize NB.
70*> \endverbatim
71*>
72*> \param[in] NNS
73*> \verbatim
74*>          NNS is INTEGER
75*>          The number of values of NRHS contained in the vector NSVAL.
76*> \endverbatim
77*>
78*> \param[in] NSVAL
79*> \verbatim
80*>          NSVAL is INTEGER array, dimension (NNS)
81*>          The values of the number of right hand sides NRHS.
82*> \endverbatim
83*>
84*> \param[in] THRESH
85*> \verbatim
86*>          THRESH is DOUBLE PRECISION
87*>          The threshold value for the test ratios.  A result is
88*>          included in the output file if RESULT >= THRESH.  To have
89*>          every test ratio printed, use THRESH = 0.
90*> \endverbatim
91*>
92*> \param[in] TSTERR
93*> \verbatim
94*>          TSTERR is LOGICAL
95*>          Flag that indicates whether error exits are to be tested.
96*> \endverbatim
97*>
98*> \param[in] NMAX
99*> \verbatim
100*>          NMAX is INTEGER
101*>          The maximum value permitted for N, used in dimensioning the
102*>          work arrays.
103*> \endverbatim
104*>
105*> \param[out] A
106*> \verbatim
107*>          A is COMPLEX*16 array, dimension (NMAX*NMAX)
108*> \endverbatim
109*>
110*> \param[out] AFAC
111*> \verbatim
112*>          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
113*> \endverbatim
114*>
115*> \param[out] AINV
116*> \verbatim
117*>          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
118*> \endverbatim
119*>
120*> \param[out] B
121*> \verbatim
122*>          B is COMPLEX*16 array, dimension (NMAX*NSMAX)
123*>          where NSMAX is the largest entry in NSVAL.
124*> \endverbatim
125*>
126*> \param[out] X
127*> \verbatim
128*>          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
129*> \endverbatim
130*>
131*> \param[out] XACT
132*> \verbatim
133*>          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
134*> \endverbatim
135*>
136*> \param[out] WORK
137*> \verbatim
138*>          WORK is COMPLEX*16 array, dimension (NMAX*max(3,NSMAX))
139*> \endverbatim
140*>
141*> \param[out] RWORK
142*> \verbatim
143*>          RWORK is DOUBLE PRECISION array, dimension (max(NMAX,2*NSMAX))
144*> \endverbatim
145*>
146*> \param[out] IWORK
147*> \verbatim
148*>          IWORK is INTEGER array, dimension (NMAX)
149*> \endverbatim
150*>
151*> \param[in] NOUT
152*> \verbatim
153*>          NOUT is INTEGER
154*>          The unit number for output.
155*> \endverbatim
156*
157*  Authors:
158*  ========
159*
160*> \author Univ. of Tennessee
161*> \author Univ. of California Berkeley
162*> \author Univ. of Colorado Denver
163*> \author NAG Ltd.
164*
165*> \ingroup complex16_lin
166*
167*  =====================================================================
168      SUBROUTINE ZCHKHE( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
169     $                   THRESH, TSTERR, NMAX, A, AFAC, AINV, B, X,
170     $                   XACT, WORK, RWORK, IWORK, NOUT )
171*
172*  -- LAPACK test routine --
173*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
174*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176*     .. Scalar Arguments ..
177      LOGICAL            TSTERR
178      INTEGER            NMAX, NN, NNB, NNS, NOUT
179      DOUBLE PRECISION   THRESH
180*     ..
181*     .. Array Arguments ..
182      LOGICAL            DOTYPE( * )
183      INTEGER            IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
184      DOUBLE PRECISION   RWORK( * )
185      COMPLEX*16         A( * ), AFAC( * ), AINV( * ), B( * ),
186     $                   WORK( * ), X( * ), XACT( * )
187*     ..
188*
189*  =====================================================================
190*
191*     .. Parameters ..
192      DOUBLE PRECISION   ZERO
193      PARAMETER          ( ZERO = 0.0D+0 )
194      COMPLEX*16         CZERO
195      PARAMETER          ( CZERO = ( 0.0D+0, 0.0D+0 ) )
196      INTEGER            NTYPES
197      PARAMETER          ( NTYPES = 10 )
198      INTEGER            NTESTS
199      PARAMETER          ( NTESTS = 9 )
200*     ..
201*     .. Local Scalars ..
202      LOGICAL            TRFCON, ZEROT
203      CHARACTER          DIST, TYPE, UPLO, XTYPE
204      CHARACTER*3        PATH
205      INTEGER            I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
206     $                   IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
207     $                   N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
208      DOUBLE PRECISION   ANORM, CNDNUM, RCOND, RCONDC
209*     ..
210*     .. Local Arrays ..
211      CHARACTER          UPLOS( 2 )
212      INTEGER            ISEED( 4 ), ISEEDY( 4 )
213      DOUBLE PRECISION   RESULT( NTESTS )
214*     ..
215*     .. External Functions ..
216      DOUBLE PRECISION   DGET06, ZLANHE
217      EXTERNAL           DGET06, ZLANHE
218*     ..
219*     .. External Subroutines ..
220      EXTERNAL           ALAERH, ALAHD, ALASUM, XLAENV, ZERRHE, ZGET04,
221     $                   ZHECON, ZHERFS, ZHET01, ZHETRF, ZHETRI2,
222     $                   ZHETRS, ZLACPY, ZLAIPD, ZLARHS, ZLATB4, ZLATMS,
223     $                   ZPOT02, ZPOT03, ZPOT05
224*     ..
225*     .. Intrinsic Functions ..
226      INTRINSIC          MAX, MIN
227*     ..
228*     .. Scalars in Common ..
229      LOGICAL            LERR, OK
230      CHARACTER*32       SRNAMT
231      INTEGER            INFOT, NUNIT
232*     ..
233*     .. Common blocks ..
234      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
235      COMMON             / SRNAMC / SRNAMT
236*     ..
237*     .. Data statements ..
238      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
239      DATA               UPLOS / 'U', 'L' /
240*     ..
241*     .. Executable Statements ..
242*
243*     Initialize constants and the random number seed.
244*
245      PATH( 1: 1 ) = 'Zomplex precision'
246      PATH( 2: 3 ) = 'HE'
247      NRUN = 0
248      NFAIL = 0
249      NERRS = 0
250      DO 10 I = 1, 4
251         ISEED( I ) = ISEEDY( I )
252   10 CONTINUE
253*
254*     Test the error exits
255*
256      IF( TSTERR )
257     $   CALL ZERRHE( PATH, NOUT )
258      INFOT = 0
259*
260*     Set the minimum block size for which the block routine should
261*     be used, which will be later returned by ILAENV
262*
263      CALL XLAENV( 2, 2 )
264*
265*     Do for each value of N in NVAL
266*
267      DO 180 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         IZERO = 0
276         DO 170 IMAT = 1, NIMAT
277*
278*           Do the tests only if DOTYPE( IMAT ) is true.
279*
280            IF( .NOT.DOTYPE( IMAT ) )
281     $         GO TO 170
282*
283*           Skip types 3, 4, 5, or 6 if the matrix size is too small.
284*
285            ZEROT = IMAT.GE.3 .AND. IMAT.LE.6
286            IF( ZEROT .AND. N.LT.IMAT-2 )
287     $         GO TO 170
288*
289*           Do first for UPLO = 'U', then for UPLO = 'L'
290*
291            DO 160 IUPLO = 1, 2
292               UPLO = UPLOS( IUPLO )
293*
294*              Set up parameters with ZLATB4 for the matrix generator
295*              based on the type of matrix to be generated.
296*
297               CALL ZLATB4( PATH, IMAT, N, N, TYPE, KL, KU, ANORM, MODE,
298     $                      CNDNUM, DIST )
299*
300*              Generate a matrix with ZLATMS.
301*
302               SRNAMT = 'ZLATMS'
303               CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
304     $                      CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
305     $                      INFO )
306*
307*              Check error code from ZLATMS and handle error.
308*
309               IF( INFO.NE.0 ) THEN
310                  CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N, -1,
311     $                         -1, -1, IMAT, NFAIL, NERRS, NOUT )
312*
313*                 Skip all tests for this generated matrix
314*
315                  GO TO 160
316               END IF
317*
318*              For types 3-6, zero one or more rows and columns of
319*              the matrix to test that INFO is returned correctly.
320*
321               IF( ZEROT ) THEN
322                  IF( IMAT.EQ.3 ) THEN
323                     IZERO = 1
324                  ELSE IF( IMAT.EQ.4 ) THEN
325                     IZERO = N
326                  ELSE
327                     IZERO = N / 2 + 1
328                  END IF
329*
330                  IF( IMAT.LT.6 ) THEN
331*
332*                    Set row and column IZERO to zero.
333*
334                     IF( IUPLO.EQ.1 ) THEN
335                        IOFF = ( IZERO-1 )*LDA
336                        DO 20 I = 1, IZERO - 1
337                           A( IOFF+I ) = CZERO
338   20                   CONTINUE
339                        IOFF = IOFF + IZERO
340                        DO 30 I = IZERO, N
341                           A( IOFF ) = CZERO
342                           IOFF = IOFF + LDA
343   30                   CONTINUE
344                     ELSE
345                        IOFF = IZERO
346                        DO 40 I = 1, IZERO - 1
347                           A( IOFF ) = CZERO
348                           IOFF = IOFF + LDA
349   40                   CONTINUE
350                        IOFF = IOFF - IZERO
351                        DO 50 I = IZERO, N
352                           A( IOFF+I ) = CZERO
353   50                   CONTINUE
354                     END IF
355                  ELSE
356                     IF( IUPLO.EQ.1 ) THEN
357*
358*                       Set the first IZERO rows and columns to zero.
359*
360                        IOFF = 0
361                        DO 70 J = 1, N
362                           I2 = MIN( J, IZERO )
363                           DO 60 I = 1, I2
364                              A( IOFF+I ) = CZERO
365   60                      CONTINUE
366                           IOFF = IOFF + LDA
367   70                   CONTINUE
368                     ELSE
369*
370*                       Set the last IZERO rows and columns to zero.
371*
372                        IOFF = 0
373                        DO 90 J = 1, N
374                           I1 = MAX( J, IZERO )
375                           DO 80 I = I1, N
376                              A( IOFF+I ) = CZERO
377   80                      CONTINUE
378                           IOFF = IOFF + LDA
379   90                   CONTINUE
380                     END IF
381                  END IF
382               ELSE
383                  IZERO = 0
384               END IF
385*
386*              End generate test matrix A.
387*
388*
389*              Set the imaginary part of the diagonals.
390*
391               CALL ZLAIPD( N, A, LDA+1, 0 )
392*
393*              Do for each value of NB in NBVAL
394*
395               DO 150 INB = 1, NNB
396*
397*                 Set the optimal blocksize, which will be later
398*                 returned by ILAENV.
399*
400                  NB = NBVAL( INB )
401                  CALL XLAENV( 1, NB )
402*
403*                 Copy the test matrix A into matrix AFAC which
404*                 will be factorized in place. This is needed to
405*                 preserve the test matrix A for subsequent tests.
406*
407                  CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
408*
409*                 Compute the L*D*L**T or U*D*U**T factorization of the
410*                 matrix. IWORK stores details of the interchanges and
411*                 the block structure of D. AINV is a work array for
412*                 block factorization, LWORK is the length of AINV.
413*
414                  LWORK = MAX( 2, NB )*LDA
415                  SRNAMT = 'ZHETRF'
416                  CALL ZHETRF( UPLO, N, AFAC, LDA, IWORK, AINV, LWORK,
417     $                         INFO )
418*
419*                 Adjust the expected value of INFO to account for
420*                 pivoting.
421*
422                  K = IZERO
423                  IF( K.GT.0 ) THEN
424  100                CONTINUE
425                     IF( IWORK( K ).LT.0 ) THEN
426                        IF( IWORK( K ).NE.-K ) THEN
427                           K = -IWORK( K )
428                           GO TO 100
429                        END IF
430                     ELSE IF( IWORK( K ).NE.K ) THEN
431                        K = IWORK( K )
432                        GO TO 100
433                     END IF
434                  END IF
435*
436*                 Check error code from ZHETRF and handle error.
437*
438                  IF( INFO.NE.K )
439     $               CALL ALAERH( PATH, 'ZHETRF', INFO, K, UPLO, N, N,
440     $                            -1, -1, NB, IMAT, NFAIL, NERRS, NOUT )
441*
442*                 Set the condition estimate flag if the INFO is not 0.
443*
444                  IF( INFO.NE.0 ) THEN
445                     TRFCON = .TRUE.
446                  ELSE
447                     TRFCON = .FALSE.
448                  END IF
449*
450*+    TEST 1
451*                 Reconstruct matrix from factors and compute residual.
452*
453                  CALL ZHET01( UPLO, N, A, LDA, AFAC, LDA, IWORK, AINV,
454     $                         LDA, RWORK, RESULT( 1 ) )
455                  NT = 1
456*
457*+    TEST 2
458*                 Form the inverse and compute the residual.
459*
460                  IF( INB.EQ.1 .AND. .NOT.TRFCON ) THEN
461                     CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
462                     SRNAMT = 'ZHETRI2'
463                     LWORK = (N+NB+1)*(NB+3)
464                     CALL ZHETRI2( UPLO, N, AINV, LDA, IWORK, WORK,
465     $                            LWORK, INFO )
466*
467*                    Check error code from ZHETRI and handle error.
468*
469                     IF( INFO.NE.0 )
470     $                  CALL ALAERH( PATH, 'ZHETRI', INFO, -1, UPLO, N,
471     $                               N, -1, -1, -1, IMAT, NFAIL, NERRS,
472     $                               NOUT )
473*
474*                    Compute the residual for a symmetric matrix times
475*                    its inverse.
476*
477                     CALL ZPOT03( UPLO, N, A, LDA, AINV, LDA, WORK, LDA,
478     $                            RWORK, RCONDC, RESULT( 2 ) )
479                     NT = 2
480                  END IF
481*
482*                 Print information about the tests that did not pass
483*                 the threshold.
484*
485                  DO 110 K = 1, NT
486                     IF( RESULT( K ).GE.THRESH ) THEN
487                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
488     $                     CALL ALAHD( NOUT, PATH )
489                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
490     $                     RESULT( K )
491                        NFAIL = NFAIL + 1
492                     END IF
493  110             CONTINUE
494                  NRUN = NRUN + NT
495*
496*                 Skip the other tests if this is not the first block
497*                 size.
498*
499                  IF( INB.GT.1 )
500     $               GO TO 150
501*
502*                 Do only the condition estimate if INFO is not 0.
503*
504                  IF( TRFCON ) THEN
505                     RCONDC = ZERO
506                     GO TO 140
507                  END IF
508*
509*                 Do for each value of NRHS in NSVAL.
510*
511                  DO 130 IRHS = 1, NNS
512                     NRHS = NSVAL( IRHS )
513*
514*+    TEST 3 (Using TRS)
515*                 Solve and compute residual for  A * X = B.
516*
517*                    Choose a set of NRHS random solution vectors
518*                    stored in XACT and set up the right hand side B
519*
520                     SRNAMT = 'ZLARHS'
521                     CALL ZLARHS( PATH, XTYPE, UPLO, ' ', N, N, KL, KU,
522     $                            NRHS, A, LDA, XACT, LDA, B, LDA,
523     $                            ISEED, INFO )
524                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
525*
526                     SRNAMT = 'ZHETRS'
527                     CALL ZHETRS( UPLO, N, NRHS, AFAC, LDA, IWORK, X,
528     $                            LDA, INFO )
529*
530*                    Check error code from ZHETRS and handle error.
531*
532                     IF( INFO.NE.0 )
533     $                  CALL ALAERH( PATH, 'ZHETRS', INFO, 0, UPLO, N,
534     $                               N, -1, -1, NRHS, IMAT, NFAIL,
535     $                               NERRS, NOUT )
536*
537                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK, LDA )
538*
539*                    Compute the residual for the solution
540*
541                     CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA, WORK,
542     $                            LDA, RWORK, RESULT( 3 ) )
543*
544*+    TEST 4 (Using TRS2)
545*                 Solve and compute residual for  A * X = B.
546*
547*                    Choose a set of NRHS random solution vectors
548*                    stored in XACT and set up the right hand side B
549*
550                     SRNAMT = 'ZLARHS'
551                     CALL ZLARHS( PATH, XTYPE, UPLO, ' ', N, N, KL, KU,
552     $                            NRHS, A, LDA, XACT, LDA, B, LDA,
553     $                            ISEED, INFO )
554                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
555*
556                     SRNAMT = 'ZHETRS2'
557                     CALL ZHETRS2( UPLO, N, NRHS, AFAC, LDA, IWORK, X,
558     $                            LDA, WORK, INFO )
559*
560*                    Check error code from ZHETRS2 and handle error.
561*
562                     IF( INFO.NE.0 )
563     $                  CALL ALAERH( PATH, 'ZHETRS2', INFO, 0, UPLO, N,
564     $                               N, -1, -1, NRHS, IMAT, NFAIL,
565     $                               NERRS, NOUT )
566*
567                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK, LDA )
568*
569*                    Compute the residual for the solution
570*
571                     CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA, WORK,
572     $                            LDA, RWORK, RESULT( 4 ) )
573*
574*+    TEST 5
575*                 Check solution from generated exact solution.
576*
577                     CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
578     $                            RESULT( 5 ) )
579*
580*+    TESTS 6, 7, and 8
581*                 Use iterative refinement to improve the solution.
582*
583                     SRNAMT = 'ZHERFS'
584                     CALL ZHERFS( UPLO, N, NRHS, A, LDA, AFAC, LDA,
585     $                            IWORK, B, LDA, X, LDA, RWORK,
586     $                            RWORK( NRHS+1 ), WORK,
587     $                            RWORK( 2*NRHS+1 ), INFO )
588*
589*                    Check error code from ZHERFS.
590*
591                     IF( INFO.NE.0 )
592     $                  CALL ALAERH( PATH, 'ZHERFS', INFO, 0, UPLO, N,
593     $                               N, -1, -1, NRHS, IMAT, NFAIL,
594     $                               NERRS, NOUT )
595*
596                     CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
597     $                            RESULT( 6 ) )
598                     CALL ZPOT05( UPLO, N, NRHS, A, LDA, B, LDA, X, LDA,
599     $                            XACT, LDA, RWORK, RWORK( NRHS+1 ),
600     $                            RESULT( 7 ) )
601*
602*                    Print information about the tests that did not pass
603*                    the threshold.
604*
605                     DO 120 K = 3, 8
606                        IF( RESULT( K ).GE.THRESH ) THEN
607                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
608     $                        CALL ALAHD( NOUT, PATH )
609                           WRITE( NOUT, FMT = 9998 )UPLO, N, NRHS,
610     $                        IMAT, K, RESULT( K )
611                           NFAIL = NFAIL + 1
612                        END IF
613  120                CONTINUE
614                     NRUN = NRUN + 6
615*
616*                 End do for each value of NRHS in NSVAL.
617*
618  130             CONTINUE
619*
620*+    TEST 9
621*                 Get an estimate of RCOND = 1/CNDNUM.
622*
623  140             CONTINUE
624                  ANORM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
625                  SRNAMT = 'ZHECON'
626                  CALL ZHECON( UPLO, N, AFAC, LDA, IWORK, ANORM, RCOND,
627     $                         WORK, INFO )
628*
629*                 Check error code from ZHECON and handle error.
630*
631                  IF( INFO.NE.0 )
632     $               CALL ALAERH( PATH, 'ZHECON', INFO, 0, UPLO, N, N,
633     $                            -1, -1, -1, IMAT, NFAIL, NERRS, NOUT )
634*
635                  RESULT( 9 ) = DGET06( RCOND, RCONDC )
636*
637*                 Print information about the tests that did not pass
638*                 the threshold.
639*
640                  IF( RESULT( 9 ).GE.THRESH ) THEN
641                     IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
642     $                  CALL ALAHD( NOUT, PATH )
643                     WRITE( NOUT, FMT = 9997 )UPLO, N, IMAT, 9,
644     $                  RESULT( 9 )
645                     NFAIL = NFAIL + 1
646                  END IF
647                  NRUN = NRUN + 1
648  150          CONTINUE
649  160       CONTINUE
650  170    CONTINUE
651  180 CONTINUE
652*
653*     Print a summary of the results.
654*
655      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
656*
657 9999 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NB =', I4, ', type ',
658     $      I2, ', test ', I2, ', ratio =', G12.5 )
659 9998 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NRHS=', I3, ', type ',
660     $      I2, ', test(', I2, ') =', G12.5 )
661 9997 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ',', 10X, ' type ', I2,
662     $      ', test(', I2, ') =', G12.5 )
663      RETURN
664*
665*     End of ZCHKHE
666*
667      END
668