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