LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssyrfsx.f
Go to the documentation of this file.
1*> \brief \b SSYRFSX
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SSYRFSX + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssyrfsx.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssyrfsx.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssyrfsx.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SSYRFSX( 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, IWORK, 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( * ), IWORK( * )
34* REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
35* $ X( LDX, * ), WORK( * )
36* REAL S( * ), PARAMS( * ), BERR( * ),
37* $ ERR_BNDS_NORM( NRHS, * ),
38* $ ERR_BNDS_COMP( NRHS, * )
39* ..
40*
41*
42*> \par Purpose:
43* =============
44*>
45*> \verbatim
46*>
47*> SSYRFSX improves the computed solution to a system of linear
48*> equations when the coefficient matrix is symmetric indefinite, and
49*> provides error bounds and backward error estimates for the
50*> solution. In addition to normwise error bound, the code provides
51*> maximum componentwise error bound if possible. See comments for
52*> ERR_BNDS_NORM and ERR_BNDS_COMP for details of the error bounds.
53*>
54*> The original system of linear equations may have been equilibrated
55*> before calling this routine, as described by arguments EQUED and S
56*> below. In this case, the solution and error bounds returned are
57*> for the original unequilibrated system.
58*> \endverbatim
59*
60* Arguments:
61* ==========
62*
63*> \verbatim
64*> Some optional parameters are bundled in the PARAMS array. These
65*> settings determine how refinement is performed, but often the
66*> defaults are acceptable. If the defaults are acceptable, users
67*> can pass NPARAMS = 0 which prevents the source code from accessing
68*> the PARAMS argument.
69*> \endverbatim
70*>
71*> \param[in] UPLO
72*> \verbatim
73*> UPLO is CHARACTER*1
74*> = 'U': Upper triangle of A is stored;
75*> = 'L': Lower triangle of A is stored.
76*> \endverbatim
77*>
78*> \param[in] EQUED
79*> \verbatim
80*> EQUED is CHARACTER*1
81*> Specifies the form of equilibration that was done to A
82*> before calling this routine. This is needed to compute
83*> the solution and error bounds correctly.
84*> = 'N': No equilibration
85*> = 'Y': Both row and column equilibration, i.e., A has been
86*> replaced by diag(S) * A * diag(S).
87*> The right hand side B has been changed accordingly.
88*> \endverbatim
89*>
90*> \param[in] N
91*> \verbatim
92*> N is INTEGER
93*> The order of the matrix A. N >= 0.
94*> \endverbatim
95*>
96*> \param[in] NRHS
97*> \verbatim
98*> NRHS is INTEGER
99*> The number of right hand sides, i.e., the number of columns
100*> of the matrices B and X. NRHS >= 0.
101*> \endverbatim
102*>
103*> \param[in] A
104*> \verbatim
105*> A is REAL array, dimension (LDA,N)
106*> The symmetric matrix A. If UPLO = 'U', the leading N-by-N
107*> upper triangular part of A contains the upper triangular
108*> part of the matrix A, and the strictly lower triangular
109*> part of A is not referenced. If UPLO = 'L', the leading
110*> N-by-N lower triangular part of A contains the lower
111*> triangular part of the matrix A, and the strictly upper
112*> triangular part of A is not referenced.
113*> \endverbatim
114*>
115*> \param[in] LDA
116*> \verbatim
117*> LDA is INTEGER
118*> The leading dimension of the array A. LDA >= max(1,N).
119*> \endverbatim
120*>
121*> \param[in] AF
122*> \verbatim
123*> AF is REAL array, dimension (LDAF,N)
124*> The factored form of the matrix A. AF contains the block
125*> diagonal matrix D and the multipliers used to obtain the
126*> factor U or L from the factorization A = U*D*U**T or A =
127*> L*D*L**T as computed by SSYTRF.
128*> \endverbatim
129*>
130*> \param[in] LDAF
131*> \verbatim
132*> LDAF is INTEGER
133*> The leading dimension of the array AF. LDAF >= max(1,N).
134*> \endverbatim
135*>
136*> \param[in] IPIV
137*> \verbatim
138*> IPIV is INTEGER array, dimension (N)
139*> Details of the interchanges and the block structure of D
140*> as determined by SSYTRF.
141*> \endverbatim
142*>
143*> \param[in,out] S
144*> \verbatim
145*> S is REAL array, dimension (N)
146*> The scale factors for A. If EQUED = 'Y', A is multiplied on
147*> the left and right by diag(S). S is an input argument if FACT =
148*> 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED
149*> = 'Y', each element of S must be positive. If S is output, each
150*> element of S is a power of the radix. If S is input, each element
151*> of S should be a power of the radix to ensure a reliable solution
152*> and error estimates. Scaling by powers of the radix does not cause
153*> rounding errors unless the result underflows or overflows.
154*> Rounding errors during scaling lead to refining with a matrix that
155*> is not equivalent to the input matrix, producing error estimates
156*> that may not be reliable.
157*> \endverbatim
158*>
159*> \param[in] B
160*> \verbatim
161*> B is REAL array, dimension (LDB,NRHS)
162*> The right hand side matrix B.
163*> \endverbatim
164*>
165*> \param[in] LDB
166*> \verbatim
167*> LDB is INTEGER
168*> The leading dimension of the array B. LDB >= max(1,N).
169*> \endverbatim
170*>
171*> \param[in,out] X
172*> \verbatim
173*> X is REAL array, dimension (LDX,NRHS)
174*> On entry, the solution matrix X, as computed by SGETRS.
175*> On exit, the improved solution matrix X.
176*> \endverbatim
177*>
178*> \param[in] LDX
179*> \verbatim
180*> LDX is INTEGER
181*> The leading dimension of the array X. LDX >= max(1,N).
182*> \endverbatim
183*>
184*> \param[out] RCOND
185*> \verbatim
186*> RCOND is REAL
187*> Reciprocal scaled condition number. This is an estimate of the
188*> reciprocal Skeel condition number of the matrix A after
189*> equilibration (if done). If this is less than the machine
190*> precision (in particular, if it is zero), the matrix is singular
191*> to working precision. Note that the error may still be small even
192*> if this number is very small and the matrix appears ill-
193*> conditioned.
194*> \endverbatim
195*>
196*> \param[out] BERR
197*> \verbatim
198*> BERR is REAL array, dimension (NRHS)
199*> Componentwise relative backward error. This is the
200*> componentwise relative backward error of each solution vector X(j)
201*> (i.e., the smallest relative change in any element of A or B that
202*> makes X(j) an exact solution).
203*> \endverbatim
204*>
205*> \param[in] N_ERR_BNDS
206*> \verbatim
207*> N_ERR_BNDS is INTEGER
208*> Number of error bounds to return for each right hand side
209*> and each type (normwise or componentwise). See ERR_BNDS_NORM and
210*> ERR_BNDS_COMP below.
211*> \endverbatim
212*>
213*> \param[out] ERR_BNDS_NORM
214*> \verbatim
215*> ERR_BNDS_NORM is REAL array, dimension (NRHS, N_ERR_BNDS)
216*> For each right-hand side, this array contains information about
217*> various error bounds and condition numbers corresponding to the
218*> normwise relative error, which is defined as follows:
219*>
220*> Normwise relative error in the ith solution vector:
221*> max_j (abs(XTRUE(j,i) - X(j,i)))
222*> ------------------------------
223*> max_j abs(X(j,i))
224*>
225*> The array is indexed by the type of error information as described
226*> below. There currently are up to three pieces of information
227*> returned.
228*>
229*> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
230*> right-hand side.
231*>
232*> The second index in ERR_BNDS_NORM(:,err) contains the following
233*> three fields:
234*> err = 1 "Trust/don't trust" boolean. Trust the answer if the
235*> reciprocal condition number is less than the threshold
236*> sqrt(n) * slamch('Epsilon').
237*>
238*> err = 2 "Guaranteed" error bound: The estimated forward error,
239*> almost certainly within a factor of 10 of the true error
240*> so long as the next entry is greater than the threshold
241*> sqrt(n) * slamch('Epsilon'). This error bound should only
242*> be trusted if the previous boolean is true.
243*>
244*> err = 3 Reciprocal condition number: Estimated normwise
245*> reciprocal condition number. Compared with the threshold
246*> sqrt(n) * slamch('Epsilon') to determine if the error
247*> estimate is "guaranteed". These reciprocal condition
248*> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
249*> appropriately scaled matrix Z.
250*> Let Z = S*A, where S scales each row by a power of the
251*> radix so all absolute row sums of Z are approximately 1.
252*>
253*> See Lapack Working Note 165 for further details and extra
254*> cautions.
255*> \endverbatim
256*>
257*> \param[out] ERR_BNDS_COMP
258*> \verbatim
259*> ERR_BNDS_COMP is REAL array, dimension (NRHS, N_ERR_BNDS)
260*> For each right-hand side, this array contains information about
261*> various error bounds and condition numbers corresponding to the
262*> componentwise relative error, which is defined as follows:
263*>
264*> Componentwise relative error in the ith solution vector:
265*> abs(XTRUE(j,i) - X(j,i))
266*> max_j ----------------------
267*> abs(X(j,i))
268*>
269*> The array is indexed by the right-hand side i (on which the
270*> componentwise relative error depends), and the type of error
271*> information as described below. There currently are up to three
272*> pieces of information returned for each right-hand side. If
273*> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
274*> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS < 3, then at most
275*> the first (:,N_ERR_BNDS) entries are returned.
276*>
277*> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
278*> right-hand side.
279*>
280*> The second index in ERR_BNDS_COMP(:,err) contains the following
281*> three fields:
282*> err = 1 "Trust/don't trust" boolean. Trust the answer if the
283*> reciprocal condition number is less than the threshold
284*> sqrt(n) * slamch('Epsilon').
285*>
286*> err = 2 "Guaranteed" error bound: The estimated forward error,
287*> almost certainly within a factor of 10 of the true error
288*> so long as the next entry is greater than the threshold
289*> sqrt(n) * slamch('Epsilon'). This error bound should only
290*> be trusted if the previous boolean is true.
291*>
292*> err = 3 Reciprocal condition number: Estimated componentwise
293*> reciprocal condition number. Compared with the threshold
294*> sqrt(n) * slamch('Epsilon') to determine if the error
295*> estimate is "guaranteed". These reciprocal condition
296*> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
297*> appropriately scaled matrix Z.
298*> Let Z = S*(A*diag(x)), where x is the solution for the
299*> current right-hand side and S scales each row of
300*> A*diag(x) by a power of the radix so all absolute row
301*> sums of Z are approximately 1.
302*>
303*> See Lapack Working Note 165 for further details and extra
304*> cautions.
305*> \endverbatim
306*>
307*> \param[in] NPARAMS
308*> \verbatim
309*> NPARAMS is INTEGER
310*> Specifies the number of parameters set in PARAMS. If <= 0, the
311*> PARAMS array is never referenced and default values are used.
312*> \endverbatim
313*>
314*> \param[in,out] PARAMS
315*> \verbatim
316*> PARAMS is REAL array, dimension NPARAMS
317*> Specifies algorithm parameters. If an entry is < 0.0, then
318*> that entry will be filled with default value used for that
319*> parameter. Only positions up to NPARAMS are accessed; defaults
320*> are used for higher-numbered parameters.
321*>
322*> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
323*> refinement or not.
324*> Default: 1.0
325*> = 0.0: No refinement is performed, and no error bounds are
326*> computed.
327*> = 1.0: Use the double-precision refinement algorithm,
328*> possibly with doubled-single computations if the
329*> compilation environment does not support DOUBLE
330*> PRECISION.
331*> (other values are reserved for future use)
332*>
333*> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
334*> computations allowed for refinement.
335*> Default: 10
336*> Aggressive: Set to 100 to permit convergence using approximate
337*> factorizations or factorizations other than LU. If
338*> the factorization uses a technique other than
339*> Gaussian elimination, the guarantees in
340*> err_bnds_norm and err_bnds_comp may no longer be
341*> trustworthy.
342*>
343*> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
344*> will attempt to find a solution with small componentwise
345*> relative error in the double-precision algorithm. Positive
346*> is true, 0.0 is false.
347*> Default: 1.0 (attempt componentwise convergence)
348*> \endverbatim
349*>
350*> \param[out] WORK
351*> \verbatim
352*> WORK is REAL array, dimension (4*N)
353*> \endverbatim
354*>
355*> \param[out] IWORK
356*> \verbatim
357*> IWORK is INTEGER array, dimension (N)
358*> \endverbatim
359*>
360*> \param[out] INFO
361*> \verbatim
362*> INFO is INTEGER
363*> = 0: Successful exit. The solution to every right-hand side is
364*> guaranteed.
365*> < 0: If INFO = -i, the i-th argument had an illegal value
366*> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
367*> has been completed, but the factor U is exactly singular, so
368*> the solution and error bounds could not be computed. RCOND = 0
369*> is returned.
370*> = N+J: The solution corresponding to the Jth right-hand side is
371*> not guaranteed. The solutions corresponding to other right-
372*> hand sides K with K > J may not be guaranteed as well, but
373*> only the first such right-hand side is reported. If a small
374*> componentwise error is not requested (PARAMS(3) = 0.0) then
375*> the Jth right-hand side is the first with a normwise error
376*> bound that is not guaranteed (the smallest J such
377*> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
378*> the Jth right-hand side is the first with either a normwise or
379*> componentwise error bound that is not guaranteed (the smallest
380*> J such that either ERR_BNDS_NORM(J,1) = 0.0 or
381*> ERR_BNDS_COMP(J,1) = 0.0). See the definition of
382*> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
383*> about all of the right-hand sides check ERR_BNDS_NORM or
384*> ERR_BNDS_COMP.
385*> \endverbatim
386*
387* Authors:
388* ========
389*
390*> \author Univ. of Tennessee
391*> \author Univ. of California Berkeley
392*> \author Univ. of Colorado Denver
393*> \author NAG Ltd.
394*
395*> \ingroup herfsx
396*
397* =====================================================================
398 SUBROUTINE ssyrfsx( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV,
399 $ S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS,
400 $ ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
401 $ WORK, IWORK, INFO )
402*
403* -- LAPACK computational routine --
404* -- LAPACK is a software package provided by Univ. of Tennessee, --
405* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
406*
407* .. Scalar Arguments ..
408 CHARACTER UPLO, EQUED
409 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
410 $ N_ERR_BNDS
411 REAL RCOND
412* ..
413* .. Array Arguments ..
414 INTEGER IPIV( * ), IWORK( * )
415 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
416 $ X( LDX, * ), WORK( * )
417 REAL S( * ), PARAMS( * ), BERR( * ),
418 $ err_bnds_norm( nrhs, * ),
419 $ err_bnds_comp( nrhs, * )
420* ..
421*
422* ==================================================================
423*
424* .. Parameters ..
425 REAL ZERO, ONE
426 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
427 REAL ITREF_DEFAULT, ITHRESH_DEFAULT,
428 $ componentwise_default
429 REAL RTHRESH_DEFAULT, DZTHRESH_DEFAULT
430 parameter( itref_default = 1.0 )
431 parameter( ithresh_default = 10.0 )
432 parameter( componentwise_default = 1.0 )
433 parameter( rthresh_default = 0.5 )
434 parameter( dzthresh_default = 0.25 )
435 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
436 $ la_linrx_cwise_i
437 parameter( la_linrx_itref_i = 1,
438 $ la_linrx_ithresh_i = 2 )
439 parameter( la_linrx_cwise_i = 3 )
440 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
441 $ LA_LINRX_RCOND_I
442 parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
443 parameter( la_linrx_rcond_i = 3 )
444* ..
445* .. Local Scalars ..
446 CHARACTER(1) NORM
447 LOGICAL RCEQU
448 INTEGER J, PREC_TYPE, REF_TYPE, 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
460* ..
461* .. External Functions ..
462 EXTERNAL lsame, ilaprec
463 EXTERNAL slamch, slansy, sla_syrcond
464 REAL SLAMCH, SLANSY, SLA_SYRCOND
465 LOGICAL LSAME
466 INTEGER ILAPREC
467* ..
468* .. Executable Statements ..
469*
470* Check the input parameters.
471*
472 info = 0
473 ref_type = int( itref_default )
474 IF ( nparams .GE. la_linrx_itref_i ) THEN
475 IF ( params( la_linrx_itref_i ) .LT. 0.0 ) THEN
476 params( la_linrx_itref_i ) = itref_default
477 ELSE
478 ref_type = params( la_linrx_itref_i )
479 END IF
480 END IF
481*
482* Set default parameters.
483*
484 illrcond_thresh = real( n )*slamch( 'Epsilon' )
485 ithresh = int( ithresh_default )
486 rthresh = rthresh_default
487 unstable_thresh = dzthresh_default
488 ignore_cwise = componentwise_default .EQ. 0.0
489*
490 IF ( nparams.GE.la_linrx_ithresh_i ) THEN
491 IF ( params( la_linrx_ithresh_i ).LT.0.0 ) THEN
492 params( la_linrx_ithresh_i ) = ithresh
493 ELSE
494 ithresh = int( params( la_linrx_ithresh_i ) )
495 END IF
496 END IF
497 IF ( nparams.GE.la_linrx_cwise_i ) THEN
498 IF ( params( la_linrx_cwise_i ).LT.0.0 ) THEN
499 IF ( ignore_cwise ) THEN
500 params( la_linrx_cwise_i ) = 0.0
501 ELSE
502 params( la_linrx_cwise_i ) = 1.0
503 END IF
504 ELSE
505 ignore_cwise = params( la_linrx_cwise_i ) .EQ. 0.0
506 END IF
507 END IF
508 IF ( ref_type .EQ. 0 .OR. n_err_bnds .EQ. 0 ) THEN
509 n_norms = 0
510 ELSE IF ( ignore_cwise ) THEN
511 n_norms = 1
512 ELSE
513 n_norms = 2
514 END IF
515*
516 rcequ = lsame( equed, 'Y' )
517*
518* Test input parameters.
519*
520 IF ( .NOT.lsame( uplo, 'U' ) .AND. .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( 'SSYRFSX', -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 = slansy( norm, uplo, n, a, lda, work )
588 CALL ssycon( uplo, n, af, ldaf, ipiv, anorm, rcond, work,
589 $ iwork, 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 sla_syrfsx_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( n+1 ), work( 1 ), work( 2*n+1 ), work( 1 ), 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 = sla_syrcond( uplo, n, a, lda, af, ldaf, ipiv,
612 $ -1, s, info, work, iwork )
613 ELSE
614 rcond_tmp = sla_syrcond( uplo, n, a, lda, af, ldaf, ipiv,
615 $ 0, s, info, work, iwork )
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 = sla_syrcond( uplo, n, a, lda, af, ldaf, ipiv,
660 $ 1, x(1,j), info, work, iwork )
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 SSYRFSX
696*
697 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssycon(uplo, n, a, lda, ipiv, anorm, rcond, work, iwork, info)
SSYCON
Definition ssycon.f:130
subroutine ssyrfsx(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, iwork, info)
SSYRFSX
Definition ssyrfsx.f:402
integer function ilaprec(prec)
ILAPREC
Definition ilaprec.f:58
real function sla_syrcond(uplo, n, a, lda, af, ldaf, ipiv, cmode, c, info, work, iwork)
SLA_SYRCOND estimates the Skeel condition number for a symmetric indefinite matrix.
subroutine sla_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)
SLA_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 slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:122
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48