LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgerfsx.f
Go to the documentation of this file.
1*> \brief \b ZGERFSX
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGERFSX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgerfsx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgerfsx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgerfsx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGERFSX( TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV,
20* R, C, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS,
21* ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
22* WORK, RWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER TRANS, EQUED
26* INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
27* $ N_ERR_BNDS
28* DOUBLE PRECISION RCOND
29* ..
30* .. Array Arguments ..
31* INTEGER IPIV( * )
32* COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
33* $ X( LDX , * ), WORK( * )
34* DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
35* $ ERR_BNDS_NORM( NRHS, * ),
36* $ ERR_BNDS_COMP( NRHS, * ), RWORK( * )
37* ..
38*
39*
40*> \par Purpose:
41* =============
42*>
43*> \verbatim
44*>
45*> ZGERFSX improves the computed solution to a system of linear
46*> equations and provides error bounds and backward error estimates
47*> for the solution. In addition to normwise error bound, the code
48*> provides maximum componentwise error bound if possible. See
49*> comments for ERR_BNDS_NORM and ERR_BNDS_COMP for details of the
50*> error bounds.
51*>
52*> The original system of linear equations may have been equilibrated
53*> before calling this routine, as described by arguments EQUED, R
54*> and C below. In this case, the solution and error bounds returned
55*> are for the original unequilibrated system.
56*> \endverbatim
57*
58* Arguments:
59* ==========
60*
61*> \verbatim
62*> Some optional parameters are bundled in the PARAMS array. These
63*> settings determine how refinement is performed, but often the
64*> defaults are acceptable. If the defaults are acceptable, users
65*> can pass NPARAMS = 0 which prevents the source code from accessing
66*> the PARAMS argument.
67*> \endverbatim
68*>
69*> \param[in] TRANS
70*> \verbatim
71*> TRANS is CHARACTER*1
72*> Specifies the form of the system of equations:
73*> = 'N': A * X = B (No transpose)
74*> = 'T': A**T * X = B (Transpose)
75*> = 'C': A**H * X = B (Conjugate transpose)
76*> \endverbatim
77*>
78*> \param[in] EQUED
79*> \verbatim
80*> EQUED is CHARACTER*1
81*> Specifies the form of equilibration that was done to A
82*> before calling this routine. This is needed to compute
83*> the solution and error bounds correctly.
84*> = 'N': No equilibration
85*> = 'R': Row equilibration, i.e., A has been premultiplied by
86*> diag(R).
87*> = 'C': Column equilibration, i.e., A has been postmultiplied
88*> by diag(C).
89*> = 'B': Both row and column equilibration, i.e., A has been
90*> replaced by diag(R) * A * diag(C).
91*> The right hand side B has been changed accordingly.
92*> \endverbatim
93*>
94*> \param[in] N
95*> \verbatim
96*> N is INTEGER
97*> The order of the matrix A. N >= 0.
98*> \endverbatim
99*>
100*> \param[in] NRHS
101*> \verbatim
102*> NRHS is INTEGER
103*> The number of right hand sides, i.e., the number of columns
104*> of the matrices B and X. NRHS >= 0.
105*> \endverbatim
106*>
107*> \param[in] A
108*> \verbatim
109*> A is COMPLEX*16 array, dimension (LDA,N)
110*> The original N-by-N matrix A.
111*> \endverbatim
112*>
113*> \param[in] LDA
114*> \verbatim
115*> LDA is INTEGER
116*> The leading dimension of the array A. LDA >= max(1,N).
117*> \endverbatim
118*>
119*> \param[in] AF
120*> \verbatim
121*> AF is COMPLEX*16 array, dimension (LDAF,N)
122*> The factors L and U from the factorization A = P*L*U
123*> as computed by ZGETRF.
124*> \endverbatim
125*>
126*> \param[in] LDAF
127*> \verbatim
128*> LDAF is INTEGER
129*> The leading dimension of the array AF. LDAF >= max(1,N).
130*> \endverbatim
131*>
132*> \param[in] IPIV
133*> \verbatim
134*> IPIV is INTEGER array, dimension (N)
135*> The pivot indices from ZGETRF; for 1<=i<=N, row i of the
136*> matrix was interchanged with row IPIV(i).
137*> \endverbatim
138*>
139*> \param[in] R
140*> \verbatim
141*> R is DOUBLE PRECISION array, dimension (N)
142*> The row scale factors for A. If EQUED = 'R' or 'B', A is
143*> multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
144*> is not accessed.
145*> If R is accessed, each element of R should be a power of the radix
146*> to ensure a reliable solution and error estimates. Scaling by
147*> powers of the radix does not cause rounding errors unless the
148*> result underflows or overflows. Rounding errors during scaling
149*> lead to refining with a matrix that is not equivalent to the
150*> input matrix, producing error estimates that may not be
151*> reliable.
152*> \endverbatim
153*>
154*> \param[in] C
155*> \verbatim
156*> C is DOUBLE PRECISION array, dimension (N)
157*> The column scale factors for A. If EQUED = 'C' or 'B', A is
158*> multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
159*> is not accessed.
160*> If C is accessed, each element of C should be a power of the radix
161*> to ensure a reliable solution and error estimates. Scaling by
162*> powers of the radix does not cause rounding errors unless the
163*> result underflows or overflows. Rounding errors during scaling
164*> lead to refining with a matrix that is not equivalent to the
165*> input matrix, producing error estimates that may not be
166*> reliable.
167*> \endverbatim
168*>
169*> \param[in] B
170*> \verbatim
171*> B is COMPLEX*16 array, dimension (LDB,NRHS)
172*> The right hand side matrix B.
173*> \endverbatim
174*>
175*> \param[in] LDB
176*> \verbatim
177*> LDB is INTEGER
178*> The leading dimension of the array B. LDB >= max(1,N).
179*> \endverbatim
180*>
181*> \param[in,out] X
182*> \verbatim
183*> X is COMPLEX*16 array, dimension (LDX,NRHS)
184*> On entry, the solution matrix X, as computed by ZGETRS.
185*> On exit, the improved solution matrix X.
186*> \endverbatim
187*>
188*> \param[in] LDX
189*> \verbatim
190*> LDX is INTEGER
191*> The leading dimension of the array X. LDX >= max(1,N).
192*> \endverbatim
193*>
194*> \param[out] RCOND
195*> \verbatim
196*> RCOND is DOUBLE PRECISION
197*> Reciprocal scaled condition number. This is an estimate of the
198*> reciprocal Skeel condition number of the matrix A after
199*> equilibration (if done). If this is less than the machine
200*> precision (in particular, if it is zero), the matrix is singular
201*> to working precision. Note that the error may still be small even
202*> if this number is very small and the matrix appears ill-
203*> conditioned.
204*> \endverbatim
205*>
206*> \param[out] BERR
207*> \verbatim
208*> BERR is DOUBLE PRECISION array, dimension (NRHS)
209*> Componentwise relative backward error. This is the
210*> componentwise relative backward error of each solution vector X(j)
211*> (i.e., the smallest relative change in any element of A or B that
212*> makes X(j) an exact solution).
213*> \endverbatim
214*>
215*> \param[in] N_ERR_BNDS
216*> \verbatim
217*> N_ERR_BNDS is INTEGER
218*> Number of error bounds to return for each right hand side
219*> and each type (normwise or componentwise). See ERR_BNDS_NORM and
220*> ERR_BNDS_COMP below.
221*> \endverbatim
222*>
223*> \param[out] ERR_BNDS_NORM
224*> \verbatim
225*> ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
226*> For each right-hand side, this array contains information about
227*> various error bounds and condition numbers corresponding to the
228*> normwise relative error, which is defined as follows:
229*>
230*> Normwise relative error in the ith solution vector:
231*> max_j (abs(XTRUE(j,i) - X(j,i)))
232*> ------------------------------
233*> max_j abs(X(j,i))
234*>
235*> The array is indexed by the type of error information as described
236*> below. There currently are up to three pieces of information
237*> returned.
238*>
239*> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
240*> right-hand side.
241*>
242*> The second index in ERR_BNDS_NORM(:,err) contains the following
243*> three fields:
244*> err = 1 "Trust/don't trust" boolean. Trust the answer if the
245*> reciprocal condition number is less than the threshold
246*> sqrt(n) * dlamch('Epsilon').
247*>
248*> err = 2 "Guaranteed" error bound: The estimated forward error,
249*> almost certainly within a factor of 10 of the true error
250*> so long as the next entry is greater than the threshold
251*> sqrt(n) * dlamch('Epsilon'). This error bound should only
252*> be trusted if the previous boolean is true.
253*>
254*> err = 3 Reciprocal condition number: Estimated normwise
255*> reciprocal condition number. Compared with the threshold
256*> sqrt(n) * dlamch('Epsilon') to determine if the error
257*> estimate is "guaranteed". These reciprocal condition
258*> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
259*> appropriately scaled matrix Z.
260*> Let Z = S*A, where S scales each row by a power of the
261*> radix so all absolute row sums of Z are approximately 1.
262*>
263*> See Lapack Working Note 165 for further details and extra
264*> cautions.
265*> \endverbatim
266*>
267*> \param[out] ERR_BNDS_COMP
268*> \verbatim
269*> ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
270*> For each right-hand side, this array contains information about
271*> various error bounds and condition numbers corresponding to the
272*> componentwise relative error, which is defined as follows:
273*>
274*> Componentwise relative error in the ith solution vector:
275*> abs(XTRUE(j,i) - X(j,i))
276*> max_j ----------------------
277*> abs(X(j,i))
278*>
279*> The array is indexed by the right-hand side i (on which the
280*> componentwise relative error depends), and the type of error
281*> information as described below. There currently are up to three
282*> pieces of information returned for each right-hand side. If
283*> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
284*> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS < 3, then at most
285*> the first (:,N_ERR_BNDS) entries are returned.
286*>
287*> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
288*> right-hand side.
289*>
290*> The second index in ERR_BNDS_COMP(:,err) contains the following
291*> three fields:
292*> err = 1 "Trust/don't trust" boolean. Trust the answer if the
293*> reciprocal condition number is less than the threshold
294*> sqrt(n) * dlamch('Epsilon').
295*>
296*> err = 2 "Guaranteed" error bound: The estimated forward error,
297*> almost certainly within a factor of 10 of the true error
298*> so long as the next entry is greater than the threshold
299*> sqrt(n) * dlamch('Epsilon'). This error bound should only
300*> be trusted if the previous boolean is true.
301*>
302*> err = 3 Reciprocal condition number: Estimated componentwise
303*> reciprocal condition number. Compared with the threshold
304*> sqrt(n) * dlamch('Epsilon') to determine if the error
305*> estimate is "guaranteed". These reciprocal condition
306*> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
307*> appropriately scaled matrix Z.
308*> Let Z = S*(A*diag(x)), where x is the solution for the
309*> current right-hand side and S scales each row of
310*> A*diag(x) by a power of the radix so all absolute row
311*> sums of Z are approximately 1.
312*>
313*> See Lapack Working Note 165 for further details and extra
314*> cautions.
315*> \endverbatim
316*>
317*> \param[in] NPARAMS
318*> \verbatim
319*> NPARAMS is INTEGER
320*> Specifies the number of parameters set in PARAMS. If <= 0, the
321*> PARAMS array is never referenced and default values are used.
322*> \endverbatim
323*>
324*> \param[in,out] PARAMS
325*> \verbatim
326*> PARAMS is DOUBLE PRECISION array, dimension NPARAMS
327*> Specifies algorithm parameters. If an entry is < 0.0, then
328*> that entry will be filled with default value used for that
329*> parameter. Only positions up to NPARAMS are accessed; defaults
330*> are used for higher-numbered parameters.
331*>
332*> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
333*> refinement or not.
334*> Default: 1.0D+0
335*> = 0.0: No refinement is performed, and no error bounds are
336*> computed.
337*> = 1.0: Use the double-precision refinement algorithm,
338*> possibly with doubled-single computations if the
339*> compilation environment does not support DOUBLE
340*> PRECISION.
341*> (other values are reserved for future use)
342*>
343*> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
344*> computations allowed for refinement.
345*> Default: 10
346*> Aggressive: Set to 100 to permit convergence using approximate
347*> factorizations or factorizations other than LU. If
348*> the factorization uses a technique other than
349*> Gaussian elimination, the guarantees in
350*> err_bnds_norm and err_bnds_comp may no longer be
351*> trustworthy.
352*>
353*> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
354*> will attempt to find a solution with small componentwise
355*> relative error in the double-precision algorithm. Positive
356*> is true, 0.0 is false.
357*> Default: 1.0 (attempt componentwise convergence)
358*> \endverbatim
359*>
360*> \param[out] WORK
361*> \verbatim
362*> WORK is COMPLEX*16 array, dimension (2*N)
363*> \endverbatim
364*>
365*> \param[out] RWORK
366*> \verbatim
367*> RWORK is DOUBLE PRECISION array, dimension (2*N)
368*> \endverbatim
369*>
370*> \param[out] INFO
371*> \verbatim
372*> INFO is INTEGER
373*> = 0: Successful exit. The solution to every right-hand side is
374*> guaranteed.
375*> < 0: If INFO = -i, the i-th argument had an illegal value
376*> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
377*> has been completed, but the factor U is exactly singular, so
378*> the solution and error bounds could not be computed. RCOND = 0
379*> is returned.
380*> = N+J: The solution corresponding to the Jth right-hand side is
381*> not guaranteed. The solutions corresponding to other right-
382*> hand sides K with K > J may not be guaranteed as well, but
383*> only the first such right-hand side is reported. If a small
384*> componentwise error is not requested (PARAMS(3) = 0.0) then
385*> the Jth right-hand side is the first with a normwise error
386*> bound that is not guaranteed (the smallest J such
387*> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
388*> the Jth right-hand side is the first with either a normwise or
389*> componentwise error bound that is not guaranteed (the smallest
390*> J such that either ERR_BNDS_NORM(J,1) = 0.0 or
391*> ERR_BNDS_COMP(J,1) = 0.0). See the definition of
392*> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
393*> about all of the right-hand sides check ERR_BNDS_NORM or
394*> ERR_BNDS_COMP.
395*> \endverbatim
396*
397* Authors:
398* ========
399*
400*> \author Univ. of Tennessee
401*> \author Univ. of California Berkeley
402*> \author Univ. of Colorado Denver
403*> \author NAG Ltd.
404*
405*> \ingroup gerfsx
406*
407* =====================================================================
408 SUBROUTINE zgerfsx( TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF,
409 $ IPIV,
410 $ R, C, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS,
411 $ ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
412 $ WORK, RWORK, INFO )
413*
414* -- LAPACK computational routine --
415* -- LAPACK is a software package provided by Univ. of Tennessee, --
416* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
417*
418* .. Scalar Arguments ..
419 CHARACTER TRANS, EQUED
420 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
421 $ N_ERR_BNDS
422 DOUBLE PRECISION RCOND
423* ..
424* .. Array Arguments ..
425 INTEGER IPIV( * )
426 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
427 $ X( LDX , * ), WORK( * )
428 DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
429 $ err_bnds_norm( nrhs, * ),
430 $ err_bnds_comp( nrhs, * ), rwork( * )
431* ..
432*
433* ==================================================================
434*
435* .. Parameters ..
436 DOUBLE PRECISION ZERO, ONE
437 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
438 DOUBLE PRECISION ITREF_DEFAULT, ITHRESH_DEFAULT
439 DOUBLE PRECISION COMPONENTWISE_DEFAULT, RTHRESH_DEFAULT
440 DOUBLE PRECISION DZTHRESH_DEFAULT
441 parameter( itref_default = 1.0d+0 )
442 parameter( ithresh_default = 10.0d+0 )
443 parameter( componentwise_default = 1.0d+0 )
444 parameter( rthresh_default = 0.5d+0 )
445 parameter( dzthresh_default = 0.25d+0 )
446 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
447 $ la_linrx_cwise_i
448 parameter( la_linrx_itref_i = 1,
449 $ la_linrx_ithresh_i = 2 )
450 parameter( la_linrx_cwise_i = 3 )
451 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
452 $ la_linrx_rcond_i
453 parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
454 parameter( la_linrx_rcond_i = 3 )
455* ..
456* .. Local Scalars ..
457 CHARACTER(1) NORM
458 LOGICAL ROWEQU, COLEQU, NOTRAN
459 INTEGER J, TRANS_TYPE, PREC_TYPE, REF_TYPE
460 INTEGER N_NORMS
461 DOUBLE PRECISION ANORM, RCOND_TMP
462 DOUBLE PRECISION ILLRCOND_THRESH, ERR_LBND, CWISE_WRONG
463 LOGICAL IGNORE_CWISE
464 INTEGER ITHRESH
465 DOUBLE PRECISION RTHRESH, UNSTABLE_THRESH
466* ..
467* .. External Subroutines ..
469* ..
470* .. Intrinsic Functions ..
471 INTRINSIC max, sqrt, transfer
472* ..
473* .. External Functions ..
474 EXTERNAL lsame, ilatrans, ilaprec
475 EXTERNAL dlamch, zlange, zla_gercond_x,
477 DOUBLE PRECISION DLAMCH, ZLANGE, ZLA_GERCOND_X, ZLA_GERCOND_C
478 LOGICAL LSAME
479 INTEGER ILATRANS, ILAPREC
480* ..
481* .. Executable Statements ..
482*
483* Check the input parameters.
484*
485 info = 0
486 trans_type = ilatrans( trans )
487 ref_type = int( itref_default )
488 IF ( nparams .GE. la_linrx_itref_i ) THEN
489 IF ( params( la_linrx_itref_i ) .LT. 0.0d+0 ) THEN
490 params( la_linrx_itref_i ) = itref_default
491 ELSE
492 ref_type = params( la_linrx_itref_i )
493 END IF
494 END IF
495*
496* Set default parameters.
497*
498 illrcond_thresh = dble( n ) * dlamch( 'Epsilon' )
499 ithresh = int( ithresh_default )
500 rthresh = rthresh_default
501 unstable_thresh = dzthresh_default
502 ignore_cwise = componentwise_default .EQ. 0.0d+0
503*
504 IF ( nparams.GE.la_linrx_ithresh_i ) THEN
505 IF ( params( la_linrx_ithresh_i ).LT.0.0d+0 ) THEN
506 params(la_linrx_ithresh_i) = ithresh
507 ELSE
508 ithresh = int( params( la_linrx_ithresh_i ) )
509 END IF
510 END IF
511 IF ( nparams.GE.la_linrx_cwise_i ) THEN
512 IF ( params( la_linrx_cwise_i ).LT.0.0d+0 ) THEN
513 IF ( ignore_cwise ) THEN
514 params( la_linrx_cwise_i ) = 0.0d+0
515 ELSE
516 params( la_linrx_cwise_i ) = 1.0d+0
517 END IF
518 ELSE
519 ignore_cwise = params( la_linrx_cwise_i ) .EQ. 0.0d+0
520 END IF
521 END IF
522 IF ( ref_type .EQ. 0 .OR. n_err_bnds .EQ. 0 ) THEN
523 n_norms = 0
524 ELSE IF ( ignore_cwise ) THEN
525 n_norms = 1
526 ELSE
527 n_norms = 2
528 END IF
529*
530 notran = lsame( trans, 'N' )
531 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
532 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
533*
534* Test input parameters.
535*
536 IF( trans_type.EQ.-1 ) THEN
537 info = -1
538 ELSE IF( .NOT.rowequ .AND. .NOT.colequ .AND.
539 $ .NOT.lsame( equed, 'N' ) ) THEN
540 info = -2
541 ELSE IF( n.LT.0 ) THEN
542 info = -3
543 ELSE IF( nrhs.LT.0 ) THEN
544 info = -4
545 ELSE IF( lda.LT.max( 1, n ) ) THEN
546 info = -6
547 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
548 info = -8
549 ELSE IF( ldb.LT.max( 1, n ) ) THEN
550 info = -13
551 ELSE IF( ldx.LT.max( 1, n ) ) THEN
552 info = -15
553 END IF
554 IF( info.NE.0 ) THEN
555 CALL xerbla( 'ZGERFSX', -info )
556 RETURN
557 END IF
558*
559* Quick return if possible.
560*
561 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
562 rcond = 1.0d+0
563 DO j = 1, nrhs
564 berr( j ) = 0.0d+0
565 IF ( n_err_bnds .GE. 1 ) THEN
566 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
567 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
568 END IF
569 IF ( n_err_bnds .GE. 2 ) THEN
570 err_bnds_norm( j, la_linrx_err_i ) = 0.0d+0
571 err_bnds_comp( j, la_linrx_err_i ) = 0.0d+0
572 END IF
573 IF ( n_err_bnds .GE. 3 ) THEN
574 err_bnds_norm( j, la_linrx_rcond_i ) = 1.0d+0
575 err_bnds_comp( j, la_linrx_rcond_i ) = 1.0d+0
576 END IF
577 END DO
578 RETURN
579 END IF
580*
581* Default to failure.
582*
583 rcond = 0.0d+0
584 DO j = 1, nrhs
585 berr( j ) = 1.0d+0
586 IF ( n_err_bnds .GE. 1 ) THEN
587 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
588 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
589 END IF
590 IF ( n_err_bnds .GE. 2 ) THEN
591 err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
592 err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
593 END IF
594 IF ( n_err_bnds .GE. 3 ) THEN
595 err_bnds_norm( j, la_linrx_rcond_i ) = 0.0d+0
596 err_bnds_comp( j, la_linrx_rcond_i ) = 0.0d+0
597 END IF
598 END DO
599*
600* Compute the norm of A and the reciprocal of the condition
601* number of A.
602*
603 IF( notran ) THEN
604 norm = 'I'
605 ELSE
606 norm = '1'
607 END IF
608 anorm = zlange( norm, n, n, a, lda, rwork )
609 CALL zgecon( norm, n, af, ldaf, anorm, rcond, work, rwork,
610 $ info )
611*
612* Perform refinement on each right-hand side
613*
614 IF ( ref_type .NE. 0 ) THEN
615
616 prec_type = ilaprec( 'E' )
617
618 IF ( notran ) THEN
619 CALL zla_gerfsx_extended( prec_type, trans_type, n,
620 $ nrhs, a, lda, af, ldaf, ipiv, colequ, c, b,
621 $ ldb, x, ldx, berr, n_norms, err_bnds_norm,
622 $ err_bnds_comp, work, rwork, work(n+1),
623 $ transfer(rwork(1:2*n), (/ (zero, zero) /), n),
624 $ rcond, ithresh, rthresh, unstable_thresh, ignore_cwise,
625 $ info )
626 ELSE
627 CALL zla_gerfsx_extended( prec_type, trans_type, n,
628 $ nrhs, a, lda, af, ldaf, ipiv, rowequ, r, b,
629 $ ldb, x, ldx, berr, n_norms, err_bnds_norm,
630 $ err_bnds_comp, work, rwork, work(n+1),
631 $ transfer(rwork(1:2*n), (/ (zero, zero) /), n),
632 $ rcond, ithresh, rthresh, unstable_thresh, ignore_cwise,
633 $ info )
634 END IF
635 END IF
636
637 err_lbnd = max( 10.0d+0,
638 $ sqrt( dble( n ) ) ) * dlamch( 'Epsilon' )
639 IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 1 ) THEN
640*
641* Compute scaled normwise condition number cond(A*C).
642*
643 IF ( colequ .AND. notran ) THEN
644 rcond_tmp = zla_gercond_c( trans, n, a, lda, af, ldaf,
645 $ ipiv,
646 $ c, .true., info, work, rwork )
647 ELSE IF ( rowequ .AND. .NOT. notran ) THEN
648 rcond_tmp = zla_gercond_c( trans, n, a, lda, af, ldaf,
649 $ ipiv,
650 $ r, .true., info, work, rwork )
651 ELSE
652 rcond_tmp = zla_gercond_c( trans, n, a, lda, af, ldaf,
653 $ ipiv,
654 $ c, .false., info, work, rwork )
655 END IF
656 DO j = 1, nrhs
657*
658* Cap the error at 1.0.
659*
660 IF ( n_err_bnds .GE. la_linrx_err_i
661 $ .AND. err_bnds_norm( j, la_linrx_err_i ) .GT. 1.0d+0 )
662 $ err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
663*
664* Threshold the error (see LAWN).
665*
666 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
667 err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
668 err_bnds_norm( j, la_linrx_trust_i ) = 0.0d+0
669 IF ( info .LE. n ) info = n + j
670 ELSE IF (err_bnds_norm( j, la_linrx_err_i ) .LT. err_lbnd)
671 $ THEN
672 err_bnds_norm( j, la_linrx_err_i ) = err_lbnd
673 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
674 END IF
675*
676* Save the condition number.
677*
678 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
679 err_bnds_norm( j, la_linrx_rcond_i ) = rcond_tmp
680 END IF
681 END DO
682 END IF
683
684 IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 2 ) THEN
685*
686* Compute componentwise condition number cond(A*diag(Y(:,J))) for
687* each right-hand side using the current solution as an estimate of
688* the true solution. If the componentwise error estimate is too
689* large, then the solution is a lousy estimate of truth and the
690* estimated RCOND may be too optimistic. To avoid misleading users,
691* the inverse condition number is set to 0.0 when the estimated
692* cwise error is at least CWISE_WRONG.
693*
694 cwise_wrong = sqrt( dlamch( 'Epsilon' ) )
695 DO j = 1, nrhs
696 IF ( err_bnds_comp( j, la_linrx_err_i ) .LT. cwise_wrong )
697 $ THEN
698 rcond_tmp = zla_gercond_x( trans, n, a, lda, af, ldaf,
699 $ ipiv, x(1,j), info, work, rwork )
700 ELSE
701 rcond_tmp = 0.0d+0
702 END IF
703*
704* Cap the error at 1.0.
705*
706 IF ( n_err_bnds .GE. la_linrx_err_i
707 $ .AND. err_bnds_comp( j, la_linrx_err_i ) .GT. 1.0d+0 )
708 $ err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
709*
710* Threshold the error (see LAWN).
711*
712 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
713 err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
714 err_bnds_comp( j, la_linrx_trust_i ) = 0.0d+0
715 IF ( params( la_linrx_cwise_i ) .EQ. 1.0d+0
716 $ .AND. info.LT.n + j ) info = n + j
717 ELSE IF ( err_bnds_comp( j, la_linrx_err_i )
718 $ .LT. err_lbnd ) THEN
719 err_bnds_comp( j, la_linrx_err_i ) = err_lbnd
720 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
721 END IF
722*
723* Save the condition number.
724*
725 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
726 err_bnds_comp( j, la_linrx_rcond_i ) = rcond_tmp
727 END IF
728
729 END DO
730 END IF
731*
732 RETURN
733*
734* End of ZGERFSX
735*
736 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgecon(norm, n, a, lda, anorm, rcond, work, rwork, info)
ZGECON
Definition zgecon.f:130
subroutine zgerfsx(trans, equed, n, nrhs, a, lda, af, ldaf, ipiv, r, c, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
ZGERFSX
Definition zgerfsx.f:413
integer function ilaprec(prec)
ILAPREC
Definition ilaprec.f:56
integer function ilatrans(trans)
ILATRANS
Definition ilatrans.f:56
double precision function zla_gercond_c(trans, n, a, lda, af, ldaf, ipiv, c, capply, info, work, rwork)
ZLA_GERCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general matrices.
double precision function zla_gercond_x(trans, n, a, lda, af, ldaf, ipiv, x, info, work, rwork)
ZLA_GERCOND_X computes the infinity norm condition number of op(A)*diag(x) for general matrices.
subroutine zla_gerfsx_extended(prec_type, trans_type, n, nrhs, a, lda, af, ldaf, ipiv, colequ, c, b, ldb, y, ldy, berr_out, n_norms, errs_n, errs_c, res, ayb, dy, y_tail, rcond, ithresh, rthresh, dz_ub, ignore_cwise, info)
ZLA_GERFSX_EXTENDED
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48