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