1*> \brief \b ZDRVPOX
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 ZDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12*                          A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
13*                          RWORK, NOUT )
14*
15*       .. Scalar Arguments ..
16*       LOGICAL            TSTERR
17*       INTEGER            NMAX, NN, NOUT, NRHS
18*       DOUBLE PRECISION   THRESH
19*       ..
20*       .. Array Arguments ..
21*       LOGICAL            DOTYPE( * )
22*       INTEGER            NVAL( * )
23*       DOUBLE PRECISION   RWORK( * ), S( * )
24*       COMPLEX*16         A( * ), AFAC( * ), ASAV( * ), B( * ),
25*      $                   BSAV( * ), WORK( * ), X( * ), XACT( * )
26*       ..
27*
28*
29*> \par Purpose:
30*  =============
31*>
32*> \verbatim
33*>
34*> ZDRVPO tests the driver routines ZPOSV, -SVX, and -SVXX.
35*>
36*> Note that this file is used only when the XBLAS are available,
37*> otherwise zdrvpo.f defines this subroutine.
38*> \endverbatim
39*
40*  Arguments:
41*  ==========
42*
43*> \param[in] DOTYPE
44*> \verbatim
45*>          DOTYPE is LOGICAL array, dimension (NTYPES)
46*>          The matrix types to be used for testing.  Matrices of type j
47*>          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
48*>          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
49*> \endverbatim
50*>
51*> \param[in] NN
52*> \verbatim
53*>          NN is INTEGER
54*>          The number of values of N contained in the vector NVAL.
55*> \endverbatim
56*>
57*> \param[in] NVAL
58*> \verbatim
59*>          NVAL is INTEGER array, dimension (NN)
60*>          The values of the matrix dimension N.
61*> \endverbatim
62*>
63*> \param[in] NRHS
64*> \verbatim
65*>          NRHS is INTEGER
66*>          The number of right hand side vectors to be generated for
67*>          each linear system.
68*> \endverbatim
69*>
70*> \param[in] THRESH
71*> \verbatim
72*>          THRESH is DOUBLE PRECISION
73*>          The threshold value for the test ratios.  A result is
74*>          included in the output file if RESULT >= THRESH.  To have
75*>          every test ratio printed, use THRESH = 0.
76*> \endverbatim
77*>
78*> \param[in] TSTERR
79*> \verbatim
80*>          TSTERR is LOGICAL
81*>          Flag that indicates whether error exits are to be tested.
82*> \endverbatim
83*>
84*> \param[in] NMAX
85*> \verbatim
86*>          NMAX is INTEGER
87*>          The maximum value permitted for N, used in dimensioning the
88*>          work arrays.
89*> \endverbatim
90*>
91*> \param[out] A
92*> \verbatim
93*>          A is COMPLEX*16 array, dimension (NMAX*NMAX)
94*> \endverbatim
95*>
96*> \param[out] AFAC
97*> \verbatim
98*>          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
99*> \endverbatim
100*>
101*> \param[out] ASAV
102*> \verbatim
103*>          ASAV is COMPLEX*16 array, dimension (NMAX*NMAX)
104*> \endverbatim
105*>
106*> \param[out] B
107*> \verbatim
108*>          B is COMPLEX*16 array, dimension (NMAX*NRHS)
109*> \endverbatim
110*>
111*> \param[out] BSAV
112*> \verbatim
113*>          BSAV is COMPLEX*16 array, dimension (NMAX*NRHS)
114*> \endverbatim
115*>
116*> \param[out] X
117*> \verbatim
118*>          X is COMPLEX*16 array, dimension (NMAX*NRHS)
119*> \endverbatim
120*>
121*> \param[out] XACT
122*> \verbatim
123*>          XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
124*> \endverbatim
125*>
126*> \param[out] S
127*> \verbatim
128*>          S is DOUBLE PRECISION array, dimension (NMAX)
129*> \endverbatim
130*>
131*> \param[out] WORK
132*> \verbatim
133*>          WORK is COMPLEX*16 array, dimension
134*>                      (NMAX*max(3,NRHS))
135*> \endverbatim
136*>
137*> \param[out] RWORK
138*> \verbatim
139*>          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
140*> \endverbatim
141*>
142*> \param[in] NOUT
143*> \verbatim
144*>          NOUT is INTEGER
145*>          The unit number for output.
146*> \endverbatim
147*
148*  Authors:
149*  ========
150*
151*> \author Univ. of Tennessee
152*> \author Univ. of California Berkeley
153*> \author Univ. of Colorado Denver
154*> \author NAG Ltd.
155*
156*> \ingroup complex16_lin
157*
158*  =====================================================================
159      SUBROUTINE ZDRVPO( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
160     $                   A, AFAC, ASAV, B, BSAV, X, XACT, S, WORK,
161     $                   RWORK, NOUT )
162*
163*  -- LAPACK test routine --
164*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
165*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
166*
167*     .. Scalar Arguments ..
168      LOGICAL            TSTERR
169      INTEGER            NMAX, NN, NOUT, NRHS
170      DOUBLE PRECISION   THRESH
171*     ..
172*     .. Array Arguments ..
173      LOGICAL            DOTYPE( * )
174      INTEGER            NVAL( * )
175      DOUBLE PRECISION   RWORK( * ), S( * )
176      COMPLEX*16         A( * ), AFAC( * ), ASAV( * ), B( * ),
177     $                   BSAV( * ), WORK( * ), X( * ), XACT( * )
178*     ..
179*
180*  =====================================================================
181*
182*     .. Parameters ..
183      DOUBLE PRECISION   ONE, ZERO
184      PARAMETER          ( ONE = 1.0D+0, ZERO = 0.0D+0 )
185      INTEGER            NTYPES
186      PARAMETER          ( NTYPES = 9 )
187      INTEGER            NTESTS
188      PARAMETER          ( NTESTS = 6 )
189*     ..
190*     .. Local Scalars ..
191      LOGICAL            EQUIL, NOFACT, PREFAC, ZEROT
192      CHARACTER          DIST, EQUED, FACT, TYPE, UPLO, XTYPE
193      CHARACTER*3        PATH
194      INTEGER            I, IEQUED, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
195     $                   IZERO, K, K1, KL, KU, LDA, MODE, N, NB, NBMIN,
196     $                   NERRS, NFACT, NFAIL, NIMAT, NRUN, NT,
197     $                   N_ERR_BNDS
198      DOUBLE PRECISION   AINVNM, AMAX, ANORM, CNDNUM, RCOND, RCONDC,
199     $                   ROLDC, SCOND, RPVGRW_SVXX
200*     ..
201*     .. Local Arrays ..
202      CHARACTER          EQUEDS( 2 ), FACTS( 3 ), UPLOS( 2 )
203      INTEGER            ISEED( 4 ), ISEEDY( 4 )
204      DOUBLE PRECISION   RESULT( NTESTS ), BERR( NRHS ),
205     $                   ERRBNDS_N( NRHS, 3 ), ERRBNDS_C( NRHS, 3 )
206*     ..
207*     .. External Functions ..
208      LOGICAL            LSAME
209      DOUBLE PRECISION   DGET06, ZLANHE
210      EXTERNAL           LSAME, DGET06, ZLANHE
211*     ..
212*     .. External Subroutines ..
213      EXTERNAL           ALADHD, ALAERH, ALASVM, XLAENV, ZERRVX, ZGET04,
214     $                   ZLACPY, ZLAIPD, ZLAQHE, ZLARHS, ZLASET, ZLATB4,
215     $                   ZLATMS, ZPOEQU, ZPOSV, ZPOSVX, ZPOT01, ZPOT02,
216     $                   ZPOT05, ZPOTRF, ZPOTRI, ZPOSVXX
217*     ..
218*     .. Scalars in Common ..
219      LOGICAL            LERR, OK
220      CHARACTER*32       SRNAMT
221      INTEGER            INFOT, NUNIT
222*     ..
223*     .. Common blocks ..
224      COMMON             / INFOC / INFOT, NUNIT, OK, LERR
225      COMMON             / SRNAMC / SRNAMT
226*     ..
227*     .. Intrinsic Functions ..
228      INTRINSIC          DCMPLX, MAX
229*     ..
230*     .. Data statements ..
231      DATA               ISEEDY / 1988, 1989, 1990, 1991 /
232      DATA               UPLOS / 'U', 'L' /
233      DATA               FACTS / 'F', 'N', 'E' /
234      DATA               EQUEDS / 'N', 'Y' /
235*     ..
236*     .. Executable Statements ..
237*
238*     Initialize constants and the random number seed.
239*
240      PATH( 1: 1 ) = 'Zomplex precision'
241      PATH( 2: 3 ) = 'PO'
242      NRUN = 0
243      NFAIL = 0
244      NERRS = 0
245      DO 10 I = 1, 4
246         ISEED( I ) = ISEEDY( I )
247   10 CONTINUE
248*
249*     Test the error exits
250*
251      IF( TSTERR )
252     $   CALL ZERRVX( PATH, NOUT )
253      INFOT = 0
254*
255*     Set the block size and minimum block size for testing.
256*
257      NB = 1
258      NBMIN = 2
259      CALL XLAENV( 1, NB )
260      CALL XLAENV( 2, NBMIN )
261*
262*     Do for each value of N in NVAL
263*
264      DO 130 IN = 1, NN
265         N = NVAL( IN )
266         LDA = MAX( N, 1 )
267         XTYPE = 'N'
268         NIMAT = NTYPES
269         IF( N.LE.0 )
270     $      NIMAT = 1
271*
272         DO 120 IMAT = 1, NIMAT
273*
274*           Do the tests only if DOTYPE( IMAT ) is true.
275*
276            IF( .NOT.DOTYPE( IMAT ) )
277     $         GO TO 120
278*
279*           Skip types 3, 4, or 5 if the matrix size is too small.
280*
281            ZEROT = IMAT.GE.3 .AND. IMAT.LE.5
282            IF( ZEROT .AND. N.LT.IMAT-2 )
283     $         GO TO 120
284*
285*           Do first for UPLO = 'U', then for UPLO = 'L'
286*
287            DO 110 IUPLO = 1, 2
288               UPLO = UPLOS( IUPLO )
289*
290*              Set up parameters with ZLATB4 and generate a test matrix
291*              with ZLATMS.
292*
293               CALL ZLATB4( PATH, IMAT, N, N, TYPE, KL, KU, ANORM, MODE,
294     $                      CNDNUM, DIST )
295*
296               SRNAMT = 'ZLATMS'
297               CALL ZLATMS( N, N, DIST, ISEED, TYPE, RWORK, MODE,
298     $                      CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
299     $                      INFO )
300*
301*              Check error code from ZLATMS.
302*
303               IF( INFO.NE.0 ) THEN
304                  CALL ALAERH( PATH, 'ZLATMS', INFO, 0, UPLO, N, N, -1,
305     $                         -1, -1, IMAT, NFAIL, NERRS, NOUT )
306                  GO TO 110
307               END IF
308*
309*              For types 3-5, zero one row and column of the matrix to
310*              test that INFO is returned correctly.
311*
312               IF( ZEROT ) THEN
313                  IF( IMAT.EQ.3 ) THEN
314                     IZERO = 1
315                  ELSE IF( IMAT.EQ.4 ) THEN
316                     IZERO = N
317                  ELSE
318                     IZERO = N / 2 + 1
319                  END IF
320                  IOFF = ( IZERO-1 )*LDA
321*
322*                 Set row and column IZERO of A to 0.
323*
324                  IF( IUPLO.EQ.1 ) THEN
325                     DO 20 I = 1, IZERO - 1
326                        A( IOFF+I ) = ZERO
327   20                CONTINUE
328                     IOFF = IOFF + IZERO
329                     DO 30 I = IZERO, N
330                        A( IOFF ) = ZERO
331                        IOFF = IOFF + LDA
332   30                CONTINUE
333                  ELSE
334                     IOFF = IZERO
335                     DO 40 I = 1, IZERO - 1
336                        A( IOFF ) = ZERO
337                        IOFF = IOFF + LDA
338   40                CONTINUE
339                     IOFF = IOFF - IZERO
340                     DO 50 I = IZERO, N
341                        A( IOFF+I ) = ZERO
342   50                CONTINUE
343                  END IF
344               ELSE
345                  IZERO = 0
346               END IF
347*
348*              Set the imaginary part of the diagonals.
349*
350               CALL ZLAIPD( N, A, LDA+1, 0 )
351*
352*              Save a copy of the matrix A in ASAV.
353*
354               CALL ZLACPY( UPLO, N, N, A, LDA, ASAV, LDA )
355*
356               DO 100 IEQUED = 1, 2
357                  EQUED = EQUEDS( IEQUED )
358                  IF( IEQUED.EQ.1 ) THEN
359                     NFACT = 3
360                  ELSE
361                     NFACT = 1
362                  END IF
363*
364                  DO 90 IFACT = 1, NFACT
365                     FACT = FACTS( IFACT )
366                     PREFAC = LSAME( FACT, 'F' )
367                     NOFACT = LSAME( FACT, 'N' )
368                     EQUIL = LSAME( FACT, 'E' )
369*
370                     IF( ZEROT ) THEN
371                        IF( PREFAC )
372     $                     GO TO 90
373                        RCONDC = ZERO
374*
375                     ELSE IF( .NOT.LSAME( FACT, 'N' ) ) THEN
376*
377*                       Compute the condition number for comparison with
378*                       the value returned by ZPOSVX (FACT = 'N' reuses
379*                       the condition number from the previous iteration
380*                       with FACT = 'F').
381*
382                        CALL ZLACPY( UPLO, N, N, ASAV, LDA, AFAC, LDA )
383                        IF( EQUIL .OR. IEQUED.GT.1 ) THEN
384*
385*                          Compute row and column scale factors to
386*                          equilibrate the matrix A.
387*
388                           CALL ZPOEQU( N, AFAC, LDA, S, SCOND, AMAX,
389     $                                  INFO )
390                           IF( INFO.EQ.0 .AND. N.GT.0 ) THEN
391                              IF( IEQUED.GT.1 )
392     $                           SCOND = ZERO
393*
394*                             Equilibrate the matrix.
395*
396                              CALL ZLAQHE( UPLO, N, AFAC, LDA, S, SCOND,
397     $                                     AMAX, EQUED )
398                           END IF
399                        END IF
400*
401*                       Save the condition number of the
402*                       non-equilibrated system for use in ZGET04.
403*
404                        IF( EQUIL )
405     $                     ROLDC = RCONDC
406*
407*                       Compute the 1-norm of A.
408*
409                        ANORM = ZLANHE( '1', UPLO, N, AFAC, LDA, RWORK )
410*
411*                       Factor the matrix A.
412*
413                        CALL ZPOTRF( UPLO, N, AFAC, LDA, INFO )
414*
415*                       Form the inverse of A.
416*
417                        CALL ZLACPY( UPLO, N, N, AFAC, LDA, A, LDA )
418                        CALL ZPOTRI( UPLO, N, A, LDA, INFO )
419*
420*                       Compute the 1-norm condition number of A.
421*
422                        AINVNM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
423                        IF( ANORM.LE.ZERO .OR. AINVNM.LE.ZERO ) THEN
424                           RCONDC = ONE
425                        ELSE
426                           RCONDC = ( ONE / ANORM ) / AINVNM
427                        END IF
428                     END IF
429*
430*                    Restore the matrix A.
431*
432                     CALL ZLACPY( UPLO, N, N, ASAV, LDA, A, LDA )
433*
434*                    Form an exact solution and set the right hand side.
435*
436                     SRNAMT = 'ZLARHS'
437                     CALL ZLARHS( PATH, XTYPE, UPLO, ' ', N, N, KL, KU,
438     $                            NRHS, A, LDA, XACT, LDA, B, LDA,
439     $                            ISEED, INFO )
440                     XTYPE = 'C'
441                     CALL ZLACPY( 'Full', N, NRHS, B, LDA, BSAV, LDA )
442*
443                     IF( NOFACT ) THEN
444*
445*                       --- Test ZPOSV  ---
446*
447*                       Compute the L*L' or U'*U factorization of the
448*                       matrix and solve the system.
449*
450                        CALL ZLACPY( UPLO, N, N, A, LDA, AFAC, LDA )
451                        CALL ZLACPY( 'Full', N, NRHS, B, LDA, X, LDA )
452*
453                        SRNAMT = 'ZPOSV '
454                        CALL ZPOSV( UPLO, N, NRHS, AFAC, LDA, X, LDA,
455     $                              INFO )
456*
457*                       Check error code from ZPOSV .
458*
459                        IF( INFO.NE.IZERO ) THEN
460                           CALL ALAERH( PATH, 'ZPOSV ', INFO, IZERO,
461     $                                  UPLO, N, N, -1, -1, NRHS, IMAT,
462     $                                  NFAIL, NERRS, NOUT )
463                           GO TO 70
464                        ELSE IF( INFO.NE.0 ) THEN
465                           GO TO 70
466                        END IF
467*
468*                       Reconstruct matrix from factors and compute
469*                       residual.
470*
471                        CALL ZPOT01( UPLO, N, A, LDA, AFAC, LDA, RWORK,
472     $                               RESULT( 1 ) )
473*
474*                       Compute residual of the computed solution.
475*
476                        CALL ZLACPY( 'Full', N, NRHS, B, LDA, WORK,
477     $                               LDA )
478                        CALL ZPOT02( UPLO, N, NRHS, A, LDA, X, LDA,
479     $                               WORK, LDA, RWORK, RESULT( 2 ) )
480*
481*                       Check solution from generated exact solution.
482*
483                        CALL ZGET04( N, NRHS, X, LDA, XACT, LDA, RCONDC,
484     $                               RESULT( 3 ) )
485                        NT = 3
486*
487*                       Print information about the tests that did not
488*                       pass the threshold.
489*
490                        DO 60 K = 1, NT
491                           IF( RESULT( K ).GE.THRESH ) THEN
492                              IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
493     $                           CALL ALADHD( NOUT, PATH )
494                              WRITE( NOUT, FMT = 9999 )'ZPOSV ', UPLO,
495     $                           N, IMAT, K, RESULT( K )
496                              NFAIL = NFAIL + 1
497                           END IF
498   60                   CONTINUE
499                        NRUN = NRUN + NT
500   70                   CONTINUE
501                     END IF
502*
503*                    --- Test ZPOSVX ---
504*
505                     IF( .NOT.PREFAC )
506     $                  CALL ZLASET( UPLO, N, N, DCMPLX( ZERO ),
507     $                               DCMPLX( ZERO ), AFAC, LDA )
508                     CALL ZLASET( 'Full', N, NRHS, DCMPLX( ZERO ),
509     $                            DCMPLX( ZERO ), X, LDA )
510                     IF( IEQUED.GT.1 .AND. N.GT.0 ) THEN
511*
512*                       Equilibrate the matrix if FACT='F' and
513*                       EQUED='Y'.
514*
515                        CALL ZLAQHE( UPLO, N, A, LDA, S, SCOND, AMAX,
516     $                               EQUED )
517                     END IF
518*
519*                    Solve the system and compute the condition number
520*                    and error bounds using ZPOSVX.
521*
522                     SRNAMT = 'ZPOSVX'
523                     CALL ZPOSVX( FACT, UPLO, N, NRHS, A, LDA, AFAC,
524     $                            LDA, EQUED, S, B, LDA, X, LDA, RCOND,
525     $                            RWORK, RWORK( NRHS+1 ), WORK,
526     $                            RWORK( 2*NRHS+1 ), INFO )
527*
528*                    Check the error code from ZPOSVX.
529*
530                     IF( INFO.NE.IZERO ) THEN
531                        CALL ALAERH( PATH, 'ZPOSVX', INFO, IZERO,
532     $                               FACT // UPLO, N, N, -1, -1, NRHS,
533     $                               IMAT, NFAIL, NERRS, NOUT )
534                        GO TO 90
535                     END IF
536*
537                     IF( INFO.EQ.0 ) THEN
538                        IF( .NOT.PREFAC ) THEN
539*
540*                          Reconstruct matrix from factors and compute
541*                          residual.
542*
543                           CALL ZPOT01( UPLO, N, A, LDA, AFAC, LDA,
544     $                                  RWORK( 2*NRHS+1 ), RESULT( 1 ) )
545                           K1 = 1
546                        ELSE
547                           K1 = 2
548                        END IF
549*
550*                       Compute residual of the computed solution.
551*
552                        CALL ZLACPY( 'Full', N, NRHS, BSAV, LDA, WORK,
553     $                               LDA )
554                        CALL ZPOT02( UPLO, N, NRHS, ASAV, LDA, X, LDA,
555     $                               WORK, LDA, RWORK( 2*NRHS+1 ),
556     $                               RESULT( 2 ) )
557*
558*                       Check solution from generated exact solution.
559*
560                        IF( NOFACT .OR. ( PREFAC .AND. LSAME( EQUED,
561     $                      'N' ) ) ) THEN
562                           CALL ZGET04( N, NRHS, X, LDA, XACT, LDA,
563     $                                  RCONDC, RESULT( 3 ) )
564                        ELSE
565                           CALL ZGET04( N, NRHS, X, LDA, XACT, LDA,
566     $                                  ROLDC, RESULT( 3 ) )
567                        END IF
568*
569*                       Check the error bounds from iterative
570*                       refinement.
571*
572                        CALL ZPOT05( UPLO, N, NRHS, ASAV, LDA, B, LDA,
573     $                               X, LDA, XACT, LDA, RWORK,
574     $                               RWORK( NRHS+1 ), RESULT( 4 ) )
575                     ELSE
576                        K1 = 6
577                     END IF
578*
579*                    Compare RCOND from ZPOSVX with the computed value
580*                    in RCONDC.
581*
582                     RESULT( 6 ) = DGET06( RCOND, RCONDC )
583*
584*                    Print information about the tests that did not pass
585*                    the threshold.
586*
587                     DO 80 K = K1, 6
588                        IF( RESULT( K ).GE.THRESH ) THEN
589                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
590     $                        CALL ALADHD( NOUT, PATH )
591                           IF( PREFAC ) THEN
592                              WRITE( NOUT, FMT = 9997 )'ZPOSVX', FACT,
593     $                           UPLO, N, EQUED, IMAT, K, RESULT( K )
594                           ELSE
595                              WRITE( NOUT, FMT = 9998 )'ZPOSVX', FACT,
596     $                           UPLO, N, IMAT, K, RESULT( K )
597                           END IF
598                           NFAIL = NFAIL + 1
599                        END IF
600   80                CONTINUE
601                     NRUN = NRUN + 7 - K1
602*
603*                    --- Test ZPOSVXX ---
604*
605*                    Restore the matrices A and B.
606*
607                     CALL ZLACPY( 'Full', N, N, ASAV, LDA, A, LDA )
608                     CALL ZLACPY( 'Full', N, NRHS, BSAV, LDA, B, LDA )
609
610                     IF( .NOT.PREFAC )
611     $                  CALL ZLASET( UPLO, N, N, DCMPLX( ZERO ),
612     $                               DCMPLX( ZERO ), AFAC, LDA )
613                     CALL ZLASET( 'Full', N, NRHS, DCMPLX( ZERO ),
614     $                            DCMPLX( ZERO ), X, LDA )
615                     IF( IEQUED.GT.1 .AND. N.GT.0 ) THEN
616*
617*                       Equilibrate the matrix if FACT='F' and
618*                       EQUED='Y'.
619*
620                        CALL ZLAQHE( UPLO, N, A, LDA, S, SCOND, AMAX,
621     $                               EQUED )
622                     END IF
623*
624*                    Solve the system and compute the condition number
625*                    and error bounds using ZPOSVXX.
626*
627                     SRNAMT = 'ZPOSVXX'
628                     N_ERR_BNDS = 3
629                     CALL ZPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AFAC,
630     $                    LDA, EQUED, S, B, LDA, X,
631     $                    LDA, RCOND, RPVGRW_SVXX, BERR, N_ERR_BNDS,
632     $                    ERRBNDS_N, ERRBNDS_C, 0, ZERO, WORK,
633     $                    RWORK( 2*NRHS+1 ), INFO )
634*
635*                    Check the error code from ZPOSVXX.
636*
637                     IF( INFO.EQ.N+1 ) GOTO 90
638                     IF( INFO.NE.IZERO ) THEN
639                        CALL ALAERH( PATH, 'ZPOSVXX', INFO, IZERO,
640     $                               FACT // UPLO, N, N, -1, -1, NRHS,
641     $                               IMAT, NFAIL, NERRS, NOUT )
642                        GO TO 90
643                     END IF
644*
645                     IF( INFO.EQ.0 ) THEN
646                        IF( .NOT.PREFAC ) THEN
647*
648*                          Reconstruct matrix from factors and compute
649*                          residual.
650*
651                           CALL ZPOT01( UPLO, N, A, LDA, AFAC, LDA,
652     $                                  RWORK( 2*NRHS+1 ), RESULT( 1 ) )
653                           K1 = 1
654                        ELSE
655                           K1 = 2
656                        END IF
657*
658*                       Compute residual of the computed solution.
659*
660                        CALL ZLACPY( 'Full', N, NRHS, BSAV, LDA, WORK,
661     $                               LDA )
662                        CALL ZPOT02( UPLO, N, NRHS, ASAV, LDA, X, LDA,
663     $                               WORK, LDA, RWORK( 2*NRHS+1 ),
664     $                               RESULT( 2 ) )
665*
666*                       Check solution from generated exact solution.
667*
668                        IF( NOFACT .OR. ( PREFAC .AND. LSAME( EQUED,
669     $                      'N' ) ) ) THEN
670                           CALL ZGET04( N, NRHS, X, LDA, XACT, LDA,
671     $                                  RCONDC, RESULT( 3 ) )
672                        ELSE
673                           CALL ZGET04( N, NRHS, X, LDA, XACT, LDA,
674     $                                  ROLDC, RESULT( 3 ) )
675                        END IF
676*
677*                       Check the error bounds from iterative
678*                       refinement.
679*
680                        CALL ZPOT05( UPLO, N, NRHS, ASAV, LDA, B, LDA,
681     $                               X, LDA, XACT, LDA, RWORK,
682     $                               RWORK( NRHS+1 ), RESULT( 4 ) )
683                     ELSE
684                        K1 = 6
685                     END IF
686*
687*                    Compare RCOND from ZPOSVXX with the computed value
688*                    in RCONDC.
689*
690                     RESULT( 6 ) = DGET06( RCOND, RCONDC )
691*
692*                    Print information about the tests that did not pass
693*                    the threshold.
694*
695                     DO 85 K = K1, 6
696                        IF( RESULT( K ).GE.THRESH ) THEN
697                           IF( NFAIL.EQ.0 .AND. NERRS.EQ.0 )
698     $                        CALL ALADHD( NOUT, PATH )
699                           IF( PREFAC ) THEN
700                              WRITE( NOUT, FMT = 9997 )'ZPOSVXX', FACT,
701     $                           UPLO, N, EQUED, IMAT, K, RESULT( K )
702                           ELSE
703                              WRITE( NOUT, FMT = 9998 )'ZPOSVXX', FACT,
704     $                           UPLO, N, IMAT, K, RESULT( K )
705                           END IF
706                           NFAIL = NFAIL + 1
707                        END IF
708   85                CONTINUE
709                     NRUN = NRUN + 7 - K1
710   90             CONTINUE
711  100          CONTINUE
712  110       CONTINUE
713  120    CONTINUE
714  130 CONTINUE
715*
716*     Print a summary of the results.
717*
718      CALL ALASVM( PATH, NOUT, NFAIL, NRUN, NERRS )
719*
720
721*     Test Error Bounds for ZGESVXX
722
723      CALL ZEBCHVXX(THRESH, PATH)
724
725 9999 FORMAT( 1X, A, ', UPLO=''', A1, ''', N =', I5, ', type ', I1,
726     $      ', test(', I1, ')=', G12.5 )
727 9998 FORMAT( 1X, A, ', FACT=''', A1, ''', UPLO=''', A1, ''', N=', I5,
728     $      ', type ', I1, ', test(', I1, ')=', G12.5 )
729 9997 FORMAT( 1X, A, ', FACT=''', A1, ''', UPLO=''', A1, ''', N=', I5,
730     $      ', EQUED=''', A1, ''', type ', I1, ', test(', I1, ') =',
731     $      G12.5 )
732      RETURN
733*
734*     End of ZDRVPOX
735*
736      END
737