1*> \brief <b> ZPOSVXX computes the solution to system of linear equations A * X = B for PO matrices</b>
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZPOSVXX + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zposvxx.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zposvxx.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zposvxx.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18*  Definition:
19*  ===========
20*
21*       SUBROUTINE ZPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
22*                           S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
23*                           N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
24*                           NPARAMS, PARAMS, WORK, RWORK, INFO )
25*
26*       .. Scalar Arguments ..
27*       CHARACTER          EQUED, FACT, UPLO
28*       INTEGER            INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
29*      $                   N_ERR_BNDS
30*       DOUBLE PRECISION   RCOND, RPVGRW
31*       ..
32*       .. Array Arguments ..
33*       COMPLEX*16         A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
34*      $                   WORK( * ), X( LDX, * )
35*       DOUBLE PRECISION   S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
36*      $                   ERR_BNDS_NORM( NRHS, * ),
37*      $                   ERR_BNDS_COMP( NRHS, * )
38*       ..
39*
40*
41*> \par Purpose:
42*  =============
43*>
44*> \verbatim
45*>
46*>    ZPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
47*>    to compute the solution to a complex*16 system of linear equations
48*>    A * X = B, where A is an N-by-N Hermitian positive definite matrix
49*>    and X and B are N-by-NRHS matrices.
50*>
51*>    If requested, both normwise and maximum componentwise error bounds
52*>    are returned. ZPOSVXX will return a solution with a tiny
53*>    guaranteed error (O(eps) where eps is the working machine
54*>    precision) unless the matrix is very ill-conditioned, in which
55*>    case a warning is returned. Relevant condition numbers also are
56*>    calculated and returned.
57*>
58*>    ZPOSVXX accepts user-provided factorizations and equilibration
59*>    factors; see the definitions of the FACT and EQUED options.
60*>    Solving with refinement and using a factorization from a previous
61*>    ZPOSVXX call will also produce a solution with either O(eps)
62*>    errors or warnings, but we cannot make that claim for general
63*>    user-provided factorizations and equilibration factors if they
64*>    differ from what ZPOSVXX would itself produce.
65*> \endverbatim
66*
67*> \par Description:
68*  =================
69*>
70*> \verbatim
71*>
72*>    The following steps are performed:
73*>
74*>    1. If FACT = 'E', double precision scaling factors are computed to equilibrate
75*>    the system:
76*>
77*>      diag(S)*A*diag(S)     *inv(diag(S))*X = diag(S)*B
78*>
79*>    Whether or not the system will be equilibrated depends on the
80*>    scaling of the matrix A, but if equilibration is used, A is
81*>    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
82*>
83*>    2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
84*>    factor the matrix A (after equilibration if FACT = 'E') as
85*>       A = U**T* U,  if UPLO = 'U', or
86*>       A = L * L**T,  if UPLO = 'L',
87*>    where U is an upper triangular matrix and L is a lower triangular
88*>    matrix.
89*>
90*>    3. If the leading i-by-i principal minor is not positive definite,
91*>    then the routine returns with INFO = i. Otherwise, the factored
92*>    form of A is used to estimate the condition number of the matrix
93*>    A (see argument RCOND).  If the reciprocal of the condition number
94*>    is less than machine precision, the routine still goes on to solve
95*>    for X and compute error bounds as described below.
96*>
97*>    4. The system of equations is solved for X using the factored form
98*>    of A.
99*>
100*>    5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
101*>    the routine will use iterative refinement to try to get a small
102*>    error and error bounds.  Refinement calculates the residual to at
103*>    least twice the working precision.
104*>
105*>    6. If equilibration was used, the matrix X is premultiplied by
106*>    diag(S) so that it solves the original system before
107*>    equilibration.
108*> \endverbatim
109*
110*  Arguments:
111*  ==========
112*
113*> \verbatim
114*>     Some optional parameters are bundled in the PARAMS array.  These
115*>     settings determine how refinement is performed, but often the
116*>     defaults are acceptable.  If the defaults are acceptable, users
117*>     can pass NPARAMS = 0 which prevents the source code from accessing
118*>     the PARAMS argument.
119*> \endverbatim
120*>
121*> \param[in] FACT
122*> \verbatim
123*>          FACT is CHARACTER*1
124*>     Specifies whether or not the factored form of the matrix A is
125*>     supplied on entry, and if not, whether the matrix A should be
126*>     equilibrated before it is factored.
127*>       = 'F':  On entry, AF contains the factored form of A.
128*>               If EQUED is not 'N', the matrix A has been
129*>               equilibrated with scaling factors given by S.
130*>               A and AF are not modified.
131*>       = 'N':  The matrix A will be copied to AF and factored.
132*>       = 'E':  The matrix A will be equilibrated if necessary, then
133*>               copied to AF and factored.
134*> \endverbatim
135*>
136*> \param[in] UPLO
137*> \verbatim
138*>          UPLO is CHARACTER*1
139*>       = 'U':  Upper triangle of A is stored;
140*>       = 'L':  Lower triangle of A is stored.
141*> \endverbatim
142*>
143*> \param[in] N
144*> \verbatim
145*>          N is INTEGER
146*>     The number of linear equations, i.e., the order of the
147*>     matrix A.  N >= 0.
148*> \endverbatim
149*>
150*> \param[in] NRHS
151*> \verbatim
152*>          NRHS is INTEGER
153*>     The number of right hand sides, i.e., the number of columns
154*>     of the matrices B and X.  NRHS >= 0.
155*> \endverbatim
156*>
157*> \param[in,out] A
158*> \verbatim
159*>          A is COMPLEX*16 array, dimension (LDA,N)
160*>     On entry, the Hermitian matrix A, except if FACT = 'F' and EQUED =
161*>     'Y', then A must contain the equilibrated matrix
162*>     diag(S)*A*diag(S).  If UPLO = 'U', the leading N-by-N upper
163*>     triangular part of A contains the upper triangular part of the
164*>     matrix A, and the strictly lower triangular part of A is not
165*>     referenced.  If UPLO = 'L', the leading N-by-N lower triangular
166*>     part of A contains the lower triangular part of the matrix A, and
167*>     the strictly upper triangular part of A is not referenced.  A is
168*>     not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
169*>     'N' on exit.
170*>
171*>     On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
172*>     diag(S)*A*diag(S).
173*> \endverbatim
174*>
175*> \param[in] LDA
176*> \verbatim
177*>          LDA is INTEGER
178*>     The leading dimension of the array A.  LDA >= max(1,N).
179*> \endverbatim
180*>
181*> \param[in,out] AF
182*> \verbatim
183*>          AF is COMPLEX*16 array, dimension (LDAF,N)
184*>     If FACT = 'F', then AF is an input argument and on entry
185*>     contains the triangular factor U or L from the Cholesky
186*>     factorization A = U**T*U or A = L*L**T, in the same storage
187*>     format as A.  If EQUED .ne. 'N', then AF is the factored
188*>     form of the equilibrated matrix diag(S)*A*diag(S).
189*>
190*>     If FACT = 'N', then AF is an output argument and on exit
191*>     returns the triangular factor U or L from the Cholesky
192*>     factorization A = U**T*U or A = L*L**T of the original
193*>     matrix A.
194*>
195*>     If FACT = 'E', then AF is an output argument and on exit
196*>     returns the triangular factor U or L from the Cholesky
197*>     factorization A = U**T*U or A = L*L**T of the equilibrated
198*>     matrix A (see the description of A for the form of the
199*>     equilibrated matrix).
200*> \endverbatim
201*>
202*> \param[in] LDAF
203*> \verbatim
204*>          LDAF is INTEGER
205*>     The leading dimension of the array AF.  LDAF >= max(1,N).
206*> \endverbatim
207*>
208*> \param[in,out] EQUED
209*> \verbatim
210*>          EQUED is CHARACTER*1
211*>     Specifies the form of equilibration that was done.
212*>       = 'N':  No equilibration (always true if FACT = 'N').
213*>       = 'Y':  Both row and column equilibration, i.e., A has been
214*>               replaced by diag(S) * A * diag(S).
215*>     EQUED is an input argument if FACT = 'F'; otherwise, it is an
216*>     output argument.
217*> \endverbatim
218*>
219*> \param[in,out] S
220*> \verbatim
221*>          S is DOUBLE PRECISION array, dimension (N)
222*>     The row scale factors for A.  If EQUED = 'Y', A is multiplied on
223*>     the left and right by diag(S).  S is an input argument if FACT =
224*>     'F'; otherwise, S is an output argument.  If FACT = 'F' and EQUED
225*>     = 'Y', each element of S must be positive.  If S is output, each
226*>     element of S is a power of the radix. If S is input, each element
227*>     of S should be a power of the radix to ensure a reliable solution
228*>     and error estimates. Scaling by powers of the radix does not cause
229*>     rounding errors unless the result underflows or overflows.
230*>     Rounding errors during scaling lead to refining with a matrix that
231*>     is not equivalent to the input matrix, producing error estimates
232*>     that may not be reliable.
233*> \endverbatim
234*>
235*> \param[in,out] B
236*> \verbatim
237*>          B is COMPLEX*16 array, dimension (LDB,NRHS)
238*>     On entry, the N-by-NRHS right hand side matrix B.
239*>     On exit,
240*>     if EQUED = 'N', B is not modified;
241*>     if EQUED = 'Y', B is overwritten by diag(S)*B;
242*> \endverbatim
243*>
244*> \param[in] LDB
245*> \verbatim
246*>          LDB is INTEGER
247*>     The leading dimension of the array B.  LDB >= max(1,N).
248*> \endverbatim
249*>
250*> \param[out] X
251*> \verbatim
252*>          X is COMPLEX*16 array, dimension (LDX,NRHS)
253*>     If INFO = 0, the N-by-NRHS solution matrix X to the original
254*>     system of equations.  Note that A and B are modified on exit if
255*>     EQUED .ne. 'N', and the solution to the equilibrated system is
256*>     inv(diag(S))*X.
257*> \endverbatim
258*>
259*> \param[in] LDX
260*> \verbatim
261*>          LDX is INTEGER
262*>     The leading dimension of the array X.  LDX >= max(1,N).
263*> \endverbatim
264*>
265*> \param[out] RCOND
266*> \verbatim
267*>          RCOND is DOUBLE PRECISION
268*>     Reciprocal scaled condition number.  This is an estimate of the
269*>     reciprocal Skeel condition number of the matrix A after
270*>     equilibration (if done).  If this is less than the machine
271*>     precision (in particular, if it is zero), the matrix is singular
272*>     to working precision.  Note that the error may still be small even
273*>     if this number is very small and the matrix appears ill-
274*>     conditioned.
275*> \endverbatim
276*>
277*> \param[out] RPVGRW
278*> \verbatim
279*>          RPVGRW is DOUBLE PRECISION
280*>     Reciprocal pivot growth.  On exit, this contains the reciprocal
281*>     pivot growth factor norm(A)/norm(U). The "max absolute element"
282*>     norm is used.  If this is much less than 1, then the stability of
283*>     the LU factorization of the (equilibrated) matrix A could be poor.
284*>     This also means that the solution X, estimated condition numbers,
285*>     and error bounds could be unreliable. If factorization fails with
286*>     0<INFO<=N, then this contains the reciprocal pivot growth factor
287*>     for the leading INFO columns of A.
288*> \endverbatim
289*>
290*> \param[out] BERR
291*> \verbatim
292*>          BERR is DOUBLE PRECISION array, dimension (NRHS)
293*>     Componentwise relative backward error.  This is the
294*>     componentwise relative backward error of each solution vector X(j)
295*>     (i.e., the smallest relative change in any element of A or B that
296*>     makes X(j) an exact solution).
297*> \endverbatim
298*>
299*> \param[in] N_ERR_BNDS
300*> \verbatim
301*>          N_ERR_BNDS is INTEGER
302*>     Number of error bounds to return for each right hand side
303*>     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
304*>     ERR_BNDS_COMP below.
305*> \endverbatim
306*>
307*> \param[out] ERR_BNDS_NORM
308*> \verbatim
309*>          ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
310*>     For each right-hand side, this array contains information about
311*>     various error bounds and condition numbers corresponding to the
312*>     normwise relative error, which is defined as follows:
313*>
314*>     Normwise relative error in the ith solution vector:
315*>             max_j (abs(XTRUE(j,i) - X(j,i)))
316*>            ------------------------------
317*>                  max_j abs(X(j,i))
318*>
319*>     The array is indexed by the type of error information as described
320*>     below. There currently are up to three pieces of information
321*>     returned.
322*>
323*>     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
324*>     right-hand side.
325*>
326*>     The second index in ERR_BNDS_NORM(:,err) contains the following
327*>     three fields:
328*>     err = 1 "Trust/don't trust" boolean. Trust the answer if the
329*>              reciprocal condition number is less than the threshold
330*>              sqrt(n) * dlamch('Epsilon').
331*>
332*>     err = 2 "Guaranteed" error bound: The estimated forward error,
333*>              almost certainly within a factor of 10 of the true error
334*>              so long as the next entry is greater than the threshold
335*>              sqrt(n) * dlamch('Epsilon'). This error bound should only
336*>              be trusted if the previous boolean is true.
337*>
338*>     err = 3  Reciprocal condition number: Estimated normwise
339*>              reciprocal condition number.  Compared with the threshold
340*>              sqrt(n) * dlamch('Epsilon') to determine if the error
341*>              estimate is "guaranteed". These reciprocal condition
342*>              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
343*>              appropriately scaled matrix Z.
344*>              Let Z = S*A, where S scales each row by a power of the
345*>              radix so all absolute row sums of Z are approximately 1.
346*>
347*>     See Lapack Working Note 165 for further details and extra
348*>     cautions.
349*> \endverbatim
350*>
351*> \param[out] ERR_BNDS_COMP
352*> \verbatim
353*>          ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
354*>     For each right-hand side, this array contains information about
355*>     various error bounds and condition numbers corresponding to the
356*>     componentwise relative error, which is defined as follows:
357*>
358*>     Componentwise relative error in the ith solution vector:
359*>                    abs(XTRUE(j,i) - X(j,i))
360*>             max_j ----------------------
361*>                         abs(X(j,i))
362*>
363*>     The array is indexed by the right-hand side i (on which the
364*>     componentwise relative error depends), and the type of error
365*>     information as described below. There currently are up to three
366*>     pieces of information returned for each right-hand side. If
367*>     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
368*>     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS < 3, then at most
369*>     the first (:,N_ERR_BNDS) entries are returned.
370*>
371*>     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
372*>     right-hand side.
373*>
374*>     The second index in ERR_BNDS_COMP(:,err) contains the following
375*>     three fields:
376*>     err = 1 "Trust/don't trust" boolean. Trust the answer if the
377*>              reciprocal condition number is less than the threshold
378*>              sqrt(n) * dlamch('Epsilon').
379*>
380*>     err = 2 "Guaranteed" error bound: The estimated forward error,
381*>              almost certainly within a factor of 10 of the true error
382*>              so long as the next entry is greater than the threshold
383*>              sqrt(n) * dlamch('Epsilon'). This error bound should only
384*>              be trusted if the previous boolean is true.
385*>
386*>     err = 3  Reciprocal condition number: Estimated componentwise
387*>              reciprocal condition number.  Compared with the threshold
388*>              sqrt(n) * dlamch('Epsilon') to determine if the error
389*>              estimate is "guaranteed". These reciprocal condition
390*>              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
391*>              appropriately scaled matrix Z.
392*>              Let Z = S*(A*diag(x)), where x is the solution for the
393*>              current right-hand side and S scales each row of
394*>              A*diag(x) by a power of the radix so all absolute row
395*>              sums of Z are approximately 1.
396*>
397*>     See Lapack Working Note 165 for further details and extra
398*>     cautions.
399*> \endverbatim
400*>
401*> \param[in] NPARAMS
402*> \verbatim
403*>          NPARAMS is INTEGER
404*>     Specifies the number of parameters set in PARAMS.  If <= 0, the
405*>     PARAMS array is never referenced and default values are used.
406*> \endverbatim
407*>
408*> \param[in,out] PARAMS
409*> \verbatim
410*>          PARAMS is DOUBLE PRECISION array, dimension NPARAMS
411*>     Specifies algorithm parameters.  If an entry is < 0.0, then
412*>     that entry will be filled with default value used for that
413*>     parameter.  Only positions up to NPARAMS are accessed; defaults
414*>     are used for higher-numbered parameters.
415*>
416*>       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
417*>            refinement or not.
418*>         Default: 1.0D+0
419*>            = 0.0:  No refinement is performed, and no error bounds are
420*>                    computed.
421*>            = 1.0:  Use the extra-precise refinement algorithm.
422*>              (other values are reserved for future use)
423*>
424*>       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
425*>            computations allowed for refinement.
426*>         Default: 10
427*>         Aggressive: Set to 100 to permit convergence using approximate
428*>                     factorizations or factorizations other than LU. If
429*>                     the factorization uses a technique other than
430*>                     Gaussian elimination, the guarantees in
431*>                     err_bnds_norm and err_bnds_comp may no longer be
432*>                     trustworthy.
433*>
434*>       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
435*>            will attempt to find a solution with small componentwise
436*>            relative error in the double-precision algorithm.  Positive
437*>            is true, 0.0 is false.
438*>         Default: 1.0 (attempt componentwise convergence)
439*> \endverbatim
440*>
441*> \param[out] WORK
442*> \verbatim
443*>          WORK is COMPLEX*16 array, dimension (2*N)
444*> \endverbatim
445*>
446*> \param[out] RWORK
447*> \verbatim
448*>          RWORK is DOUBLE PRECISION array, dimension (2*N)
449*> \endverbatim
450*>
451*> \param[out] INFO
452*> \verbatim
453*>          INFO is INTEGER
454*>       = 0:  Successful exit. The solution to every right-hand side is
455*>         guaranteed.
456*>       < 0:  If INFO = -i, the i-th argument had an illegal value
457*>       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
458*>         has been completed, but the factor U is exactly singular, so
459*>         the solution and error bounds could not be computed. RCOND = 0
460*>         is returned.
461*>       = N+J: The solution corresponding to the Jth right-hand side is
462*>         not guaranteed. The solutions corresponding to other right-
463*>         hand sides K with K > J may not be guaranteed as well, but
464*>         only the first such right-hand side is reported. If a small
465*>         componentwise error is not requested (PARAMS(3) = 0.0) then
466*>         the Jth right-hand side is the first with a normwise error
467*>         bound that is not guaranteed (the smallest J such
468*>         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
469*>         the Jth right-hand side is the first with either a normwise or
470*>         componentwise error bound that is not guaranteed (the smallest
471*>         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
472*>         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
473*>         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
474*>         about all of the right-hand sides check ERR_BNDS_NORM or
475*>         ERR_BNDS_COMP.
476*> \endverbatim
477*
478*  Authors:
479*  ========
480*
481*> \author Univ. of Tennessee
482*> \author Univ. of California Berkeley
483*> \author Univ. of Colorado Denver
484*> \author NAG Ltd.
485*
486*> \ingroup complex16POsolve
487*
488*  =====================================================================
489      SUBROUTINE ZPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
490     $                    S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
491     $                    N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
492     $                    NPARAMS, PARAMS, WORK, RWORK, INFO )
493*
494*  -- LAPACK driver routine --
495*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
496*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
497*
498*     .. Scalar Arguments ..
499      CHARACTER          EQUED, FACT, UPLO
500      INTEGER            INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
501     $                   N_ERR_BNDS
502      DOUBLE PRECISION   RCOND, RPVGRW
503*     ..
504*     .. Array Arguments ..
505      COMPLEX*16         A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
506     $                   WORK( * ), X( LDX, * )
507      DOUBLE PRECISION   S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
508     $                   ERR_BNDS_NORM( NRHS, * ),
509     $                   ERR_BNDS_COMP( NRHS, * )
510*     ..
511*
512*  ==================================================================
513*
514*     .. Parameters ..
515      DOUBLE PRECISION   ZERO, ONE
516      PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
517      INTEGER            FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
518      INTEGER            RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
519      INTEGER            CMP_ERR_I, PIV_GROWTH_I
520      PARAMETER          ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
521     $                   BERR_I = 3 )
522      PARAMETER          ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
523      PARAMETER          ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
524     $                   PIV_GROWTH_I = 9 )
525*     ..
526*     .. Local Scalars ..
527      LOGICAL            EQUIL, NOFACT, RCEQU
528      INTEGER            INFEQU, J
529      DOUBLE PRECISION   AMAX, BIGNUM, SMIN, SMAX, SCOND, SMLNUM
530*     ..
531*     .. External Functions ..
532      EXTERNAL           LSAME, DLAMCH, ZLA_PORPVGRW
533      LOGICAL            LSAME
534      DOUBLE PRECISION   DLAMCH, ZLA_PORPVGRW
535*     ..
536*     .. External Subroutines ..
537      EXTERNAL           ZPOEQUB, ZPOTRF, ZPOTRS, ZLACPY,
538     $                   ZLAQHE, XERBLA, ZLASCL2, ZPORFSX
539*     ..
540*     .. Intrinsic Functions ..
541      INTRINSIC          MAX, MIN
542*     ..
543*     .. Executable Statements ..
544*
545      INFO = 0
546      NOFACT = LSAME( FACT, 'N' )
547      EQUIL = LSAME( FACT, 'E' )
548      SMLNUM = DLAMCH( 'Safe minimum' )
549      BIGNUM = ONE / SMLNUM
550      IF( NOFACT .OR. EQUIL ) THEN
551         EQUED = 'N'
552         RCEQU = .FALSE.
553      ELSE
554         RCEQU = LSAME( EQUED, 'Y' )
555      ENDIF
556*
557*     Default is failure.  If an input parameter is wrong or
558*     factorization fails, make everything look horrible.  Only the
559*     pivot growth is set here, the rest is initialized in ZPORFSX.
560*
561      RPVGRW = ZERO
562*
563*     Test the input parameters.  PARAMS is not tested until ZPORFSX.
564*
565      IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
566     $     LSAME( FACT, 'F' ) ) THEN
567         INFO = -1
568      ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND.
569     $         .NOT.LSAME( UPLO, 'L' ) ) THEN
570         INFO = -2
571      ELSE IF( N.LT.0 ) THEN
572         INFO = -3
573      ELSE IF( NRHS.LT.0 ) THEN
574         INFO = -4
575      ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
576         INFO = -6
577      ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
578         INFO = -8
579      ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
580     $        ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
581         INFO = -9
582      ELSE
583         IF ( RCEQU ) THEN
584            SMIN = BIGNUM
585            SMAX = ZERO
586            DO 10 J = 1, N
587               SMIN = MIN( SMIN, S( J ) )
588               SMAX = MAX( SMAX, S( J ) )
589 10         CONTINUE
590            IF( SMIN.LE.ZERO ) THEN
591               INFO = -10
592            ELSE IF( N.GT.0 ) THEN
593               SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
594            ELSE
595               SCOND = ONE
596            END IF
597         END IF
598         IF( INFO.EQ.0 ) THEN
599            IF( LDB.LT.MAX( 1, N ) ) THEN
600               INFO = -12
601            ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
602               INFO = -14
603            END IF
604         END IF
605      END IF
606*
607      IF( INFO.NE.0 ) THEN
608         CALL XERBLA( 'ZPOSVXX', -INFO )
609         RETURN
610      END IF
611*
612      IF( EQUIL ) THEN
613*
614*     Compute row and column scalings to equilibrate the matrix A.
615*
616         CALL ZPOEQUB( N, A, LDA, S, SCOND, AMAX, INFEQU )
617         IF( INFEQU.EQ.0 ) THEN
618*
619*     Equilibrate the matrix.
620*
621            CALL ZLAQHE( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED )
622            RCEQU = LSAME( EQUED, 'Y' )
623         END IF
624      END IF
625*
626*     Scale the right-hand side.
627*
628      IF( RCEQU ) CALL ZLASCL2( N, NRHS, S, B, LDB )
629*
630      IF( NOFACT .OR. EQUIL ) THEN
631*
632*        Compute the Cholesky factorization of A.
633*
634         CALL ZLACPY( UPLO, N, N, A, LDA, AF, LDAF )
635         CALL ZPOTRF( UPLO, N, AF, LDAF, INFO )
636*
637*        Return if INFO is non-zero.
638*
639         IF( INFO.GT.0 ) THEN
640*
641*           Pivot in column INFO is exactly 0
642*           Compute the reciprocal pivot growth factor of the
643*           leading rank-deficient INFO columns of A.
644*
645            RPVGRW = ZLA_PORPVGRW( UPLO, N, A, LDA, AF, LDAF, RWORK )
646            RETURN
647         END IF
648      END IF
649*
650*     Compute the reciprocal pivot growth factor RPVGRW.
651*
652      RPVGRW = ZLA_PORPVGRW( UPLO, N, A, LDA, AF, LDAF, RWORK )
653*
654*     Compute the solution matrix X.
655*
656      CALL ZLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
657      CALL ZPOTRS( UPLO, N, NRHS, AF, LDAF, X, LDX, INFO )
658*
659*     Use iterative refinement to improve the computed solution and
660*     compute error bounds and backward error estimates for it.
661*
662      CALL ZPORFSX( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF,
663     $     S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM,
664     $     ERR_BNDS_COMP,  NPARAMS, PARAMS, WORK, RWORK, INFO )
665
666*
667*     Scale solutions.
668*
669      IF ( RCEQU ) THEN
670         CALL ZLASCL2( N, NRHS, S, X, LDX )
671      END IF
672*
673      RETURN
674*
675*     End of ZPOSVXX
676*
677      END
678