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