1 *> \brief <b> ZHESVXX computes the solution to system of linear equations A * X = B for HE matrices</b>
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download ZHESVXX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhesvxx.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhesvxx.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhesvxx.f">
21 * SUBROUTINE ZHESVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV,
22 * EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
23 * N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
24 * NPARAMS, PARAMS, WORK, RWORK, INFO )
26 * .. Scalar Arguments ..
27 * CHARACTER EQUED, FACT, UPLO
28 * INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
30 * DOUBLE PRECISION RCOND, RPVGRW
32 * .. Array Arguments ..
34 * COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
35 * $ WORK( * ), X( LDX, * )
36 * DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
37 * $ ERR_BNDS_NORM( NRHS, * ),
38 * $ ERR_BNDS_COMP( NRHS, * )
47 *> ZHESVXX uses the diagonal pivoting factorization to compute the
48 *> solution to a complex*16 system of linear equations A * X = B, where
49 *> A is an N-by-N symmetric matrix and X and B are N-by-NRHS
52 *> If requested, both normwise and maximum componentwise error bounds
53 *> are returned. ZHESVXX will return a solution with a tiny
54 *> guaranteed error (O(eps) where eps is the working machine
55 *> precision) unless the matrix is very ill-conditioned, in which
56 *> case a warning is returned. Relevant condition numbers also are
57 *> calculated and returned.
59 *> ZHESVXX accepts user-provided factorizations and equilibration
60 *> factors; see the definitions of the FACT and EQUED options.
61 *> Solving with refinement and using a factorization from a previous
62 *> ZHESVXX call will also produce a solution with either O(eps)
63 *> errors or warnings, but we cannot make that claim for general
64 *> user-provided factorizations and equilibration factors if they
65 *> differ from what ZHESVXX would itself produce.
73 *> The following steps are performed:
75 *> 1. If FACT = 'E', double precision scaling factors are computed to equilibrate
78 *> diag(S)*A*diag(S) *inv(diag(S))*X = diag(S)*B
80 *> Whether or not the system will be equilibrated depends on the
81 *> scaling of the matrix A, but if equilibration is used, A is
82 *> overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
84 *> 2. If FACT = 'N' or 'E', the LU decomposition is used to factor
85 *> the matrix A (after equilibration if FACT = 'E') as
87 *> A = U * D * U**T, if UPLO = 'U', or
88 *> A = L * D * L**T, if UPLO = 'L',
90 *> where U (or L) is a product of permutation and unit upper (lower)
91 *> triangular matrices, and D is symmetric and block diagonal with
92 *> 1-by-1 and 2-by-2 diagonal blocks.
94 *> 3. If some D(i,i)=0, so that D is exactly singular, then the
95 *> routine returns with INFO = i. Otherwise, the factored form of A
96 *> is used to estimate the condition number of the matrix A (see
97 *> argument RCOND). If the reciprocal of the condition number is
98 *> less than machine precision, the routine still goes on to solve
99 *> for X and compute error bounds as described below.
101 *> 4. The system of equations is solved for X using the factored form
104 *> 5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
105 *> the routine will use iterative refinement to try to get a small
106 *> error and error bounds. Refinement calculates the residual to at
107 *> least twice the working precision.
109 *> 6. If equilibration was used, the matrix X is premultiplied by
110 *> diag(R) so that it solves the original system before
118 *> Some optional parameters are bundled in the PARAMS array. These
119 *> settings determine how refinement is performed, but often the
120 *> defaults are acceptable. If the defaults are acceptable, users
121 *> can pass NPARAMS = 0 which prevents the source code from accessing
122 *> the PARAMS argument.
127 *> FACT is CHARACTER*1
128 *> Specifies whether or not the factored form of the matrix A is
129 *> supplied on entry, and if not, whether the matrix A should be
130 *> equilibrated before it is factored.
131 *> = 'F': On entry, AF and IPIV contain the factored form of A.
132 *> If EQUED is not 'N', the matrix A has been
133 *> equilibrated with scaling factors given by S.
134 *> A, AF, and IPIV are not modified.
135 *> = 'N': The matrix A will be copied to AF and factored.
136 *> = 'E': The matrix A will be equilibrated if necessary, then
137 *> copied to AF and factored.
142 *> UPLO is CHARACTER*1
143 *> = 'U': Upper triangle of A is stored;
144 *> = 'L': Lower triangle of A is stored.
150 *> The number of linear equations, i.e., the order of the
157 *> The number of right hand sides, i.e., the number of columns
158 *> of the matrices B and X. NRHS >= 0.
163 *> A is COMPLEX*16 array, dimension (LDA,N)
164 *> The symmetric matrix A. If UPLO = 'U', the leading N-by-N
165 *> upper triangular part of A contains the upper triangular
166 *> part of the matrix A, and the strictly lower triangular
167 *> part of A is not referenced. If UPLO = 'L', the leading
168 *> N-by-N lower triangular part of A contains the lower
169 *> triangular part of the matrix A, and the strictly upper
170 *> triangular part of A is not referenced.
172 *> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
173 *> diag(S)*A*diag(S).
179 *> The leading dimension of the array A. LDA >= max(1,N).
184 *> AF is COMPLEX*16 array, dimension (LDAF,N)
185 *> If FACT = 'F', then AF is an input argument and on entry
186 *> contains the block diagonal matrix D and the multipliers
187 *> used to obtain the factor U or L from the factorization A =
188 *> U*D*U**T or A = L*D*L**T as computed by DSYTRF.
190 *> If FACT = 'N', then AF is an output argument and on exit
191 *> returns the block diagonal matrix D and the multipliers
192 *> used to obtain the factor U or L from the factorization A =
193 *> U*D*U**T or A = L*D*L**T.
199 *> The leading dimension of the array AF. LDAF >= max(1,N).
202 *> \param[in,out] IPIV
204 *> IPIV is INTEGER array, dimension (N)
205 *> If FACT = 'F', then IPIV is an input argument and on entry
206 *> contains details of the interchanges and the block
207 *> structure of D, as determined by ZHETRF. If IPIV(k) > 0,
208 *> then rows and columns k and IPIV(k) were interchanged and
209 *> D(k,k) is a 1-by-1 diagonal block. If UPLO = 'U' and
210 *> IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and
211 *> -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2
212 *> diagonal block. If UPLO = 'L' and IPIV(k) = IPIV(k+1) < 0,
213 *> then rows and columns k+1 and -IPIV(k) were interchanged
214 *> and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
216 *> If FACT = 'N', then IPIV is an output argument and on exit
217 *> contains details of the interchanges and the block
218 *> structure of D, as determined by ZHETRF.
221 *> \param[in,out] EQUED
223 *> EQUED is CHARACTER*1
224 *> Specifies the form of equilibration that was done.
225 *> = 'N': No equilibration (always true if FACT = 'N').
226 *> = 'Y': Both row and column equilibration, i.e., A has been
227 *> replaced by diag(S) * A * diag(S).
228 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an
234 *> S is DOUBLE PRECISION array, dimension (N)
235 *> The scale factors for A. If EQUED = 'Y', A is multiplied on
236 *> the left and right by diag(S). S is an input argument if FACT =
237 *> 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED
238 *> = 'Y', each element of S must be positive. If S is output, each
239 *> element of S is a power of the radix. If S is input, each element
240 *> of S should be a power of the radix to ensure a reliable solution
241 *> and error estimates. Scaling by powers of the radix does not cause
242 *> rounding errors unless the result underflows or overflows.
243 *> Rounding errors during scaling lead to refining with a matrix that
244 *> is not equivalent to the input matrix, producing error estimates
245 *> that may not be reliable.
250 *> B is COMPLEX*16 array, dimension (LDB,NRHS)
251 *> On entry, the N-by-NRHS right hand side matrix B.
253 *> if EQUED = 'N', B is not modified;
254 *> if EQUED = 'Y', B is overwritten by diag(S)*B;
260 *> The leading dimension of the array B. LDB >= max(1,N).
265 *> X is COMPLEX*16 array, dimension (LDX,NRHS)
266 *> If INFO = 0, the N-by-NRHS solution matrix X to the original
267 *> system of equations. Note that A and B are modified on exit if
268 *> EQUED .ne. 'N', and the solution to the equilibrated system is
275 *> The leading dimension of the array X. LDX >= max(1,N).
280 *> RCOND is DOUBLE PRECISION
281 *> Reciprocal scaled condition number. This is an estimate of the
282 *> reciprocal Skeel condition number of the matrix A after
283 *> equilibration (if done). If this is less than the machine
284 *> precision (in particular, if it is zero), the matrix is singular
285 *> to working precision. Note that the error may still be small even
286 *> if this number is very small and the matrix appears ill-
290 *> \param[out] RPVGRW
292 *> RPVGRW is DOUBLE PRECISION
293 *> Reciprocal pivot growth. On exit, this contains the reciprocal
294 *> pivot growth factor norm(A)/norm(U). The "max absolute element"
295 *> norm is used. If this is much less than 1, then the stability of
296 *> the LU factorization of the (equilibrated) matrix A could be poor.
297 *> This also means that the solution X, estimated condition numbers,
298 *> and error bounds could be unreliable. If factorization fails with
299 *> 0<INFO<=N, then this contains the reciprocal pivot growth factor
300 *> for the leading INFO columns of A.
305 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
306 *> Componentwise relative backward error. This is the
307 *> componentwise relative backward error of each solution vector X(j)
308 *> (i.e., the smallest relative change in any element of A or B that
309 *> makes X(j) an exact solution).
312 *> \param[in] N_ERR_BNDS
314 *> N_ERR_BNDS is INTEGER
315 *> Number of error bounds to return for each right hand side
316 *> and each type (normwise or componentwise). See ERR_BNDS_NORM and
317 *> ERR_BNDS_COMP below.
320 *> \param[out] ERR_BNDS_NORM
322 *> ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
323 *> For each right-hand side, this array contains information about
324 *> various error bounds and condition numbers corresponding to the
325 *> normwise relative error, which is defined as follows:
327 *> Normwise relative error in the ith solution vector:
328 *> max_j (abs(XTRUE(j,i) - X(j,i)))
329 *> ------------------------------
332 *> The array is indexed by the type of error information as described
333 *> below. There currently are up to three pieces of information
336 *> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
339 *> The second index in ERR_BNDS_NORM(:,err) contains the following
341 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
342 *> reciprocal condition number is less than the threshold
343 *> sqrt(n) * dlamch('Epsilon').
345 *> err = 2 "Guaranteed" error bound: The estimated forward error,
346 *> almost certainly within a factor of 10 of the true error
347 *> so long as the next entry is greater than the threshold
348 *> sqrt(n) * dlamch('Epsilon'). This error bound should only
349 *> be trusted if the previous boolean is true.
351 *> err = 3 Reciprocal condition number: Estimated normwise
352 *> reciprocal condition number. Compared with the threshold
353 *> sqrt(n) * dlamch('Epsilon') to determine if the error
354 *> estimate is "guaranteed". These reciprocal condition
355 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
356 *> appropriately scaled matrix Z.
357 *> Let Z = S*A, where S scales each row by a power of the
358 *> radix so all absolute row sums of Z are approximately 1.
360 *> See Lapack Working Note 165 for further details and extra
364 *> \param[out] ERR_BNDS_COMP
366 *> ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
367 *> For each right-hand side, this array contains information about
368 *> various error bounds and condition numbers corresponding to the
369 *> componentwise relative error, which is defined as follows:
371 *> Componentwise relative error in the ith solution vector:
372 *> abs(XTRUE(j,i) - X(j,i))
373 *> max_j ----------------------
376 *> The array is indexed by the right-hand side i (on which the
377 *> componentwise relative error depends), and the type of error
378 *> information as described below. There currently are up to three
379 *> pieces of information returned for each right-hand side. If
380 *> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
381 *> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most
382 *> the first (:,N_ERR_BNDS) entries are returned.
384 *> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
387 *> The second index in ERR_BNDS_COMP(:,err) contains the following
389 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
390 *> reciprocal condition number is less than the threshold
391 *> sqrt(n) * dlamch('Epsilon').
393 *> err = 2 "Guaranteed" error bound: The estimated forward error,
394 *> almost certainly within a factor of 10 of the true error
395 *> so long as the next entry is greater than the threshold
396 *> sqrt(n) * dlamch('Epsilon'). This error bound should only
397 *> be trusted if the previous boolean is true.
399 *> err = 3 Reciprocal condition number: Estimated componentwise
400 *> reciprocal condition number. Compared with the threshold
401 *> sqrt(n) * dlamch('Epsilon') to determine if the error
402 *> estimate is "guaranteed". These reciprocal condition
403 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
404 *> appropriately scaled matrix Z.
405 *> Let Z = S*(A*diag(x)), where x is the solution for the
406 *> current right-hand side and S scales each row of
407 *> A*diag(x) by a power of the radix so all absolute row
408 *> sums of Z are approximately 1.
410 *> See Lapack Working Note 165 for further details and extra
414 *> \param[in] NPARAMS
416 *> NPARAMS is INTEGER
417 *> Specifies the number of parameters set in PARAMS. If .LE. 0, the
418 *> PARAMS array is never referenced and default values are used.
421 *> \param[in,out] PARAMS
423 *> PARAMS is DOUBLE PRECISION array, dimension NPARAMS
424 *> Specifies algorithm parameters. If an entry is .LT. 0.0, then
425 *> that entry will be filled with default value used for that
426 *> parameter. Only positions up to NPARAMS are accessed; defaults
427 *> are used for higher-numbered parameters.
429 *> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
430 *> refinement or not.
432 *> = 0.0 : No refinement is performed, and no error bounds are
434 *> = 1.0 : Use the extra-precise refinement algorithm.
435 *> (other values are reserved for future use)
437 *> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
438 *> computations allowed for refinement.
440 *> Aggressive: Set to 100 to permit convergence using approximate
441 *> factorizations or factorizations other than LU. If
442 *> the factorization uses a technique other than
443 *> Gaussian elimination, the guarantees in
444 *> err_bnds_norm and err_bnds_comp may no longer be
447 *> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
448 *> will attempt to find a solution with small componentwise
449 *> relative error in the double-precision algorithm. Positive
450 *> is true, 0.0 is false.
451 *> Default: 1.0 (attempt componentwise convergence)
456 *> WORK is COMPLEX*16 array, dimension (5*N)
461 *> RWORK is DOUBLE PRECISION array, dimension (2*N)
467 *> = 0: Successful exit. The solution to every right-hand side is
469 *> < 0: If INFO = -i, the i-th argument had an illegal value
470 *> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
471 *> has been completed, but the factor U is exactly singular, so
472 *> the solution and error bounds could not be computed. RCOND = 0
474 *> = N+J: The solution corresponding to the Jth right-hand side is
475 *> not guaranteed. The solutions corresponding to other right-
476 *> hand sides K with K > J may not be guaranteed as well, but
477 *> only the first such right-hand side is reported. If a small
478 *> componentwise error is not requested (PARAMS(3) = 0.0) then
479 *> the Jth right-hand side is the first with a normwise error
480 *> bound that is not guaranteed (the smallest J such
481 *> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
482 *> the Jth right-hand side is the first with either a normwise or
483 *> componentwise error bound that is not guaranteed (the smallest
484 *> J such that either ERR_BNDS_NORM(J,1) = 0.0 or
485 *> ERR_BNDS_COMP(J,1) = 0.0). See the definition of
486 *> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
487 *> about all of the right-hand sides check ERR_BNDS_NORM or
494 *> \author Univ. of Tennessee
495 *> \author Univ. of California Berkeley
496 *> \author Univ. of Colorado Denver
501 *> \ingroup complex16HEsolve
503 * =====================================================================
504 SUBROUTINE ZHESVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV,
505 $ EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
506 $ N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
507 $ NPARAMS, PARAMS, WORK, RWORK, INFO )
509 * -- LAPACK driver routine (version 3.4.1) --
510 * -- LAPACK is a software package provided by Univ. of Tennessee, --
511 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
514 * .. Scalar Arguments ..
515 CHARACTER EQUED, FACT, UPLO
516 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
518 DOUBLE PRECISION RCOND, RPVGRW
520 * .. Array Arguments ..
522 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
523 $ WORK( * ), X( LDX, * )
524 DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
525 $ ERR_BNDS_NORM( NRHS, * ),
526 $ ERR_BNDS_COMP( NRHS, * )
529 * ==================================================================
532 DOUBLE PRECISION ZERO, ONE
533 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
534 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
535 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
536 INTEGER CMP_ERR_I, PIV_GROWTH_I
537 PARAMETER ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
539 PARAMETER ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
540 PARAMETER ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
543 * .. Local Scalars ..
544 LOGICAL EQUIL, NOFACT, RCEQU
546 DOUBLE PRECISION AMAX, BIGNUM, SMIN, SMAX, SCOND, SMLNUM
548 * .. External Functions ..
549 EXTERNAL LSAME, DLAMCH, ZLA_HERPVGRW
551 DOUBLE PRECISION DLAMCH, ZLA_HERPVGRW
553 * .. External Subroutines ..
554 EXTERNAL ZHEEQUB, ZHETRF, ZHETRS, ZLACPY,
555 $ ZLAQHE, XERBLA, ZLASCL2, ZHERFSX
557 * .. Intrinsic Functions ..
560 * .. Executable Statements ..
563 NOFACT = LSAME( FACT, 'N' )
564 EQUIL = LSAME( FACT, 'E' )
565 SMLNUM = DLAMCH( 'Safe minimum' )
566 BIGNUM = ONE / SMLNUM
567 IF( NOFACT .OR. EQUIL ) THEN
571 RCEQU = LSAME( EQUED, 'Y' )
574 * Default is failure. If an input parameter is wrong or
575 * factorization fails, make everything look horrible. Only the
576 * pivot growth is set here, the rest is initialized in ZHERFSX.
580 * Test the input parameters. PARAMS is not tested until ZHERFSX.
582 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
583 $ LSAME( FACT, 'F' ) ) THEN
585 ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND.
586 $ .NOT.LSAME( UPLO, 'L' ) ) THEN
588 ELSE IF( N.LT.0 ) THEN
590 ELSE IF( NRHS.LT.0 ) THEN
592 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
594 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
596 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
597 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
604 SMIN = MIN( SMIN, S( J ) )
605 SMAX = MAX( SMAX, S( J ) )
607 IF( SMIN.LE.ZERO ) THEN
609 ELSE IF( N.GT.0 ) THEN
610 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
616 IF( LDB.LT.MAX( 1, N ) ) THEN
618 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
625 CALL XERBLA( 'ZHESVXX', -INFO )
631 * Compute row and column scalings to equilibrate the matrix A.
633 CALL ZHEEQUB( UPLO, N, A, LDA, S, SCOND, AMAX, WORK, INFEQU )
634 IF( INFEQU.EQ.0 ) THEN
636 * Equilibrate the matrix.
638 CALL ZLAQHE( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED )
639 RCEQU = LSAME( EQUED, 'Y' )
643 * Scale the right-hand side.
645 IF( RCEQU ) CALL ZLASCL2( N, NRHS, S, B, LDB )
647 IF( NOFACT .OR. EQUIL ) THEN
649 * Compute the LDL^T or UDU^T factorization of A.
651 CALL ZLACPY( UPLO, N, N, A, LDA, AF, LDAF )
652 CALL ZHETRF( UPLO, N, AF, LDAF, IPIV, WORK, 5*MAX(1,N), INFO )
654 * Return if INFO is non-zero.
658 * Pivot in column INFO is exactly 0
659 * Compute the reciprocal pivot growth factor of the
660 * leading rank-deficient INFO columns of A.
663 $ RPVGRW = ZLA_HERPVGRW( UPLO, N, INFO, A, LDA, AF, LDAF,
669 * Compute the reciprocal pivot growth factor RPVGRW.
672 $ RPVGRW = ZLA_HERPVGRW( UPLO, N, INFO, A, LDA, AF, LDAF, IPIV,
675 * Compute the solution matrix X.
677 CALL ZLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
678 CALL ZHETRS( UPLO, N, NRHS, AF, LDAF, IPIV, X, LDX, INFO )
680 * Use iterative refinement to improve the computed solution and
681 * compute error bounds and backward error estimates for it.
683 CALL ZHERFSX( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV,
684 $ S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM,
685 $ ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, RWORK, INFO )
690 CALL ZLASCL2 ( N, NRHS, S, X, LDX )