001:       SUBROUTINE DPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
002:      $                    S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
003:      $                    N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
004:      $                    NPARAMS, PARAMS, WORK, IWORK, INFO )
005: *
006: *     -- LAPACK driver routine (version 3.2)                          --
007: *     -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and --
008: *     -- Jason Riedy of Univ. of California Berkeley.                 --
009: *     -- November 2008                                                --
010: *
011: *     -- LAPACK is a software package provided by Univ. of Tennessee, --
012: *     -- Univ. of California Berkeley and NAG Ltd.                    --
013: *
014:       IMPLICIT NONE
015: *     ..
016: *     .. Scalar Arguments ..
017:       CHARACTER          EQUED, FACT, UPLO
018:       INTEGER            INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
019:      $                   N_ERR_BNDS
020:       DOUBLE PRECISION   RCOND, RPVGRW
021: *     ..
022: *     .. Array Arguments ..
023:       INTEGER            IWORK( * )
024:       DOUBLE PRECISION   A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
025:      $                   X( LDX, * ), WORK( * )
026:       DOUBLE PRECISION   S( * ), PARAMS( * ), BERR( * ),
027:      $                   ERR_BNDS_NORM( NRHS, * ),
028:      $                   ERR_BNDS_COMP( NRHS, * )
029: *     ..
030: *
031: *     Purpose
032: *     =======
033: *
034: *     DPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
035: *     to compute the solution to a double precision system of linear equations
036: *     A * X = B, where A is an N-by-N symmetric positive definite matrix
037: *     and X and B are N-by-NRHS matrices.
038: *
039: *     If requested, both normwise and maximum componentwise error bounds
040: *     are returned. DPOSVXX will return a solution with a tiny
041: *     guaranteed error (O(eps) where eps is the working machine
042: *     precision) unless the matrix is very ill-conditioned, in which
043: *     case a warning is returned. Relevant condition numbers also are
044: *     calculated and returned.
045: *
046: *     DPOSVXX accepts user-provided factorizations and equilibration
047: *     factors; see the definitions of the FACT and EQUED options.
048: *     Solving with refinement and using a factorization from a previous
049: *     DPOSVXX call will also produce a solution with either O(eps)
050: *     errors or warnings, but we cannot make that claim for general
051: *     user-provided factorizations and equilibration factors if they
052: *     differ from what DPOSVXX would itself produce.
053: *
054: *     Description
055: *     ===========
056: *
057: *     The following steps are performed:
058: *
059: *     1. If FACT = 'E', double precision scaling factors are computed to equilibrate
060: *     the system:
061: *
062: *       diag(S)*A*diag(S)     *inv(diag(S))*X = diag(S)*B
063: *
064: *     Whether or not the system will be equilibrated depends on the
065: *     scaling of the matrix A, but if equilibration is used, A is
066: *     overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
067: *
068: *     2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
069: *     factor the matrix A (after equilibration if FACT = 'E') as
070: *        A = U**T* U,  if UPLO = 'U', or
071: *        A = L * L**T,  if UPLO = 'L',
072: *     where U is an upper triangular matrix and L is a lower triangular
073: *     matrix.
074: *
075: *     3. If the leading i-by-i principal minor is not positive definite,
076: *     then the routine returns with INFO = i. Otherwise, the factored
077: *     form of A is used to estimate the condition number of the matrix
078: *     A (see argument RCOND).  If the reciprocal of the condition number
079: *     is less than machine precision, the routine still goes on to solve
080: *     for X and compute error bounds as described below.
081: *
082: *     4. The system of equations is solved for X using the factored form
083: *     of A.
084: *
085: *     5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
086: *     the routine will use iterative refinement to try to get a small
087: *     error and error bounds.  Refinement calculates the residual to at
088: *     least twice the working precision.
089: *
090: *     6. If equilibration was used, the matrix X is premultiplied by
091: *     diag(S) so that it solves the original system before
092: *     equilibration.
093: *
094: *     Arguments
095: *     =========
096: *
097: *     Some optional parameters are bundled in the PARAMS array.  These
098: *     settings determine how refinement is performed, but often the
099: *     defaults are acceptable.  If the defaults are acceptable, users
100: *     can pass NPARAMS = 0 which prevents the source code from accessing
101: *     the PARAMS argument.
102: *
103: *     FACT    (input) CHARACTER*1
104: *     Specifies whether or not the factored form of the matrix A is
105: *     supplied on entry, and if not, whether the matrix A should be
106: *     equilibrated before it is factored.
107: *       = 'F':  On entry, AF contains the factored form of A.
108: *               If EQUED is not 'N', the matrix A has been
109: *               equilibrated with scaling factors given by S.
110: *               A and AF are not modified.
111: *       = 'N':  The matrix A will be copied to AF and factored.
112: *       = 'E':  The matrix A will be equilibrated if necessary, then
113: *               copied to AF and factored.
114: *
115: *     UPLO    (input) CHARACTER*1
116: *       = 'U':  Upper triangle of A is stored;
117: *       = 'L':  Lower triangle of A is stored.
118: *
119: *     N       (input) INTEGER
120: *     The number of linear equations, i.e., the order of the
121: *     matrix A.  N >= 0.
122: *
123: *     NRHS    (input) INTEGER
124: *     The number of right hand sides, i.e., the number of columns
125: *     of the matrices B and X.  NRHS >= 0.
126: *
127: *     A       (input/output) DOUBLE PRECISION array, dimension (LDA,N)
128: *     On entry, the symmetric matrix A, except if FACT = 'F' and EQUED =
129: *     'Y', then A must contain the equilibrated matrix
130: *     diag(S)*A*diag(S).  If UPLO = 'U', the leading N-by-N upper
131: *     triangular part of A contains the upper triangular part of the
132: *     matrix A, and the strictly lower triangular part of A is not
133: *     referenced.  If UPLO = 'L', the leading N-by-N lower triangular
134: *     part of A contains the lower triangular part of the matrix A, and
135: *     the strictly upper triangular part of A is not referenced.  A is
136: *     not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
137: *     'N' on exit.
138: *
139: *     On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
140: *     diag(S)*A*diag(S).
141: *
142: *     LDA     (input) INTEGER
143: *     The leading dimension of the array A.  LDA >= max(1,N).
144: *
145: *     AF      (input or output) DOUBLE PRECISION array, dimension (LDAF,N)
146: *     If FACT = 'F', then AF is an input argument and on entry
147: *     contains the triangular factor U or L from the Cholesky
148: *     factorization A = U**T*U or A = L*L**T, in the same storage
149: *     format as A.  If EQUED .ne. 'N', then AF is the factored
150: *     form of the equilibrated matrix diag(S)*A*diag(S).
151: *
152: *     If FACT = 'N', then AF is an output argument and on exit
153: *     returns the triangular factor U or L from the Cholesky
154: *     factorization A = U**T*U or A = L*L**T of the original
155: *     matrix A.
156: *
157: *     If FACT = 'E', then AF is an output argument and on exit
158: *     returns the triangular factor U or L from the Cholesky
159: *     factorization A = U**T*U or A = L*L**T of the equilibrated
160: *     matrix A (see the description of A for the form of the
161: *     equilibrated matrix).
162: *
163: *     LDAF    (input) INTEGER
164: *     The leading dimension of the array AF.  LDAF >= max(1,N).
165: *
166: *     EQUED   (input or output) CHARACTER*1
167: *     Specifies the form of equilibration that was done.
168: *       = 'N':  No equilibration (always true if FACT = 'N').
169: *       = 'Y':  Both row and column equilibration, i.e., A has been
170: *               replaced by diag(S) * A * diag(S).
171: *     EQUED is an input argument if FACT = 'F'; otherwise, it is an
172: *     output argument.
173: *
174: *     S       (input or output) DOUBLE PRECISION array, dimension (N)
175: *     The row scale factors for A.  If EQUED = 'Y', A is multiplied on
176: *     the left and right by diag(S).  S is an input argument if FACT =
177: *     'F'; otherwise, S is an output argument.  If FACT = 'F' and EQUED
178: *     = 'Y', each element of S must be positive.  If S is output, each
179: *     element of S is a power of the radix. If S is input, each element
180: *     of S should be a power of the radix to ensure a reliable solution
181: *     and error estimates. Scaling by powers of the radix does not cause
182: *     rounding errors unless the result underflows or overflows.
183: *     Rounding errors during scaling lead to refining with a matrix that
184: *     is not equivalent to the input matrix, producing error estimates
185: *     that may not be reliable.
186: *
187: *     B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
188: *     On entry, the N-by-NRHS right hand side matrix B.
189: *     On exit,
190: *     if EQUED = 'N', B is not modified;
191: *     if EQUED = 'Y', B is overwritten by diag(S)*B;
192: *
193: *     LDB     (input) INTEGER
194: *     The leading dimension of the array B.  LDB >= max(1,N).
195: *
196: *     X       (output) DOUBLE PRECISION array, dimension (LDX,NRHS)
197: *     If INFO = 0, the N-by-NRHS solution matrix X to the original
198: *     system of equations.  Note that A and B are modified on exit if
199: *     EQUED .ne. 'N', and the solution to the equilibrated system is
200: *     inv(diag(S))*X.
201: *
202: *     LDX     (input) INTEGER
203: *     The leading dimension of the array X.  LDX >= max(1,N).
204: *
205: *     RCOND   (output) DOUBLE PRECISION
206: *     Reciprocal scaled condition number.  This is an estimate of the
207: *     reciprocal Skeel condition number of the matrix A after
208: *     equilibration (if done).  If this is less than the machine
209: *     precision (in particular, if it is zero), the matrix is singular
210: *     to working precision.  Note that the error may still be small even
211: *     if this number is very small and the matrix appears ill-
212: *     conditioned.
213: *
214: *     RPVGRW  (output) DOUBLE PRECISION
215: *     Reciprocal pivot growth.  On exit, this contains the reciprocal
216: *     pivot growth factor norm(A)/norm(U). The "max absolute element"
217: *     norm is used.  If this is much less than 1, then the stability of
218: *     the LU factorization of the (equilibrated) matrix A could be poor.
219: *     This also means that the solution X, estimated condition numbers,
220: *     and error bounds could be unreliable. If factorization fails with
221: *     0<INFO<=N, then this contains the reciprocal pivot growth factor
222: *     for the leading INFO columns of A.
223: *
224: *     BERR    (output) DOUBLE PRECISION array, dimension (NRHS)
225: *     Componentwise relative backward error.  This is the
226: *     componentwise relative backward error of each solution vector X(j)
227: *     (i.e., the smallest relative change in any element of A or B that
228: *     makes X(j) an exact solution).
229: *
230: *     N_ERR_BNDS (input) INTEGER
231: *     Number of error bounds to return for each right hand side
232: *     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
233: *     ERR_BNDS_COMP below.
234: *
235: *     ERR_BNDS_NORM  (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
236: *     For each right-hand side, this array contains information about
237: *     various error bounds and condition numbers corresponding to the
238: *     normwise relative error, which is defined as follows:
239: *
240: *     Normwise relative error in the ith solution vector:
241: *             max_j (abs(XTRUE(j,i) - X(j,i)))
242: *            ------------------------------
243: *                  max_j abs(X(j,i))
244: *
245: *     The array is indexed by the type of error information as described
246: *     below. There currently are up to three pieces of information
247: *     returned.
248: *
249: *     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
250: *     right-hand side.
251: *
252: *     The second index in ERR_BNDS_NORM(:,err) contains the following
253: *     three fields:
254: *     err = 1 "Trust/don't trust" boolean. Trust the answer if the
255: *              reciprocal condition number is less than the threshold
256: *              sqrt(n) * dlamch('Epsilon').
257: *
258: *     err = 2 "Guaranteed" error bound: The estimated forward error,
259: *              almost certainly within a factor of 10 of the true error
260: *              so long as the next entry is greater than the threshold
261: *              sqrt(n) * dlamch('Epsilon'). This error bound should only
262: *              be trusted if the previous boolean is true.
263: *
264: *     err = 3  Reciprocal condition number: Estimated normwise
265: *              reciprocal condition number.  Compared with the threshold
266: *              sqrt(n) * dlamch('Epsilon') to determine if the error
267: *              estimate is "guaranteed". These reciprocal condition
268: *              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
269: *              appropriately scaled matrix Z.
270: *              Let Z = S*A, where S scales each row by a power of the
271: *              radix so all absolute row sums of Z are approximately 1.
272: *
273: *     See Lapack Working Note 165 for further details and extra
274: *     cautions.
275: *
276: *     ERR_BNDS_COMP  (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
277: *     For each right-hand side, this array contains information about
278: *     various error bounds and condition numbers corresponding to the
279: *     componentwise relative error, which is defined as follows:
280: *
281: *     Componentwise relative error in the ith solution vector:
282: *                    abs(XTRUE(j,i) - X(j,i))
283: *             max_j ----------------------
284: *                         abs(X(j,i))
285: *
286: *     The array is indexed by the right-hand side i (on which the
287: *     componentwise relative error depends), and the type of error
288: *     information as described below. There currently are up to three
289: *     pieces of information returned for each right-hand side. If
290: *     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
291: *     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS .LT. 3, then at most
292: *     the first (:,N_ERR_BNDS) entries are returned.
293: *
294: *     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
295: *     right-hand side.
296: *
297: *     The second index in ERR_BNDS_COMP(:,err) contains the following
298: *     three fields:
299: *     err = 1 "Trust/don't trust" boolean. Trust the answer if the
300: *              reciprocal condition number is less than the threshold
301: *              sqrt(n) * dlamch('Epsilon').
302: *
303: *     err = 2 "Guaranteed" error bound: The estimated forward error,
304: *              almost certainly within a factor of 10 of the true error
305: *              so long as the next entry is greater than the threshold
306: *              sqrt(n) * dlamch('Epsilon'). This error bound should only
307: *              be trusted if the previous boolean is true.
308: *
309: *     err = 3  Reciprocal condition number: Estimated componentwise
310: *              reciprocal condition number.  Compared with the threshold
311: *              sqrt(n) * dlamch('Epsilon') to determine if the error
312: *              estimate is "guaranteed". These reciprocal condition
313: *              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
314: *              appropriately scaled matrix Z.
315: *              Let Z = S*(A*diag(x)), where x is the solution for the
316: *              current right-hand side and S scales each row of
317: *              A*diag(x) by a power of the radix so all absolute row
318: *              sums of Z are approximately 1.
319: *
320: *     See Lapack Working Note 165 for further details and extra
321: *     cautions.
322: *
323: *     NPARAMS (input) INTEGER
324: *     Specifies the number of parameters set in PARAMS.  If .LE. 0, the
325: *     PARAMS array is never referenced and default values are used.
326: *
327: *     PARAMS  (input / output) DOUBLE PRECISION array, dimension NPARAMS
328: *     Specifies algorithm parameters.  If an entry is .LT. 0.0, then
329: *     that entry will be filled with default value used for that
330: *     parameter.  Only positions up to NPARAMS are accessed; defaults
331: *     are used for higher-numbered parameters.
332: *
333: *       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
334: *            refinement or not.
335: *         Default: 1.0D+0
336: *            = 0.0 : No refinement is performed, and no error bounds are
337: *                    computed.
338: *            = 1.0 : Use the extra-precise refinement algorithm.
339: *              (other values are reserved for future use)
340: *
341: *       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
342: *            computations allowed for refinement.
343: *         Default: 10
344: *         Aggressive: Set to 100 to permit convergence using approximate
345: *                     factorizations or factorizations other than LU. If
346: *                     the factorization uses a technique other than
347: *                     Gaussian elimination, the guarantees in
348: *                     err_bnds_norm and err_bnds_comp may no longer be
349: *                     trustworthy.
350: *
351: *       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
352: *            will attempt to find a solution with small componentwise
353: *            relative error in the double-precision algorithm.  Positive
354: *            is true, 0.0 is false.
355: *         Default: 1.0 (attempt componentwise convergence)
356: *
357: *     WORK    (workspace) DOUBLE PRECISION array, dimension (4*N)
358: *
359: *     IWORK   (workspace) INTEGER array, dimension (N)
360: *
361: *     INFO    (output) INTEGER
362: *       = 0:  Successful exit. The solution to every right-hand side is
363: *         guaranteed.
364: *       < 0:  If INFO = -i, the i-th argument had an illegal value
365: *       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
366: *         has been completed, but the factor U is exactly singular, so
367: *         the solution and error bounds could not be computed. RCOND = 0
368: *         is returned.
369: *       = N+J: The solution corresponding to the Jth right-hand side is
370: *         not guaranteed. The solutions corresponding to other right-
371: *         hand sides K with K > J may not be guaranteed as well, but
372: *         only the first such right-hand side is reported. If a small
373: *         componentwise error is not requested (PARAMS(3) = 0.0) then
374: *         the Jth right-hand side is the first with a normwise error
375: *         bound that is not guaranteed (the smallest J such
376: *         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
377: *         the Jth right-hand side is the first with either a normwise or
378: *         componentwise error bound that is not guaranteed (the smallest
379: *         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
380: *         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
381: *         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
382: *         about all of the right-hand sides check ERR_BNDS_NORM or
383: *         ERR_BNDS_COMP.
384: *
385: *     ==================================================================
386: *
387: *     .. Parameters ..
388:       DOUBLE PRECISION   ZERO, ONE
389:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
390:       INTEGER            FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
391:       INTEGER            RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
392:       INTEGER            CMP_ERR_I, PIV_GROWTH_I
393:       PARAMETER          ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
394:      $                   BERR_I = 3 )
395:       PARAMETER          ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
396:       PARAMETER          ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
397:      $                   PIV_GROWTH_I = 9 )
398: *     ..
399: *     .. Local Scalars ..
400:       LOGICAL            EQUIL, NOFACT, RCEQU
401:       INTEGER            INFEQU, J
402:       DOUBLE PRECISION   AMAX, BIGNUM, SMIN, SMAX,
403:      $                   SCOND, SMLNUM
404: *     ..
405: *     .. External Functions ..
406:       EXTERNAL           LSAME, DLAMCH, DLA_PORPVGRW
407:       LOGICAL            LSAME
408:       DOUBLE PRECISION   DLAMCH, DLA_PORPVGRW
409: *     ..
410: *     .. External Subroutines ..
411:       EXTERNAL           DPOEQUB, DPOTRF, DPOTRS, DLACPY, DLAQSY,
412:      $                   XERBLA, DLASCL2, DPORFSX
413: *     ..
414: *     .. Intrinsic Functions ..
415:       INTRINSIC          MAX, MIN
416: *     ..
417: *     .. Executable Statements ..
418: *
419:       INFO = 0
420:       NOFACT = LSAME( FACT, 'N' )
421:       EQUIL = LSAME( FACT, 'E' )
422:       SMLNUM = DLAMCH( 'Safe minimum' )
423:       BIGNUM = ONE / SMLNUM
424:       IF( NOFACT .OR. EQUIL ) THEN
425:          EQUED = 'N'
426:          RCEQU = .FALSE.
427:       ELSE
428:          RCEQU = LSAME( EQUED, 'Y' )
429:       ENDIF
430: *
431: *     Default is failure.  If an input parameter is wrong or
432: *     factorization fails, make everything look horrible.  Only the
433: *     pivot growth is set here, the rest is initialized in DPORFSX.
434: *
435:       RPVGRW = ZERO
436: *
437: *     Test the input parameters.  PARAMS is not tested until DPORFSX.
438: *
439:       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
440:      $     LSAME( FACT, 'F' ) ) THEN
441:          INFO = -1
442:       ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND.
443:      $         .NOT.LSAME( UPLO, 'L' ) ) THEN
444:          INFO = -2
445:       ELSE IF( N.LT.0 ) THEN
446:          INFO = -3
447:       ELSE IF( NRHS.LT.0 ) THEN
448:          INFO = -4
449:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
450:          INFO = -6
451:       ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
452:          INFO = -8
453:       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
454:      $        ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
455:          INFO = -9
456:       ELSE
457:          IF ( RCEQU ) THEN
458:             SMIN = BIGNUM
459:             SMAX = ZERO
460:             DO 10 J = 1, N
461:                SMIN = MIN( SMIN, S( J ) )
462:                SMAX = MAX( SMAX, S( J ) )
463:  10         CONTINUE
464:             IF( SMIN.LE.ZERO ) THEN
465:                INFO = -10
466:             ELSE IF( N.GT.0 ) THEN
467:                SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
468:             ELSE
469:                SCOND = ONE
470:             END IF
471:          END IF
472:          IF( INFO.EQ.0 ) THEN
473:             IF( LDB.LT.MAX( 1, N ) ) THEN
474:                INFO = -12
475:             ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
476:                INFO = -14
477:             END IF
478:          END IF
479:       END IF
480: *
481:       IF( INFO.NE.0 ) THEN
482:          CALL XERBLA( 'DPOSVXX', -INFO )
483:          RETURN
484:       END IF
485: *
486:       IF( EQUIL ) THEN
487: *
488: *     Compute row and column scalings to equilibrate the matrix A.
489: *
490:          CALL DPOEQUB( N, A, LDA, S, SCOND, AMAX, INFEQU )
491:          IF( INFEQU.EQ.0 ) THEN
492: *
493: *     Equilibrate the matrix.
494: *
495:             CALL DLAQSY( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED )
496:             RCEQU = LSAME( EQUED, 'Y' )
497:          END IF
498:       END IF
499: *
500: *     Scale the right-hand side.
501: *
502:       IF( RCEQU ) CALL DLASCL2( N, NRHS, S, B, LDB )
503: *
504:       IF( NOFACT .OR. EQUIL ) THEN
505: *
506: *        Compute the LU factorization of A.
507: *
508:          CALL DLACPY( UPLO, N, N, A, LDA, AF, LDAF )
509:          CALL DPOTRF( UPLO, N, AF, LDAF, INFO )
510: *
511: *        Return if INFO is non-zero.
512: *
513:          IF( INFO.NE.0 ) THEN
514: *
515: *           Pivot in column INFO is exactly 0
516: *           Compute the reciprocal pivot growth factor of the
517: *           leading rank-deficient INFO columns of A.
518: *
519:             RPVGRW = DLA_PORPVGRW( UPLO, INFO, A, LDA, AF, LDAF, WORK )
520:             RETURN
521:          ENDIF
522:       END IF
523: *
524: *     Compute the reciprocal growth factor RPVGRW.
525: *
526:       RPVGRW = DLA_PORPVGRW( UPLO, N, A, LDA, AF, LDAF, WORK )
527: *
528: *     Compute the solution matrix X.
529: *
530:       CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
531:       CALL DPOTRS( UPLO, N, NRHS, AF, LDAF, X, LDX, INFO )
532: *
533: *     Use iterative refinement to improve the computed solution and
534: *     compute error bounds and backward error estimates for it.
535: *
536:       CALL DPORFSX( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF,
537:      $     S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM,
538:      $     ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO )
539: 
540: *
541: *     Scale solutions.
542: *
543:       IF ( RCEQU ) THEN
544:          CALL DLASCL2 ( N, NRHS, S, X, LDX )
545:       END IF
546: *
547:       RETURN
548: *
549: *     End of DPOSVXX
550: *
551:       END
552: