1*> \brief \b ZCHKHE_RK
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_RK( 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( * ), E( * ),
25*      $                   WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> ZCHKHE_RK tests ZHETRF_RK, -TRI_3, -TRS_3,
35*> and -CON_3.
36*> \endverbatim
37*
38*  Arguments:
39*  ==========
40*
41*> \param[in] DOTYPE
42*> \verbatim
43*>          DOTYPE is LOGICAL array, dimension (NTYPES)
44*>          The matrix types to be used for testing.  Matrices of type j
45*>          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
46*>          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
47*> \endverbatim
48*>
49*> \param[in] NN
50*> \verbatim
51*>          NN is INTEGER
52*>          The number of values of N contained in the vector NVAL.
53*> \endverbatim
54*>
55*> \param[in] NVAL
56*> \verbatim
57*>          NVAL is INTEGER array, dimension (NN)
58*>          The values of the matrix dimension N.
59*> \endverbatim
60*>
61*> \param[in] NNB
62*> \verbatim
63*>          NNB is INTEGER
64*>          The number of values of NB contained in the vector NBVAL.
65*> \endverbatim
66*>
67*> \param[in] NBVAL
68*> \verbatim
69*>          NBVAL is INTEGER array, dimension (NNB)
70*>          The values of the blocksize NB.
71*> \endverbatim
72*>
73*> \param[in] NNS
74*> \verbatim
75*>          NNS is INTEGER
76*>          The number of values of NRHS contained in the vector NSVAL.
77*> \endverbatim
78*>
79*> \param[in] NSVAL
80*> \verbatim
81*>          NSVAL is INTEGER array, dimension (NNS)
82*>          The values of the number of right hand sides NRHS.
83*> \endverbatim
84*>
85*> \param[in] THRESH
86*> \verbatim
87*>          THRESH is DOUBLE PRECISION
88*>          The threshold value for the test ratios.  A result is
89*>          included in the output file if RESULT >= THRESH.  To have
90*>          every test ratio printed, use THRESH = 0.
91*> \endverbatim
92*>
93*> \param[in] TSTERR
94*> \verbatim
95*>          TSTERR is LOGICAL
96*>          Flag that indicates whether error exits are to be tested.
97*> \endverbatim
98*>
99*> \param[in] NMAX
100*> \verbatim
101*>          NMAX is INTEGER
102*>          The maximum value permitted for N, used in dimensioning the
103*>          work arrays.
104*> \endverbatim
105*>
106*> \param[out] A
107*> \verbatim
108*>          A is CCOMPLEX*16 array, dimension (NMAX*NMAX)
109*> \endverbatim
110*>
111*> \param[out] AFAC
112*> \verbatim
113*>          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
114*> \endverbatim
115*>
116*> \param[out] E
117*> \verbatim
118*>          E is COMPLEX*16 array, dimension (NMAX)
119*> \endverbatim
120*>
121*> \param[out] AINV
122*> \verbatim
123*>          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
124*> \endverbatim
125*>
126*> \param[out] B
127*> \verbatim
128*>          B is CCOMPLEX*16 array, dimension (NMAX*NSMAX)
129*>          where NSMAX is the largest entry in NSVAL.
130*> \endverbatim
131*>
132*> \param[out] X
133*> \verbatim
134*>          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
135*> \endverbatim
136*>
137*> \param[out] XACT
138*> \verbatim
139*>          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
140*> \endverbatim
141*>
142*> \param[out] WORK
143*> \verbatim
144*>          WORK is COMPLEX*16 array, dimension (NMAX*max(3,NSMAX))
145*> \endverbatim
146*>
147*> \param[out] RWORK
148*> \verbatim
149*>          RWORK is DOUBLE PRECISION array, dimension (max(NMAX,2*NSMAX)
150*> \endverbatim
151*>
152*> \param[out] IWORK
153*> \verbatim
154*>          IWORK is INTEGER array, dimension (2*NMAX)
155*> \endverbatim
156*>
157*> \param[in] NOUT
158*> \verbatim
159*>          NOUT is INTEGER
160*>          The unit number for output.
161*> \endverbatim
162*
163*  Authors:
164*  ========
165*
166*> \author Univ. of Tennessee
167*> \author Univ. of California Berkeley
168*> \author Univ. of Colorado Denver
169*> \author NAG Ltd.
170*
171*> \ingroup complex16_lin
172*
173*  =====================================================================
174      SUBROUTINE ZCHKHE_RK( DOTYPE, NN, NVAL, NNB, NBVAL, NNS, NSVAL,
175     $                      THRESH, TSTERR, NMAX, A, AFAC, E, AINV, B,
176     $                      X, XACT, WORK, RWORK, IWORK, NOUT )
177*
178*  -- LAPACK test routine --
179*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
180*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181*
182*     .. Scalar Arguments ..
183      LOGICAL            TSTERR
184      INTEGER            NMAX, NN, NNB, NNS, NOUT
185      DOUBLE PRECISION   THRESH
186*     ..
187*     .. Array Arguments ..
188      LOGICAL            DOTYPE( * )
189      INTEGER            IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
190      DOUBLE PRECISION   RWORK( * )
191      COMPLEX*16         A( * ), AFAC( * ), AINV( * ), B( * ), E( * ),
192     $                   WORK( * ), X( * ), XACT( * )
193*     ..
194*
195*  =====================================================================
196*
197*     .. Parameters ..
198      DOUBLE PRECISION   ZERO, ONE
199      PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
200      DOUBLE PRECISION   ONEHALF
201      PARAMETER          ( ONEHALF = 0.5D+0 )
202      DOUBLE PRECISION   EIGHT, SEVTEN
203      PARAMETER          ( EIGHT = 8.0D+0, SEVTEN = 17.0D+0 )
204      COMPLEX*16         CZERO
205      PARAMETER          ( CZERO = ( 0.0D+0, 0.0D+0 ) )
206      INTEGER            NTYPES
207      PARAMETER          ( NTYPES = 10 )
208      INTEGER            NTESTS
209      PARAMETER          ( NTESTS = 7 )
210*     ..
211*     .. Local Scalars ..
212      LOGICAL            TRFCON, ZEROT
213      CHARACTER          DIST, TYPE, UPLO, XTYPE
214      CHARACTER*3        PATH, MATPATH
215      INTEGER            I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
216     $                   ITEMP, ITEMP2, IUPLO, IZERO, J, K, KL, KU, LDA,
217     $                   LWORK, MODE, N, NB, NERRS, NFAIL, NIMAT, NRHS,
218     $                   NRUN, NT
219      DOUBLE PRECISION   ALPHA, ANORM, CNDNUM, CONST, SING_MAX,
220     $                   SING_MIN, RCOND, RCONDC, DTEMP
221*     ..
222*     .. Local Arrays ..
223      CHARACTER          UPLOS( 2 )
224      INTEGER            ISEED( 4 ), ISEEDY( 4 ), IDUMMY( 1 )
225      DOUBLE PRECISION   RESULT( NTESTS )
226      COMPLEX*16         BLOCK( 2, 2 ), ZDUMMY( 1 )
227*     ..
228*     .. External Functions ..
229      DOUBLE PRECISION   DGET06, ZLANGE, ZLANHE
230      EXTERNAL           DGET06, ZLANGE, ZLANHE
231*     ..
232*     .. External Subroutines ..
233      EXTERNAL           ALAERH, ALAHD, ALASUM, ZERRHE, ZGESVD, ZGET04,
234     $                   ZLACPY, ZLARHS, ZLATB4, ZLATMS, ZPOT02, ZPOT03,
235     $                   ZHECON_3, ZHET01_3, ZHETRF_RK, ZHETRI_3,
236     $                   ZHETRS_3, XLAENV
237*     ..
238*     .. Intrinsic Functions ..
239      INTRINSIC          DCONJG, MAX, MIN, SQRT
240*     ..
241*     .. Scalars in Common ..
242      LOGICAL            LERR, OK
243      CHARACTER*32       SRNAMT
244      INTEGER            INFOT, NUNIT
245*     ..
246*     .. Common blocks ..
247      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
248      COMMON             / SRNAMC / SRNAMT
249*     ..
250*     .. Data statements ..
251      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
252      DATA               UPLOS / 'U', 'L' /
253*     ..
254*     .. Executable Statements ..
255*
256*     Initialize constants and the random number seed.
257*
258      ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT
259*
260*     Test path
261*
262      PATH( 1: 1 ) = 'Zomplex precision'
263      PATH( 2: 3 ) = 'HK'
264*
265*     Path to generate matrices
266*
267      MATPATH( 1: 1 ) = 'Zomplex precision'
268      MATPATH( 2: 3 ) = 'HE'
269*
270      NRUN = 0
271      NFAIL = 0
272      NERRS = 0
273      DO 10 I = 1, 4
274         ISEED( I ) = ISEEDY( I )
275   10 CONTINUE
276*
277*     Test the error exits
278*
279      IF( TSTERR )
280     $   CALL ZERRHE( PATH, NOUT )
281      INFOT = 0
282*
283*     Set the minimum block size for which the block routine should
284*     be used, which will be later returned by ILAENV
285*
286      CALL XLAENV( 2, 2 )
287*
288*     Do for each value of N in NVAL
289*
290      DO 270 IN = 1, NN
291         N = NVAL( IN )
292         LDA = MAX( N, 1 )
293         XTYPE = 'N'
294         NIMAT = NTYPES
295         IF( N.LE.0 )
296     $      NIMAT = 1
297*
298         IZERO = 0
299*
300*        Do for each value of matrix type IMAT
301*
302         DO 260 IMAT = 1, NIMAT
303*
304*           Do the tests only if DOTYPE( IMAT ) is true.
305*
306            IF( .NOT.DOTYPE( IMAT ) )
307     $         GO TO 260
308*
309*           Skip types 3, 4, 5, or 6 if the matrix size is too small.
310*
311            ZEROT = IMAT.GE.3 .AND. IMAT.LE.6
312            IF( ZEROT .AND. N.LT.IMAT-2 )
313     $         GO TO 260
314*
315*           Do first for UPLO = 'U', then for UPLO = 'L'
316*
317            DO 250 IUPLO = 1, 2
318               UPLO = UPLOS( IUPLO )
319*
320*                 Begin generate the test matrix A.
321*
322*                 Set up parameters with ZLATB4 for the matrix generator
323*                 based on the type of matrix to be generated.
324*
325                  CALL ZLATB4( MATPATH, IMAT, N, N, TYPE, KL, KU, ANORM,
326     $                         MODE, CNDNUM, DIST )
327*
328*                 Generate a matrix with ZLATMS.
329*
330                  SRNAMT = 'ZLATMS'
331                  CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
332     $                         CNDNUM, ANORM, KL, KU, UPLO, A, LDA,
333     $                         WORK, INFO )
334*
335*                 Check error code from ZLATMS and handle error.
336*
337                  IF( INFO.NE.0 ) THEN
338                     CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N,
339     $                            -1, -1, -1, IMAT, NFAIL, NERRS, NOUT )
340*
341*                    Skip all tests for this generated matrix
342*
343                     GO TO 250
344                  END IF
345*
346*                 For matrix types 3-6, zero one or more rows and
347*                 columns of the matrix to test that INFO is returned
348*                 correctly.
349*
350                  IF( ZEROT ) THEN
351                     IF( IMAT.EQ.3 ) THEN
352                        IZERO = 1
353                     ELSE IF( IMAT.EQ.4 ) THEN
354                        IZERO = N
355                     ELSE
356                        IZERO = N / 2 + 1
357                     END IF
358*
359                     IF( IMAT.LT.6 ) THEN
360*
361*                       Set row and column IZERO to zero.
362*
363                        IF( IUPLO.EQ.1 ) THEN
364                           IOFF = ( IZERO-1 )*LDA
365                           DO 20 I = 1, IZERO - 1
366                              A( IOFF+I ) = CZERO
367   20                      CONTINUE
368                           IOFF = IOFF + IZERO
369                           DO 30 I = IZERO, N
370                              A( IOFF ) = CZERO
371                              IOFF = IOFF + LDA
372   30                      CONTINUE
373                        ELSE
374                           IOFF = IZERO
375                           DO 40 I = 1, IZERO - 1
376                              A( IOFF ) = CZERO
377                              IOFF = IOFF + LDA
378   40                      CONTINUE
379                           IOFF = IOFF - IZERO
380                           DO 50 I = IZERO, N
381                              A( IOFF+I ) = CZERO
382   50                      CONTINUE
383                        END IF
384                     ELSE
385                        IF( IUPLO.EQ.1 ) THEN
386*
387*                          Set the first IZERO rows and columns to zero.
388*
389                           IOFF = 0
390                           DO 70 J = 1, N
391                              I2 = MIN( J, IZERO )
392                              DO 60 I = 1, I2
393                                 A( IOFF+I ) = CZERO
394   60                         CONTINUE
395                              IOFF = IOFF + LDA
396   70                      CONTINUE
397                        ELSE
398*
399*                          Set the last IZERO rows and columns to zero.
400*
401                           IOFF = 0
402                           DO 90 J = 1, N
403                              I1 = MAX( J, IZERO )
404                              DO 80 I = I1, N
405                                 A( IOFF+I ) = CZERO
406   80                         CONTINUE
407                              IOFF = IOFF + LDA
408   90                      CONTINUE
409                        END IF
410                     END IF
411                  ELSE
412                     IZERO = 0
413                  END IF
414*
415*                 End generate the test matrix A.
416*
417*
418*              Do for each value of NB in NBVAL
419*
420               DO 240 INB = 1, NNB
421*
422*                 Set the optimal blocksize, which will be later
423*                 returned by ILAENV.
424*
425                  NB = NBVAL( INB )
426                  CALL XLAENV( 1, NB )
427*
428*                 Copy the test matrix A into matrix AFAC which
429*                 will be factorized in place. This is needed to
430*                 preserve the test matrix A for subsequent tests.
431*
432                  CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
433*
434*                 Compute the L*D*L**T or U*D*U**T factorization of the
435*                 matrix. IWORK stores details of the interchanges and
436*                 the block structure of D. AINV is a work array for
437*                 block factorization, LWORK is the length of AINV.
438*
439                  LWORK = MAX( 2, NB )*LDA
440                  SRNAMT = 'ZHETRF_RK'
441                  CALL ZHETRF_RK( UPLO, N, AFAC, LDA, E, IWORK, AINV,
442     $                            LWORK, INFO )
443*
444*                 Adjust the expected value of INFO to account for
445*                 pivoting.
446*
447                  K = IZERO
448                  IF( K.GT.0 ) THEN
449  100                CONTINUE
450                     IF( IWORK( K ).LT.0 ) THEN
451                        IF( IWORK( K ).NE.-K ) THEN
452                           K = -IWORK( K )
453                           GO TO 100
454                        END IF
455                     ELSE IF( IWORK( K ).NE.K ) THEN
456                        K = IWORK( K )
457                        GO TO 100
458                     END IF
459                  END IF
460*
461*                 Check error code from ZHETRF_RK and handle error.
462*
463                  IF( INFO.NE.K)
464     $               CALL ALAERH( PATH, 'ZHETRF_RK', INFO, K,
465     $                            UPLO, N, N, -1, -1, NB, IMAT,
466     $                            NFAIL, NERRS, NOUT )
467*
468*                 Set the condition estimate flag if the INFO is not 0.
469*
470                  IF( INFO.NE.0 ) THEN
471                     TRFCON = .TRUE.
472                  ELSE
473                     TRFCON = .FALSE.
474                  END IF
475*
476*+    TEST 1
477*                 Reconstruct matrix from factors and compute residual.
478*
479                  CALL ZHET01_3( UPLO, N, A, LDA, AFAC, LDA, E, IWORK,
480     $                           AINV, LDA, RWORK, RESULT( 1 ) )
481                  NT = 1
482*
483*+    TEST 2
484*                 Form the inverse and compute the residual,
485*                 if the factorization was competed without INFO > 0
486*                 (i.e. there is no zero rows and columns).
487*                 Do it only for the first block size.
488*
489                  IF( INB.EQ.1 .AND. .NOT.TRFCON ) THEN
490                     CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
491                     SRNAMT = 'ZHETRI_3'
492*
493*                    Another reason that we need to compute the inverse
494*                    is that ZPOT03 produces RCONDC which is used later
495*                    in TEST6 and TEST7.
496*
497                     LWORK = (N+NB+1)*(NB+3)
498                     CALL ZHETRI_3( UPLO, N, AINV, LDA, E, IWORK, WORK,
499     $                              LWORK, INFO )
500*
501*                    Check error code from ZHETRI_3 and handle error.
502*
503                     IF( INFO.NE.0 )
504     $                  CALL ALAERH( PATH, 'ZHETRI_3', INFO, -1,
505     $                               UPLO, N, N, -1, -1, -1, IMAT,
506     $                               NFAIL, NERRS, NOUT )
507*
508*                    Compute the residual for a Hermitian matrix times
509*                    its inverse.
510*
511                     CALL ZPOT03( UPLO, N, A, LDA, AINV, LDA, WORK, LDA,
512     $                            RWORK, RCONDC, RESULT( 2 ) )
513                     NT = 2
514                  END IF
515*
516*                 Print information about the tests that did not pass
517*                 the threshold.
518*
519                  DO 110 K = 1, NT
520                     IF( RESULT( K ).GE.THRESH ) THEN
521                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
522     $                     CALL ALAHD( NOUT, PATH )
523                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
524     $                     RESULT( K )
525                        NFAIL = NFAIL + 1
526                     END IF
527  110             CONTINUE
528                  NRUN = NRUN + NT
529*
530*+    TEST 3
531*                 Compute largest element in U or L
532*
533                  RESULT( 3 ) = ZERO
534                  DTEMP = ZERO
535*
536                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) ) /
537     $                    ( ONE-ALPHA )
538*
539                  IF( IUPLO.EQ.1 ) THEN
540*
541*                 Compute largest element in U
542*
543                     K = N
544  120                CONTINUE
545                     IF( K.LE.1 )
546     $                  GO TO 130
547*
548                     IF( IWORK( K ).GT.ZERO ) THEN
549*
550*                       Get max absolute value from elements
551*                       in column k in U
552*
553                        DTEMP = ZLANGE( 'M', K-1, 1,
554     $                          AFAC( ( K-1 )*LDA+1 ), LDA, RWORK )
555                     ELSE
556*
557*                       Get max absolute value from elements
558*                       in columns k and k-1 in U
559*
560                        DTEMP = ZLANGE( 'M', K-2, 2,
561     $                          AFAC( ( K-2 )*LDA+1 ), LDA, RWORK )
562                        K = K - 1
563*
564                     END IF
565*
566*                    DTEMP should be bounded by CONST
567*
568                     DTEMP = DTEMP - CONST + THRESH
569                     IF( DTEMP.GT.RESULT( 3 ) )
570     $                  RESULT( 3 ) = DTEMP
571*
572                     K = K - 1
573*
574                     GO TO 120
575  130                CONTINUE
576*
577                  ELSE
578*
579*                 Compute largest element in L
580*
581                     K = 1
582  140                CONTINUE
583                     IF( K.GE.N )
584     $                  GO TO 150
585*
586                     IF( IWORK( K ).GT.ZERO ) THEN
587*
588*                       Get max absolute value from elements
589*                       in column k in L
590*
591                        DTEMP = ZLANGE( 'M', N-K, 1,
592     $                          AFAC( ( K-1 )*LDA+K+1 ), LDA, RWORK )
593                     ELSE
594*
595*                       Get max absolute value from elements
596*                       in columns k and k+1 in L
597*
598                        DTEMP = ZLANGE( 'M', N-K-1, 2,
599     $                          AFAC( ( K-1 )*LDA+K+2 ), LDA, RWORK )
600                        K = K + 1
601*
602                     END IF
603*
604*                    DTEMP should be bounded by CONST
605*
606                     DTEMP = DTEMP - CONST + THRESH
607                     IF( DTEMP.GT.RESULT( 3 ) )
608     $                  RESULT( 3 ) = DTEMP
609*
610                     K = K + 1
611*
612                     GO TO 140
613  150                CONTINUE
614                  END IF
615*
616*
617*+    TEST 4
618*                 Compute largest 2-Norm (condition number)
619*                 of 2-by-2 diag blocks
620*
621                  RESULT( 4 ) = ZERO
622                  DTEMP = ZERO
623*
624                  CONST = ( ( ALPHA**2-ONE ) / ( ALPHA**2-ONEHALF ) )*
625     $                    ( ( ONE + ALPHA ) / ( ONE - ALPHA ) )
626                  CALL ZLACPY( UPLO, N, N, AFAC, LDA, AINV, LDA )
627*
628                  IF( IUPLO.EQ.1 ) THEN
629*
630*                    Loop backward for UPLO = 'U'
631*
632                     K = N
633  160                CONTINUE
634                     IF( K.LE.1 )
635     $                  GO TO 170
636*
637                     IF( IWORK( K ).LT.ZERO ) THEN
638*
639*                       Get the two singular values
640*                       (real and non-negative) of a 2-by-2 block,
641*                       store them in RWORK array
642*
643                        BLOCK( 1, 1 ) = AFAC( ( K-2 )*LDA+K-1 )
644                        BLOCK( 1, 2 ) = E( K )
645                        BLOCK( 2, 1 ) = DCONJG( BLOCK( 1, 2 ) )
646                        BLOCK( 2, 2 ) = AFAC( (K-1)*LDA+K )
647*
648                        CALL ZGESVD( 'N', 'N', 2, 2, BLOCK, 2, RWORK,
649     $                               ZDUMMY, 1, ZDUMMY, 1,
650     $                               WORK, 6, RWORK( 3 ), INFO )
651*
652*
653                        SING_MAX = RWORK( 1 )
654                        SING_MIN = RWORK( 2 )
655*
656                        DTEMP = SING_MAX / SING_MIN
657*
658*                       DTEMP should be bounded by CONST
659*
660                        DTEMP = DTEMP - CONST + THRESH
661                        IF( DTEMP.GT.RESULT( 4 ) )
662     $                     RESULT( 4 ) = DTEMP
663                        K = K - 1
664*
665                     END IF
666*
667                     K = K - 1
668*
669                     GO TO 160
670  170                CONTINUE
671*
672                  ELSE
673*
674*                    Loop forward for UPLO = 'L'
675*
676                     K = 1
677  180                CONTINUE
678                     IF( K.GE.N )
679     $                  GO TO 190
680*
681                     IF( IWORK( K ).LT.ZERO ) THEN
682*
683*                       Get the two singular values
684*                       (real and non-negative) of a 2-by-2 block,
685*                       store them in RWORK array
686*
687                        BLOCK( 1, 1 ) = AFAC( ( K-1 )*LDA+K )
688                        BLOCK( 2, 1 ) = E( K )
689                        BLOCK( 1, 2 ) = DCONJG( BLOCK( 2, 1 ) )
690                        BLOCK( 2, 2 ) = AFAC( K*LDA+K+1 )
691*
692                        CALL ZGESVD( 'N', 'N', 2, 2, BLOCK, 2, RWORK,
693     $                               ZDUMMY, 1, ZDUMMY, 1,
694     $                               WORK, 6, RWORK(3), INFO )
695*
696                        SING_MAX = RWORK( 1 )
697                        SING_MIN = RWORK( 2 )
698*
699                        DTEMP = SING_MAX / SING_MIN
700*
701*                       DTEMP should be bounded by CONST
702*
703                        DTEMP = DTEMP - CONST + THRESH
704                        IF( DTEMP.GT.RESULT( 4 ) )
705     $                     RESULT( 4 ) = DTEMP
706                        K = K + 1
707*
708                     END IF
709*
710                     K = K + 1
711*
712                     GO TO 180
713  190                CONTINUE
714                  END IF
715*
716*                 Print information about the tests that did not pass
717*                 the threshold.
718*
719                  DO 200 K = 3, 4
720                     IF( RESULT( K ).GE.THRESH ) THEN
721                        IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
722     $                     CALL ALAHD( NOUT, PATH )
723                        WRITE( NOUT, FMT = 9999 )UPLO, N, NB, IMAT, K,
724     $                     RESULT( K )
725                        NFAIL = NFAIL + 1
726                     END IF
727  200             CONTINUE
728                  NRUN = NRUN + 2
729*
730*                 Skip the other tests if this is not the first block
731*                 size.
732*
733                  IF( INB.GT.1 )
734     $               GO TO 240
735*
736*                 Do only the condition estimate if INFO is not 0.
737*
738                  IF( TRFCON ) THEN
739                     RCONDC = ZERO
740                     GO TO 230
741                  END IF
742*
743*                 Do for each value of NRHS in NSVAL.
744*
745                  DO 220 IRHS = 1, NNS
746                     NRHS = NSVAL( IRHS )
747*
748*                    Begin loop over NRHS values
749*
750*
751*+    TEST 5 ( Using TRS_3)
752*                 Solve and compute residual for  A * X = B.
753*
754*                    Choose a set of NRHS random solution vectors
755*                    stored in XACT and set up the right hand side B
756*
757                     SRNAMT = 'ZLARHS'
758                     CALL ZLARHS( MATPATH, XTYPE, UPLO, ' ', N, N,
759     $                            KL, KU, NRHS, A, LDA, XACT, LDA,
760     $                            B, LDA, ISEED, INFO )
761                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
762*
763                     SRNAMT = 'ZHETRS_3'
764                     CALL ZHETRS_3( UPLO, N, NRHS, AFAC, LDA, E, IWORK,
765     $                              X, LDA, INFO )
766*
767*                    Check error code from ZHETRS_3 and handle error.
768*
769                     IF( INFO.NE.0 )
770     $                  CALL ALAERH( PATH, 'ZHETRS_3', INFO, 0,
771     $                               UPLO, N, N, -1, -1, NRHS, IMAT,
772     $                               NFAIL, NERRS, NOUT )
773*
774                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK, LDA )
775*
776*                    Compute the residual for the solution
777*
778                     CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA, WORK,
779     $                            LDA, RWORK, RESULT( 5 ) )
780*
781*+    TEST 6
782*                 Check solution from generated exact solution.
783*
784                     CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
785     $                            RESULT( 6 ) )
786*
787*                    Print information about the tests that did not pass
788*                    the threshold.
789*
790                     DO 210 K = 5, 6
791                        IF( RESULT( K ).GE.THRESH ) THEN
792                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
793     $                        CALL ALAHD( NOUT, PATH )
794                           WRITE( NOUT, FMT = 9998 )UPLO, N, NRHS,
795     $                        IMAT, K, RESULT( K )
796                           NFAIL = NFAIL + 1
797                        END IF
798  210                CONTINUE
799                     NRUN = NRUN + 2
800*
801*                 End do for each value of NRHS in NSVAL.
802*
803  220             CONTINUE
804*
805*+    TEST 7
806*                 Get an estimate of RCOND = 1/CNDNUM.
807*
808  230             CONTINUE
809                  ANORM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
810                  SRNAMT = 'ZHECON_3'
811                  CALL ZHECON_3( UPLO, N, AFAC, LDA, E, IWORK, ANORM,
812     $                           RCOND, WORK, INFO )
813*
814*                 Check error code from ZHECON_3 and handle error.
815*
816                  IF( INFO.NE.0 )
817     $               CALL ALAERH( PATH, 'ZHECON_3', INFO, 0,
818     $                            UPLO, N, N, -1, -1, -1, IMAT,
819     $                            NFAIL, NERRS, NOUT )
820*
821*                 Compute the test ratio to compare values of RCOND
822*
823                  RESULT( 7 ) = DGET06( RCOND, RCONDC )
824*
825*                 Print information about the tests that did not pass
826*                 the threshold.
827*
828                  IF( RESULT( 7 ).GE.THRESH ) THEN
829                     IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
830     $                  CALL ALAHD( NOUT, PATH )
831                     WRITE( NOUT, FMT = 9997 )UPLO, N, IMAT, 7,
832     $                  RESULT( 7 )
833                     NFAIL = NFAIL + 1
834                  END IF
835                  NRUN = NRUN + 1
836  240          CONTINUE
837*
838  250       CONTINUE
839  260    CONTINUE
840  270 CONTINUE
841*
842*     Print a summary of the results.
843*
844      CALL ALASUM( PATH, NOUT, NFAIL, NRUN, NERRS )
845*
846 9999 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NB =', I4, ', type ',
847     $      I2, ', test ', I2, ', ratio =', G12.5 )
848 9998 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ', NRHS=', I3, ', type ',
849     $      I2, ', test ', I2, ', ratio =', G12.5 )
850 9997 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ',', 10X, ' type ', I2,
851     $      ', test ', I2, ', ratio =', G12.5 )
852      RETURN
853*
854*     End of ZCHKHE_RK
855*
856      END
857