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