3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DGERFSX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgerfsx.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgerfsx.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgerfsx.f">
21 * SUBROUTINE DGERFSX( TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV,
22 * R, C, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS,
23 * ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
26 * .. Scalar Arguments ..
27 * CHARACTER TRANS, EQUED
28 * INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
30 * DOUBLE PRECISION RCOND
32 * .. Array Arguments ..
33 * INTEGER IPIV( * ), IWORK( * )
34 * DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
35 * $ X( LDX , * ), WORK( * )
36 * DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
37 * $ ERR_BNDS_NORM( NRHS, * ),
38 * $ ERR_BNDS_COMP( NRHS, * )
47 *> DGERFSX improves the computed solution to a system of linear
48 *> equations and provides error bounds and backward error estimates
49 *> for the solution. In addition to normwise error bound, the code
50 *> provides maximum componentwise error bound if possible. See
51 *> comments for ERR_BNDS_NORM and ERR_BNDS_COMP for details of the
54 *> The original system of linear equations may have been equilibrated
55 *> before calling this routine, as described by arguments EQUED, R
56 *> and C below. In this case, the solution and error bounds returned
57 *> are for the original unequilibrated system.
64 *> Some optional parameters are bundled in the PARAMS array. These
65 *> settings determine how refinement is performed, but often the
66 *> defaults are acceptable. If the defaults are acceptable, users
67 *> can pass NPARAMS = 0 which prevents the source code from accessing
68 *> the PARAMS argument.
73 *> TRANS is CHARACTER*1
74 *> Specifies the form of the system of equations:
75 *> = 'N': A * X = B (No transpose)
76 *> = 'T': A**T * X = B (Transpose)
77 *> = 'C': A**H * X = B (Conjugate transpose = Transpose)
82 *> EQUED is CHARACTER*1
83 *> Specifies the form of equilibration that was done to A
84 *> before calling this routine. This is needed to compute
85 *> the solution and error bounds correctly.
86 *> = 'N': No equilibration
87 *> = 'R': Row equilibration, i.e., A has been premultiplied by
89 *> = 'C': Column equilibration, i.e., A has been postmultiplied
91 *> = 'B': Both row and column equilibration, i.e., A has been
92 *> replaced by diag(R) * A * diag(C).
93 *> The right hand side B has been changed accordingly.
99 *> The order of the matrix A. N >= 0.
105 *> The number of right hand sides, i.e., the number of columns
106 *> of the matrices B and X. NRHS >= 0.
111 *> A is DOUBLE PRECISION array, dimension (LDA,N)
112 *> The original N-by-N matrix A.
118 *> The leading dimension of the array A. LDA >= max(1,N).
123 *> AF is DOUBLE PRECISION array, dimension (LDAF,N)
124 *> The factors L and U from the factorization A = P*L*U
125 *> as computed by DGETRF.
131 *> The leading dimension of the array AF. LDAF >= max(1,N).
136 *> IPIV is INTEGER array, dimension (N)
137 *> The pivot indices from DGETRF; for 1<=i<=N, row i of the
138 *> matrix was interchanged with row IPIV(i).
143 *> R is DOUBLE PRECISION array, dimension (N)
144 *> The row scale factors for A. If EQUED = 'R' or 'B', A is
145 *> multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
147 *> If R is accessed, each element of R should be a power of the radix
148 *> to ensure a reliable solution and error estimates. Scaling by
149 *> powers of the radix does not cause rounding errors unless the
150 *> result underflows or overflows. Rounding errors during scaling
151 *> lead to refining with a matrix that is not equivalent to the
152 *> input matrix, producing error estimates that may not be
158 *> C is DOUBLE PRECISION array, dimension (N)
159 *> The column scale factors for A. If EQUED = 'C' or 'B', A is
160 *> multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
162 *> If C is accessed, each element of C should be a power of the radix
163 *> to ensure a reliable solution and error estimates. Scaling by
164 *> powers of the radix does not cause rounding errors unless the
165 *> result underflows or overflows. Rounding errors during scaling
166 *> lead to refining with a matrix that is not equivalent to the
167 *> input matrix, producing error estimates that may not be
173 *> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
174 *> The right hand side matrix B.
180 *> The leading dimension of the array B. LDB >= max(1,N).
185 *> X is DOUBLE PRECISION array, dimension (LDX,NRHS)
186 *> On entry, the solution matrix X, as computed by DGETRS.
187 *> On exit, the improved solution matrix X.
193 *> The leading dimension of the array X. LDX >= max(1,N).
198 *> RCOND is DOUBLE PRECISION
199 *> Reciprocal scaled condition number. This is an estimate of the
200 *> reciprocal Skeel condition number of the matrix A after
201 *> equilibration (if done). If this is less than the machine
202 *> precision (in particular, if it is zero), the matrix is singular
203 *> to working precision. Note that the error may still be small even
204 *> if this number is very small and the matrix appears ill-
210 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
211 *> Componentwise relative backward error. This is the
212 *> componentwise relative backward error of each solution vector X(j)
213 *> (i.e., the smallest relative change in any element of A or B that
214 *> makes X(j) an exact solution).
217 *> \param[in] N_ERR_BNDS
219 *> N_ERR_BNDS is INTEGER
220 *> Number of error bounds to return for each right hand side
221 *> and each type (normwise or componentwise). See ERR_BNDS_NORM and
222 *> ERR_BNDS_COMP below.
225 *> \param[out] ERR_BNDS_NORM
227 *> ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
228 *> For each right-hand side, this array contains information about
229 *> various error bounds and condition numbers corresponding to the
230 *> normwise relative error, which is defined as follows:
232 *> Normwise relative error in the ith solution vector:
233 *> max_j (abs(XTRUE(j,i) - X(j,i)))
234 *> ------------------------------
237 *> The array is indexed by the type of error information as described
238 *> below. There currently are up to three pieces of information
241 *> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
244 *> The second index in ERR_BNDS_NORM(:,err) contains the following
246 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
247 *> reciprocal condition number is less than the threshold
248 *> sqrt(n) * dlamch('Epsilon').
250 *> err = 2 "Guaranteed" error bound: The estimated forward error,
251 *> almost certainly within a factor of 10 of the true error
252 *> so long as the next entry is greater than the threshold
253 *> sqrt(n) * dlamch('Epsilon'). This error bound should only
254 *> be trusted if the previous boolean is true.
256 *> err = 3 Reciprocal condition number: Estimated normwise
257 *> reciprocal condition number. Compared with the threshold
258 *> sqrt(n) * dlamch('Epsilon') to determine if the error
259 *> estimate is "guaranteed". These reciprocal condition
260 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
261 *> appropriately scaled matrix Z.
262 *> Let Z = S*A, where S scales each row by a power of the
263 *> radix so all absolute row sums of Z are approximately 1.
265 *> See Lapack Working Note 165 for further details and extra
269 *> \param[out] ERR_BNDS_COMP
271 *> ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
272 *> For each right-hand side, this array contains information about
273 *> various error bounds and condition numbers corresponding to the
274 *> componentwise relative error, which is defined as follows:
276 *> Componentwise relative error in the ith solution vector:
277 *> abs(XTRUE(j,i) - X(j,i))
278 *> max_j ----------------------
281 *> The array is indexed by the right-hand side i (on which the
282 *> componentwise relative error depends), and the type of error
283 *> information as described below. There currently are up to three
284 *> pieces of information returned for each right-hand side. If
285 *> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
286 *> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most
287 *> the first (:,N_ERR_BNDS) entries are returned.
289 *> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
292 *> The second index in ERR_BNDS_COMP(:,err) contains the following
294 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
295 *> reciprocal condition number is less than the threshold
296 *> sqrt(n) * dlamch('Epsilon').
298 *> err = 2 "Guaranteed" error bound: The estimated forward error,
299 *> almost certainly within a factor of 10 of the true error
300 *> so long as the next entry is greater than the threshold
301 *> sqrt(n) * dlamch('Epsilon'). This error bound should only
302 *> be trusted if the previous boolean is true.
304 *> err = 3 Reciprocal condition number: Estimated componentwise
305 *> reciprocal condition number. Compared with the threshold
306 *> sqrt(n) * dlamch('Epsilon') to determine if the error
307 *> estimate is "guaranteed". These reciprocal condition
308 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
309 *> appropriately scaled matrix Z.
310 *> Let Z = S*(A*diag(x)), where x is the solution for the
311 *> current right-hand side and S scales each row of
312 *> A*diag(x) by a power of the radix so all absolute row
313 *> sums of Z are approximately 1.
315 *> See Lapack Working Note 165 for further details and extra
319 *> \param[in] NPARAMS
321 *> NPARAMS is INTEGER
322 *> Specifies the number of parameters set in PARAMS. If .LE. 0, the
323 *> PARAMS array is never referenced and default values are used.
326 *> \param[in,out] PARAMS
328 *> PARAMS is DOUBLE PRECISION array, dimension (NPARAMS)
329 *> Specifies algorithm parameters. If an entry is .LT. 0.0, then
330 *> that entry will be filled with default value used for that
331 *> parameter. Only positions up to NPARAMS are accessed; defaults
332 *> are used for higher-numbered parameters.
334 *> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
335 *> refinement or not.
337 *> = 0.0 : No refinement is performed, and no error bounds are
339 *> = 1.0 : Use the double-precision refinement algorithm,
340 *> possibly with doubled-single computations if the
341 *> compilation environment does not support DOUBLE
343 *> (other values are reserved for future use)
345 *> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
346 *> computations allowed for refinement.
348 *> Aggressive: Set to 100 to permit convergence using approximate
349 *> factorizations or factorizations other than LU. If
350 *> the factorization uses a technique other than
351 *> Gaussian elimination, the guarantees in
352 *> err_bnds_norm and err_bnds_comp may no longer be
355 *> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
356 *> will attempt to find a solution with small componentwise
357 *> relative error in the double-precision algorithm. Positive
358 *> is true, 0.0 is false.
359 *> Default: 1.0 (attempt componentwise convergence)
364 *> WORK is DOUBLE PRECISION array, dimension (4*N)
369 *> IWORK is INTEGER array, dimension (N)
375 *> = 0: Successful exit. The solution to every right-hand side is
377 *> < 0: If INFO = -i, the i-th argument had an illegal value
378 *> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
379 *> has been completed, but the factor U is exactly singular, so
380 *> the solution and error bounds could not be computed. RCOND = 0
382 *> = N+J: The solution corresponding to the Jth right-hand side is
383 *> not guaranteed. The solutions corresponding to other right-
384 *> hand sides K with K > J may not be guaranteed as well, but
385 *> only the first such right-hand side is reported. If a small
386 *> componentwise error is not requested (PARAMS(3) = 0.0) then
387 *> the Jth right-hand side is the first with a normwise error
388 *> bound that is not guaranteed (the smallest J such
389 *> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
390 *> the Jth right-hand side is the first with either a normwise or
391 *> componentwise error bound that is not guaranteed (the smallest
392 *> J such that either ERR_BNDS_NORM(J,1) = 0.0 or
393 *> ERR_BNDS_COMP(J,1) = 0.0). See the definition of
394 *> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
395 *> about all of the right-hand sides check ERR_BNDS_NORM or
402 *> \author Univ. of Tennessee
403 *> \author Univ. of California Berkeley
404 *> \author Univ. of Colorado Denver
407 *> \date November 2011
409 *> \ingroup doubleGEcomputational
411 * =====================================================================
412 SUBROUTINE DGERFSX( TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV,
413 $ R, C, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS,
414 $ ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
415 $ WORK, IWORK, INFO )
417 * -- LAPACK computational routine (version 3.4.0) --
418 * -- LAPACK is a software package provided by Univ. of Tennessee, --
419 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
422 * .. Scalar Arguments ..
423 CHARACTER TRANS, EQUED
424 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
426 DOUBLE PRECISION RCOND
428 * .. Array Arguments ..
429 INTEGER IPIV( * ), IWORK( * )
430 DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
431 $ X( LDX , * ), WORK( * )
432 DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
433 $ ERR_BNDS_NORM( NRHS, * ),
434 $ ERR_BNDS_COMP( NRHS, * )
437 * ==================================================================
440 DOUBLE PRECISION ZERO, ONE
441 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
442 DOUBLE PRECISION ITREF_DEFAULT, ITHRESH_DEFAULT
443 DOUBLE PRECISION COMPONENTWISE_DEFAULT, RTHRESH_DEFAULT
444 DOUBLE PRECISION DZTHRESH_DEFAULT
445 PARAMETER ( ITREF_DEFAULT = 1.0D+0 )
446 PARAMETER ( ITHRESH_DEFAULT = 10.0D+0 )
447 PARAMETER ( COMPONENTWISE_DEFAULT = 1.0D+0 )
448 PARAMETER ( RTHRESH_DEFAULT = 0.5D+0 )
449 PARAMETER ( DZTHRESH_DEFAULT = 0.25D+0 )
450 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
452 PARAMETER ( LA_LINRX_ITREF_I = 1,
453 $ LA_LINRX_ITHRESH_I = 2 )
454 PARAMETER ( LA_LINRX_CWISE_I = 3 )
455 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
457 PARAMETER ( LA_LINRX_TRUST_I = 1, LA_LINRX_ERR_I = 2 )
458 PARAMETER ( LA_LINRX_RCOND_I = 3 )
460 * .. Local Scalars ..
462 LOGICAL ROWEQU, COLEQU, NOTRAN
463 INTEGER J, TRANS_TYPE, PREC_TYPE, REF_TYPE
465 DOUBLE PRECISION ANORM, RCOND_TMP
466 DOUBLE PRECISION ILLRCOND_THRESH, ERR_LBND, CWISE_WRONG
469 DOUBLE PRECISION RTHRESH, UNSTABLE_THRESH
471 * .. External Subroutines ..
472 EXTERNAL XERBLA, DGECON, DLA_GERFSX_EXTENDED
474 * .. Intrinsic Functions ..
477 * .. External Functions ..
478 EXTERNAL LSAME, ILATRANS, ILAPREC
479 EXTERNAL DLAMCH, DLANGE, DLA_GERCOND
480 DOUBLE PRECISION DLAMCH, DLANGE, DLA_GERCOND
482 INTEGER ILATRANS, ILAPREC
484 * .. Executable Statements ..
486 * Check the input parameters.
489 TRANS_TYPE = ILATRANS( TRANS )
490 REF_TYPE = INT( ITREF_DEFAULT )
491 IF ( NPARAMS .GE. LA_LINRX_ITREF_I ) THEN
492 IF ( PARAMS( LA_LINRX_ITREF_I ) .LT. 0.0D+0 ) THEN
493 PARAMS( LA_LINRX_ITREF_I ) = ITREF_DEFAULT
495 REF_TYPE = PARAMS( LA_LINRX_ITREF_I )
499 * Set default parameters.
501 ILLRCOND_THRESH = DBLE( N ) * DLAMCH( 'Epsilon' )
502 ITHRESH = INT( ITHRESH_DEFAULT )
503 RTHRESH = RTHRESH_DEFAULT
504 UNSTABLE_THRESH = DZTHRESH_DEFAULT
505 IGNORE_CWISE = COMPONENTWISE_DEFAULT .EQ. 0.0D+0
507 IF ( NPARAMS.GE.LA_LINRX_ITHRESH_I ) THEN
508 IF ( PARAMS( LA_LINRX_ITHRESH_I ).LT.0.0D+0 ) THEN
509 PARAMS( LA_LINRX_ITHRESH_I ) = ITHRESH
511 ITHRESH = INT( PARAMS( LA_LINRX_ITHRESH_I ) )
514 IF ( NPARAMS.GE.LA_LINRX_CWISE_I ) THEN
515 IF ( PARAMS( LA_LINRX_CWISE_I ).LT.0.0D+0 ) THEN
516 IF ( IGNORE_CWISE ) THEN
517 PARAMS( LA_LINRX_CWISE_I ) = 0.0D+0
519 PARAMS( LA_LINRX_CWISE_I ) = 1.0D+0
522 IGNORE_CWISE = PARAMS( LA_LINRX_CWISE_I ) .EQ. 0.0D+0
525 IF ( REF_TYPE .EQ. 0 .OR. N_ERR_BNDS .EQ. 0 ) THEN
527 ELSE IF ( IGNORE_CWISE ) THEN
533 NOTRAN = LSAME( TRANS, 'N' )
534 ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
535 COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
537 * Test input parameters.
539 IF( TRANS_TYPE.EQ.-1 ) THEN
541 ELSE IF( .NOT.ROWEQU .AND. .NOT.COLEQU .AND.
542 $ .NOT.LSAME( EQUED, 'N' ) ) THEN
544 ELSE IF( N.LT.0 ) THEN
546 ELSE IF( NRHS.LT.0 ) THEN
548 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
550 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
552 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
554 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
558 CALL XERBLA( 'DGERFSX', -INFO )
562 * Quick return if possible.
564 IF( N.EQ.0 .OR. NRHS.EQ.0 ) THEN
568 IF ( N_ERR_BNDS .GE. 1 ) THEN
569 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I) = 1.0D+0
570 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 1.0D+0
572 IF ( N_ERR_BNDS .GE. 2 ) THEN
573 ERR_BNDS_NORM( J, LA_LINRX_ERR_I) = 0.0D+0
574 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 0.0D+0
576 IF ( N_ERR_BNDS .GE. 3 ) THEN
577 ERR_BNDS_NORM( J, LA_LINRX_RCOND_I) = 1.0D+0
578 ERR_BNDS_COMP( J, LA_LINRX_RCOND_I ) = 1.0D+0
584 * Default to failure.
589 IF ( N_ERR_BNDS .GE. 1 ) THEN
590 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 1.0D+0
591 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 1.0D+0
593 IF ( N_ERR_BNDS .GE. 2 ) THEN
594 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 1.0D+0
595 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 1.0D+0
597 IF ( N_ERR_BNDS .GE. 3 ) THEN
598 ERR_BNDS_NORM( J, LA_LINRX_RCOND_I ) = 0.0D+0
599 ERR_BNDS_COMP( J, LA_LINRX_RCOND_I ) = 0.0D+0
603 * Compute the norm of A and the reciprocal of the condition
611 ANORM = DLANGE( NORM, N, N, A, LDA, WORK )
612 CALL DGECON( NORM, N, AF, LDAF, ANORM, RCOND, WORK, IWORK, INFO )
614 * Perform refinement on each right-hand side
616 IF ( REF_TYPE .NE. 0 ) THEN
618 PREC_TYPE = ILAPREC( 'E' )
621 CALL DLA_GERFSX_EXTENDED( PREC_TYPE, TRANS_TYPE, N,
622 $ NRHS, A, LDA, AF, LDAF, IPIV, COLEQU, C, B,
623 $ LDB, X, LDX, BERR, N_NORMS, ERR_BNDS_NORM,
624 $ ERR_BNDS_COMP, WORK(N+1), WORK(1), WORK(2*N+1),
625 $ WORK(1), RCOND, ITHRESH, RTHRESH, UNSTABLE_THRESH,
626 $ IGNORE_CWISE, INFO )
628 CALL DLA_GERFSX_EXTENDED( PREC_TYPE, TRANS_TYPE, N,
629 $ NRHS, A, LDA, AF, LDAF, IPIV, ROWEQU, R, B,
630 $ LDB, X, LDX, BERR, N_NORMS, ERR_BNDS_NORM,
631 $ ERR_BNDS_COMP, WORK(N+1), WORK(1), WORK(2*N+1),
632 $ WORK(1), RCOND, ITHRESH, RTHRESH, UNSTABLE_THRESH,
633 $ IGNORE_CWISE, INFO )
637 ERR_LBND = MAX( 10.0D+0, SQRT( DBLE( N ) ) ) * DLAMCH( 'Epsilon' )
638 IF ( N_ERR_BNDS .GE. 1 .AND. N_NORMS .GE. 1 ) THEN
640 * Compute scaled normwise condition number cond(A*C).
642 IF ( COLEQU .AND. NOTRAN ) THEN
643 RCOND_TMP = DLA_GERCOND( TRANS, N, A, LDA, AF, LDAF, IPIV,
644 $ -1, C, INFO, WORK, IWORK )
645 ELSE IF ( ROWEQU .AND. .NOT. NOTRAN ) THEN
646 RCOND_TMP = DLA_GERCOND( TRANS, N, A, LDA, AF, LDAF, IPIV,
647 $ -1, R, INFO, WORK, IWORK )
649 RCOND_TMP = DLA_GERCOND( TRANS, N, A, LDA, AF, LDAF, IPIV,
650 $ 0, R, INFO, WORK, IWORK )
654 * Cap the error at 1.0.
656 IF ( N_ERR_BNDS .GE. LA_LINRX_ERR_I
657 $ .AND. ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) .GT. 1.0D+0 )
658 $ ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 1.0D+0
660 * Threshold the error (see LAWN).
662 IF ( RCOND_TMP .LT. ILLRCOND_THRESH ) THEN
663 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 1.0D+0
664 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 0.0D+0
665 IF ( INFO .LE. N ) INFO = N + J
666 ELSE IF ( ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) .LT. ERR_LBND )
668 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = ERR_LBND
669 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 1.0D+0
672 * Save the condition number.
674 IF ( N_ERR_BNDS .GE. LA_LINRX_RCOND_I ) THEN
675 ERR_BNDS_NORM( J, LA_LINRX_RCOND_I ) = RCOND_TMP
680 IF ( N_ERR_BNDS .GE. 1 .AND. N_NORMS .GE. 2 ) THEN
682 * Compute componentwise condition number cond(A*diag(Y(:,J))) for
683 * each right-hand side using the current solution as an estimate of
684 * the true solution. If the componentwise error estimate is too
685 * large, then the solution is a lousy estimate of truth and the
686 * estimated RCOND may be too optimistic. To avoid misleading users,
687 * the inverse condition number is set to 0.0 when the estimated
688 * cwise error is at least CWISE_WRONG.
690 CWISE_WRONG = SQRT( DLAMCH( 'Epsilon' ) )
692 IF ( ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) .LT. CWISE_WRONG )
694 RCOND_TMP = DLA_GERCOND( TRANS, N, A, LDA, AF, LDAF,
695 $ IPIV, 1, X(1,J), INFO, WORK, IWORK )
700 * Cap the error at 1.0.
702 IF ( N_ERR_BNDS .GE. LA_LINRX_ERR_I
703 $ .AND. ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) .GT. 1.0D+0 )
704 $ ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 1.0D+0
706 * Threshold the error (see LAWN).
708 IF ( RCOND_TMP .LT. ILLRCOND_THRESH ) THEN
709 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 1.0D+0
710 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 0.0D+0
711 IF ( PARAMS( LA_LINRX_CWISE_I ) .EQ. 1.0D+0
712 $ .AND. INFO.LT.N + J ) INFO = N + J
713 ELSE IF ( ERR_BNDS_COMP( J, LA_LINRX_ERR_I )
714 $ .LT. ERR_LBND ) THEN
715 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = ERR_LBND
716 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 1.0D+0
719 * Save the condition number.
721 IF ( N_ERR_BNDS .GE. LA_LINRX_RCOND_I ) THEN
722 ERR_BNDS_COMP( J, LA_LINRX_RCOND_I ) = RCOND_TMP