LAPACK  3.9.0
LAPACK: Linear Algebra PACKage
cla_syrfsx_extended.f
Go to the documentation of this file.
1 *> \brief \b CLA_SYRFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric indefinite matrices by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution.
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download CLA_SYRFSX_EXTENDED + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cla_syrfsx_extended.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cla_syrfsx_extended.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cla_syrfsx_extended.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE CLA_SYRFSX_EXTENDED( PREC_TYPE, UPLO, N, NRHS, A, LDA,
22 * AF, LDAF, IPIV, COLEQU, C, B, LDB,
23 * Y, LDY, BERR_OUT, N_NORMS,
24 * ERR_BNDS_NORM, ERR_BNDS_COMP, RES,
25 * AYB, DY, Y_TAIL, RCOND, ITHRESH,
26 * RTHRESH, DZ_UB, IGNORE_CWISE,
27 * INFO )
28 *
29 * .. Scalar Arguments ..
30 * INTEGER INFO, LDA, LDAF, LDB, LDY, N, NRHS, PREC_TYPE,
31 * $ N_NORMS, ITHRESH
32 * CHARACTER UPLO
33 * LOGICAL COLEQU, IGNORE_CWISE
34 * REAL RTHRESH, DZ_UB
35 * ..
36 * .. Array Arguments ..
37 * INTEGER IPIV( * )
38 * COMPLEX A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
39 * $ Y( LDY, * ), RES( * ), DY( * ), Y_TAIL( * )
40 * REAL C( * ), AYB( * ), RCOND, BERR_OUT( * ),
41 * $ ERR_BNDS_NORM( NRHS, * ),
42 * $ ERR_BNDS_COMP( NRHS, * )
43 * ..
44 *
45 *
46 *> \par Purpose:
47 * =============
48 *>
49 *> \verbatim
50 *>
51 *> CLA_SYRFSX_EXTENDED improves the computed solution to a system of
52 *> linear equations by performing extra-precise iterative refinement
53 *> and provides error bounds and backward error estimates for the solution.
54 *> This subroutine is called by CSYRFSX to perform iterative refinement.
55 *> In addition to normwise error bound, the code provides maximum
56 *> componentwise error bound if possible. See comments for ERR_BNDS_NORM
57 *> and ERR_BNDS_COMP for details of the error bounds. Note that this
58 *> subroutine is only resonsible for setting the second fields of
59 *> ERR_BNDS_NORM and ERR_BNDS_COMP.
60 *> \endverbatim
61 *
62 * Arguments:
63 * ==========
64 *
65 *> \param[in] PREC_TYPE
66 *> \verbatim
67 *> PREC_TYPE is INTEGER
68 *> Specifies the intermediate precision to be used in refinement.
69 *> The value is defined by ILAPREC(P) where P is a CHARACTER and P
70 *> = 'S': Single
71 *> = 'D': Double
72 *> = 'I': Indigenous
73 *> = 'X' or 'E': Extra
74 *> \endverbatim
75 *>
76 *> \param[in] UPLO
77 *> \verbatim
78 *> UPLO is CHARACTER*1
79 *> = 'U': Upper triangle of A is stored;
80 *> = 'L': Lower triangle of A is stored.
81 *> \endverbatim
82 *>
83 *> \param[in] N
84 *> \verbatim
85 *> N is INTEGER
86 *> The number of linear equations, i.e., the order of the
87 *> matrix A. N >= 0.
88 *> \endverbatim
89 *>
90 *> \param[in] NRHS
91 *> \verbatim
92 *> NRHS is INTEGER
93 *> The number of right-hand-sides, i.e., the number of columns of the
94 *> matrix B.
95 *> \endverbatim
96 *>
97 *> \param[in] A
98 *> \verbatim
99 *> A is COMPLEX array, dimension (LDA,N)
100 *> On entry, the N-by-N matrix A.
101 *> \endverbatim
102 *>
103 *> \param[in] LDA
104 *> \verbatim
105 *> LDA is INTEGER
106 *> The leading dimension of the array A. LDA >= max(1,N).
107 *> \endverbatim
108 *>
109 *> \param[in] AF
110 *> \verbatim
111 *> AF is COMPLEX array, dimension (LDAF,N)
112 *> The block diagonal matrix D and the multipliers used to
113 *> obtain the factor U or L as computed by CSYTRF.
114 *> \endverbatim
115 *>
116 *> \param[in] LDAF
117 *> \verbatim
118 *> LDAF is INTEGER
119 *> The leading dimension of the array AF. LDAF >= max(1,N).
120 *> \endverbatim
121 *>
122 *> \param[in] IPIV
123 *> \verbatim
124 *> IPIV is INTEGER array, dimension (N)
125 *> Details of the interchanges and the block structure of D
126 *> as determined by CSYTRF.
127 *> \endverbatim
128 *>
129 *> \param[in] COLEQU
130 *> \verbatim
131 *> COLEQU is LOGICAL
132 *> If .TRUE. then column equilibration was done to A before calling
133 *> this routine. This is needed to compute the solution and error
134 *> bounds correctly.
135 *> \endverbatim
136 *>
137 *> \param[in] C
138 *> \verbatim
139 *> C is REAL array, dimension (N)
140 *> The column scale factors for A. If COLEQU = .FALSE., C
141 *> is not accessed. If C is input, each element of C should be a power
142 *> of the radix to ensure a reliable solution and error estimates.
143 *> Scaling by powers of the radix does not cause rounding errors unless
144 *> the result underflows or overflows. Rounding errors during scaling
145 *> lead to refining with a matrix that is not equivalent to the
146 *> input matrix, producing error estimates that may not be
147 *> reliable.
148 *> \endverbatim
149 *>
150 *> \param[in] B
151 *> \verbatim
152 *> B is COMPLEX array, dimension (LDB,NRHS)
153 *> The right-hand-side matrix B.
154 *> \endverbatim
155 *>
156 *> \param[in] LDB
157 *> \verbatim
158 *> LDB is INTEGER
159 *> The leading dimension of the array B. LDB >= max(1,N).
160 *> \endverbatim
161 *>
162 *> \param[in,out] Y
163 *> \verbatim
164 *> Y is COMPLEX array, dimension (LDY,NRHS)
165 *> On entry, the solution matrix X, as computed by CSYTRS.
166 *> On exit, the improved solution matrix Y.
167 *> \endverbatim
168 *>
169 *> \param[in] LDY
170 *> \verbatim
171 *> LDY is INTEGER
172 *> The leading dimension of the array Y. LDY >= max(1,N).
173 *> \endverbatim
174 *>
175 *> \param[out] BERR_OUT
176 *> \verbatim
177 *> BERR_OUT is REAL array, dimension (NRHS)
178 *> On exit, BERR_OUT(j) contains the componentwise relative backward
179 *> error for right-hand-side j from the formula
180 *> max(i) ( abs(RES(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
181 *> where abs(Z) is the componentwise absolute value of the matrix
182 *> or vector Z. This is computed by CLA_LIN_BERR.
183 *> \endverbatim
184 *>
185 *> \param[in] N_NORMS
186 *> \verbatim
187 *> N_NORMS is INTEGER
188 *> Determines which error bounds to return (see ERR_BNDS_NORM
189 *> and ERR_BNDS_COMP).
190 *> If N_NORMS >= 1 return normwise error bounds.
191 *> If N_NORMS >= 2 return componentwise error bounds.
192 *> \endverbatim
193 *>
194 *> \param[in,out] ERR_BNDS_NORM
195 *> \verbatim
196 *> ERR_BNDS_NORM is REAL array, dimension (NRHS, N_ERR_BNDS)
197 *> For each right-hand side, this array contains information about
198 *> various error bounds and condition numbers corresponding to the
199 *> normwise relative error, which is defined as follows:
200 *>
201 *> Normwise relative error in the ith solution vector:
202 *> max_j (abs(XTRUE(j,i) - X(j,i)))
203 *> ------------------------------
204 *> max_j abs(X(j,i))
205 *>
206 *> The array is indexed by the type of error information as described
207 *> below. There currently are up to three pieces of information
208 *> returned.
209 *>
210 *> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
211 *> right-hand side.
212 *>
213 *> The second index in ERR_BNDS_NORM(:,err) contains the following
214 *> three fields:
215 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
216 *> reciprocal condition number is less than the threshold
217 *> sqrt(n) * slamch('Epsilon').
218 *>
219 *> err = 2 "Guaranteed" error bound: The estimated forward error,
220 *> almost certainly within a factor of 10 of the true error
221 *> so long as the next entry is greater than the threshold
222 *> sqrt(n) * slamch('Epsilon'). This error bound should only
223 *> be trusted if the previous boolean is true.
224 *>
225 *> err = 3 Reciprocal condition number: Estimated normwise
226 *> reciprocal condition number. Compared with the threshold
227 *> sqrt(n) * slamch('Epsilon') to determine if the error
228 *> estimate is "guaranteed". These reciprocal condition
229 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
230 *> appropriately scaled matrix Z.
231 *> Let Z = S*A, where S scales each row by a power of the
232 *> radix so all absolute row sums of Z are approximately 1.
233 *>
234 *> This subroutine is only responsible for setting the second field
235 *> above.
236 *> See Lapack Working Note 165 for further details and extra
237 *> cautions.
238 *> \endverbatim
239 *>
240 *> \param[in,out] ERR_BNDS_COMP
241 *> \verbatim
242 *> ERR_BNDS_COMP is REAL array, dimension (NRHS, N_ERR_BNDS)
243 *> For each right-hand side, this array contains information about
244 *> various error bounds and condition numbers corresponding to the
245 *> componentwise relative error, which is defined as follows:
246 *>
247 *> Componentwise relative error in the ith solution vector:
248 *> abs(XTRUE(j,i) - X(j,i))
249 *> max_j ----------------------
250 *> abs(X(j,i))
251 *>
252 *> The array is indexed by the right-hand side i (on which the
253 *> componentwise relative error depends), and the type of error
254 *> information as described below. There currently are up to three
255 *> pieces of information returned for each right-hand side. If
256 *> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
257 *> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS < 3, then at most
258 *> the first (:,N_ERR_BNDS) entries are returned.
259 *>
260 *> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
261 *> right-hand side.
262 *>
263 *> The second index in ERR_BNDS_COMP(:,err) contains the following
264 *> three fields:
265 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
266 *> reciprocal condition number is less than the threshold
267 *> sqrt(n) * slamch('Epsilon').
268 *>
269 *> err = 2 "Guaranteed" error bound: The estimated forward error,
270 *> almost certainly within a factor of 10 of the true error
271 *> so long as the next entry is greater than the threshold
272 *> sqrt(n) * slamch('Epsilon'). This error bound should only
273 *> be trusted if the previous boolean is true.
274 *>
275 *> err = 3 Reciprocal condition number: Estimated componentwise
276 *> reciprocal condition number. Compared with the threshold
277 *> sqrt(n) * slamch('Epsilon') to determine if the error
278 *> estimate is "guaranteed". These reciprocal condition
279 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
280 *> appropriately scaled matrix Z.
281 *> Let Z = S*(A*diag(x)), where x is the solution for the
282 *> current right-hand side and S scales each row of
283 *> A*diag(x) by a power of the radix so all absolute row
284 *> sums of Z are approximately 1.
285 *>
286 *> This subroutine is only responsible for setting the second field
287 *> above.
288 *> See Lapack Working Note 165 for further details and extra
289 *> cautions.
290 *> \endverbatim
291 *>
292 *> \param[in] RES
293 *> \verbatim
294 *> RES is COMPLEX array, dimension (N)
295 *> Workspace to hold the intermediate residual.
296 *> \endverbatim
297 *>
298 *> \param[in] AYB
299 *> \verbatim
300 *> AYB is REAL array, dimension (N)
301 *> Workspace.
302 *> \endverbatim
303 *>
304 *> \param[in] DY
305 *> \verbatim
306 *> DY is COMPLEX array, dimension (N)
307 *> Workspace to hold the intermediate solution.
308 *> \endverbatim
309 *>
310 *> \param[in] Y_TAIL
311 *> \verbatim
312 *> Y_TAIL is COMPLEX array, dimension (N)
313 *> Workspace to hold the trailing bits of the intermediate solution.
314 *> \endverbatim
315 *>
316 *> \param[in] RCOND
317 *> \verbatim
318 *> RCOND is REAL
319 *> Reciprocal scaled condition number. This is an estimate of the
320 *> reciprocal Skeel condition number of the matrix A after
321 *> equilibration (if done). If this is less than the machine
322 *> precision (in particular, if it is zero), the matrix is singular
323 *> to working precision. Note that the error may still be small even
324 *> if this number is very small and the matrix appears ill-
325 *> conditioned.
326 *> \endverbatim
327 *>
328 *> \param[in] ITHRESH
329 *> \verbatim
330 *> ITHRESH is INTEGER
331 *> The maximum number of residual computations allowed for
332 *> refinement. The default is 10. For 'aggressive' set to 100 to
333 *> permit convergence using approximate factorizations or
334 *> factorizations other than LU. If the factorization uses a
335 *> technique other than Gaussian elimination, the guarantees in
336 *> ERR_BNDS_NORM and ERR_BNDS_COMP may no longer be trustworthy.
337 *> \endverbatim
338 *>
339 *> \param[in] RTHRESH
340 *> \verbatim
341 *> RTHRESH is REAL
342 *> Determines when to stop refinement if the error estimate stops
343 *> decreasing. Refinement will stop when the next solution no longer
344 *> satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is
345 *> the infinity norm of Z. RTHRESH satisfies 0 < RTHRESH <= 1. The
346 *> default value is 0.5. For 'aggressive' set to 0.9 to permit
347 *> convergence on extremely ill-conditioned matrices. See LAWN 165
348 *> for more details.
349 *> \endverbatim
350 *>
351 *> \param[in] DZ_UB
352 *> \verbatim
353 *> DZ_UB is REAL
354 *> Determines when to start considering componentwise convergence.
355 *> Componentwise convergence is only considered after each component
356 *> of the solution Y is stable, which we definte as the relative
357 *> change in each component being less than DZ_UB. The default value
358 *> is 0.25, requiring the first bit to be stable. See LAWN 165 for
359 *> more details.
360 *> \endverbatim
361 *>
362 *> \param[in] IGNORE_CWISE
363 *> \verbatim
364 *> IGNORE_CWISE is LOGICAL
365 *> If .TRUE. then ignore componentwise convergence. Default value
366 *> is .FALSE..
367 *> \endverbatim
368 *>
369 *> \param[out] INFO
370 *> \verbatim
371 *> INFO is INTEGER
372 *> = 0: Successful exit.
373 *> < 0: if INFO = -i, the ith argument to CLA_SYRFSX_EXTENDED had an illegal
374 *> value
375 *> \endverbatim
376 *
377 * Authors:
378 * ========
379 *
380 *> \author Univ. of Tennessee
381 *> \author Univ. of California Berkeley
382 *> \author Univ. of Colorado Denver
383 *> \author NAG Ltd.
384 *
385 *> \date June 2017
386 *
387 *> \ingroup complexSYcomputational
388 *
389 * =====================================================================
390  SUBROUTINE cla_syrfsx_extended( PREC_TYPE, UPLO, N, NRHS, A, LDA,
391  $ AF, LDAF, IPIV, COLEQU, C, B, LDB,
392  $ Y, LDY, BERR_OUT, N_NORMS,
393  $ ERR_BNDS_NORM, ERR_BNDS_COMP, RES,
394  $ AYB, DY, Y_TAIL, RCOND, ITHRESH,
395  $ RTHRESH, DZ_UB, IGNORE_CWISE,
396  $ INFO )
397 *
398 * -- LAPACK computational routine (version 3.7.1) --
399 * -- LAPACK is a software package provided by Univ. of Tennessee, --
400 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
401 * June 2017
402 *
403 * .. Scalar Arguments ..
404  INTEGER INFO, LDA, LDAF, LDB, LDY, N, NRHS, PREC_TYPE,
405  $ N_NORMS, ITHRESH
406  CHARACTER UPLO
407  LOGICAL COLEQU, IGNORE_CWISE
408  REAL RTHRESH, DZ_UB
409 * ..
410 * .. Array Arguments ..
411  INTEGER IPIV( * )
412  COMPLEX A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
413  $ y( ldy, * ), res( * ), dy( * ), y_tail( * )
414  REAL C( * ), AYB( * ), RCOND, BERR_OUT( * ),
415  $ ERR_BNDS_NORM( NRHS, * ),
416  $ ERR_BNDS_COMP( NRHS, * )
417 * ..
418 *
419 * =====================================================================
420 *
421 * .. Local Scalars ..
422  INTEGER UPLO2, CNT, I, J, X_STATE, Z_STATE,
423  $ Y_PREC_STATE
424  REAL YK, DYK, YMIN, NORMY, NORMX, NORMDX, DXRAT,
425  $ DZRAT, PREVNORMDX, PREV_DZ_Z, DXRATMAX,
426  $ DZRATMAX, DX_X, DZ_Z, FINAL_DX_X, FINAL_DZ_Z,
427  $ EPS, HUGEVAL, INCR_THRESH
428  LOGICAL INCR_PREC, UPPER
429  COMPLEX ZDUM
430 * ..
431 * .. Parameters ..
432  INTEGER UNSTABLE_STATE, WORKING_STATE, CONV_STATE,
433  $ NOPROG_STATE, BASE_RESIDUAL, EXTRA_RESIDUAL,
434  $ extra_y
435  parameter( unstable_state = 0, working_state = 1,
436  $ conv_state = 2, noprog_state = 3 )
437  parameter( base_residual = 0, extra_residual = 1,
438  $ extra_y = 2 )
439  INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
440  INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
441  INTEGER CMP_ERR_I, PIV_GROWTH_I
442  PARAMETER ( FINAL_NRM_ERR_I = 1, final_cmp_err_i = 2,
443  $ berr_i = 3 )
444  parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
445  parameter( cmp_rcond_i = 7, cmp_err_i = 8,
446  $ piv_growth_i = 9 )
447  INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
448  $ la_linrx_cwise_i
449  parameter( la_linrx_itref_i = 1,
450  $ la_linrx_ithresh_i = 2 )
451  parameter( la_linrx_cwise_i = 3 )
452  INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
453  $ la_linrx_rcond_i
454  parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
455  parameter( la_linrx_rcond_i = 3 )
456 * ..
457 * .. External Functions ..
458  LOGICAL LSAME
459  EXTERNAL ILAUPLO
460  INTEGER ILAUPLO
461 * ..
462 * .. External Subroutines ..
463  EXTERNAL caxpy, ccopy, csytrs, csymv, blas_csymv_x,
464  $ blas_csymv2_x, cla_syamv, cla_wwaddw,
465  $ cla_lin_berr
466  REAL SLAMCH
467 * ..
468 * .. Intrinsic Functions ..
469  INTRINSIC abs, real, aimag, max, min
470 * ..
471 * .. Statement Functions ..
472  REAL CABS1
473 * ..
474 * .. Statement Function Definitions ..
475  cabs1( zdum ) = abs( real( zdum ) ) + abs( aimag( zdum ) )
476 * ..
477 * .. Executable Statements ..
478 *
479  info = 0
480  upper = lsame( uplo, 'U' )
481  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
482  info = -2
483  ELSE IF( n.LT.0 ) THEN
484  info = -3
485  ELSE IF( nrhs.LT.0 ) THEN
486  info = -4
487  ELSE IF( lda.LT.max( 1, n ) ) THEN
488  info = -6
489  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
490  info = -8
491  ELSE IF( ldb.LT.max( 1, n ) ) THEN
492  info = -13
493  ELSE IF( ldy.LT.max( 1, n ) ) THEN
494  info = -15
495  END IF
496  IF( info.NE.0 ) THEN
497  CALL xerbla( 'CLA_SYRFSX_EXTENDED', -info )
498  RETURN
499  END IF
500  eps = slamch( 'Epsilon' )
501  hugeval = slamch( 'Overflow' )
502 * Force HUGEVAL to Inf
503  hugeval = hugeval * hugeval
504 * Using HUGEVAL may lead to spurious underflows.
505  incr_thresh = real( n ) * eps
506 
507  IF ( lsame( uplo, 'L' ) ) THEN
508  uplo2 = ilauplo( 'L' )
509  ELSE
510  uplo2 = ilauplo( 'U' )
511  ENDIF
512 
513  DO j = 1, nrhs
514  y_prec_state = extra_residual
515  IF ( y_prec_state .EQ. extra_y ) THEN
516  DO i = 1, n
517  y_tail( i ) = 0.0
518  END DO
519  END IF
520 
521  dxrat = 0.0
522  dxratmax = 0.0
523  dzrat = 0.0
524  dzratmax = 0.0
525  final_dx_x = hugeval
526  final_dz_z = hugeval
527  prevnormdx = hugeval
528  prev_dz_z = hugeval
529  dz_z = hugeval
530  dx_x = hugeval
531 
532  x_state = working_state
533  z_state = unstable_state
534  incr_prec = .false.
535 
536  DO cnt = 1, ithresh
537 *
538 * Compute residual RES = B_s - op(A_s) * Y,
539 * op(A) = A, A**T, or A**H depending on TRANS (and type).
540 *
541  CALL ccopy( n, b( 1, j ), 1, res, 1 )
542  IF ( y_prec_state .EQ. base_residual ) THEN
543  CALL csymv( uplo, n, cmplx(-1.0), a, lda, y(1,j), 1,
544  $ cmplx(1.0), res, 1 )
545  ELSE IF ( y_prec_state .EQ. extra_residual ) THEN
546  CALL blas_csymv_x( uplo2, n, cmplx(-1.0), a, lda,
547  $ y( 1, j ), 1, cmplx(1.0), res, 1, prec_type )
548  ELSE
549  CALL blas_csymv2_x(uplo2, n, cmplx(-1.0), a, lda,
550  $ y(1, j), y_tail, 1, cmplx(1.0), res, 1, prec_type)
551  END IF
552 
553 ! XXX: RES is no longer needed.
554  CALL ccopy( n, res, 1, dy, 1 )
555  CALL csytrs( uplo, n, 1, af, ldaf, ipiv, dy, n, info )
556 *
557 * Calculate relative changes DX_X, DZ_Z and ratios DXRAT, DZRAT.
558 *
559  normx = 0.0
560  normy = 0.0
561  normdx = 0.0
562  dz_z = 0.0
563  ymin = hugeval
564 
565  DO i = 1, n
566  yk = cabs1( y( i, j ) )
567  dyk = cabs1( dy( i ) )
568 
569  IF ( yk .NE. 0.0 ) THEN
570  dz_z = max( dz_z, dyk / yk )
571  ELSE IF ( dyk .NE. 0.0 ) THEN
572  dz_z = hugeval
573  END IF
574 
575  ymin = min( ymin, yk )
576 
577  normy = max( normy, yk )
578 
579  IF ( colequ ) THEN
580  normx = max( normx, yk * c( i ) )
581  normdx = max( normdx, dyk * c( i ) )
582  ELSE
583  normx = normy
584  normdx = max( normdx, dyk )
585  END IF
586  END DO
587 
588  IF ( normx .NE. 0.0 ) THEN
589  dx_x = normdx / normx
590  ELSE IF ( normdx .EQ. 0.0 ) THEN
591  dx_x = 0.0
592  ELSE
593  dx_x = hugeval
594  END IF
595 
596  dxrat = normdx / prevnormdx
597  dzrat = dz_z / prev_dz_z
598 *
599 * Check termination criteria.
600 *
601  IF ( ymin*rcond .LT. incr_thresh*normy
602  $ .AND. y_prec_state .LT. extra_y )
603  $ incr_prec = .true.
604 
605  IF ( x_state .EQ. noprog_state .AND. dxrat .LE. rthresh )
606  $ x_state = working_state
607  IF ( x_state .EQ. working_state ) THEN
608  IF ( dx_x .LE. eps ) THEN
609  x_state = conv_state
610  ELSE IF ( dxrat .GT. rthresh ) THEN
611  IF ( y_prec_state .NE. extra_y ) THEN
612  incr_prec = .true.
613  ELSE
614  x_state = noprog_state
615  END IF
616  ELSE
617  IF (dxrat .GT. dxratmax) dxratmax = dxrat
618  END IF
619  IF ( x_state .GT. working_state ) final_dx_x = dx_x
620  END IF
621 
622  IF ( z_state .EQ. unstable_state .AND. dz_z .LE. dz_ub )
623  $ z_state = working_state
624  IF ( z_state .EQ. noprog_state .AND. dzrat .LE. rthresh )
625  $ z_state = working_state
626  IF ( z_state .EQ. working_state ) THEN
627  IF ( dz_z .LE. eps ) THEN
628  z_state = conv_state
629  ELSE IF ( dz_z .GT. dz_ub ) THEN
630  z_state = unstable_state
631  dzratmax = 0.0
632  final_dz_z = hugeval
633  ELSE IF ( dzrat .GT. rthresh ) THEN
634  IF ( y_prec_state .NE. extra_y ) THEN
635  incr_prec = .true.
636  ELSE
637  z_state = noprog_state
638  END IF
639  ELSE
640  IF ( dzrat .GT. dzratmax ) dzratmax = dzrat
641  END IF
642  IF ( z_state .GT. working_state ) final_dz_z = dz_z
643  END IF
644 
645  IF ( x_state.NE.working_state.AND.
646  $ ( ignore_cwise.OR.z_state.NE.working_state ) )
647  $ GOTO 666
648 
649  IF ( incr_prec ) THEN
650  incr_prec = .false.
651  y_prec_state = y_prec_state + 1
652  DO i = 1, n
653  y_tail( i ) = 0.0
654  END DO
655  END IF
656 
657  prevnormdx = normdx
658  prev_dz_z = dz_z
659 *
660 * Update soluton.
661 *
662  IF ( y_prec_state .LT. extra_y ) THEN
663  CALL caxpy( n, cmplx(1.0), dy, 1, y(1,j), 1 )
664  ELSE
665  CALL cla_wwaddw( n, y(1,j), y_tail, dy )
666  END IF
667 
668  END DO
669 * Target of "IF (Z_STOP .AND. X_STOP)". Sun's f77 won't EXIT.
670  666 CONTINUE
671 *
672 * Set final_* when cnt hits ithresh.
673 *
674  IF ( x_state .EQ. working_state ) final_dx_x = dx_x
675  IF ( z_state .EQ. working_state ) final_dz_z = dz_z
676 *
677 * Compute error bounds.
678 *
679  IF ( n_norms .GE. 1 ) THEN
680  err_bnds_norm( j, la_linrx_err_i ) =
681  $ final_dx_x / (1 - dxratmax)
682  END IF
683  IF ( n_norms .GE. 2 ) THEN
684  err_bnds_comp( j, la_linrx_err_i ) =
685  $ final_dz_z / (1 - dzratmax)
686  END IF
687 *
688 * Compute componentwise relative backward error from formula
689 * max(i) ( abs(R(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
690 * where abs(Z) is the componentwise absolute value of the matrix
691 * or vector Z.
692 *
693 * Compute residual RES = B_s - op(A_s) * Y,
694 * op(A) = A, A**T, or A**H depending on TRANS (and type).
695 *
696  CALL ccopy( n, b( 1, j ), 1, res, 1 )
697  CALL csymv( uplo, n, cmplx(-1.0), a, lda, y(1,j), 1,
698  $ cmplx(1.0), res, 1 )
699 
700  DO i = 1, n
701  ayb( i ) = cabs1( b( i, j ) )
702  END DO
703 *
704 * Compute abs(op(A_s))*abs(Y) + abs(B_s).
705 *
706  CALL cla_syamv ( uplo2, n, 1.0,
707  $ a, lda, y(1, j), 1, 1.0, ayb, 1 )
708 
709  CALL cla_lin_berr ( n, n, 1, res, ayb, berr_out( j ) )
710 *
711 * End of loop for each RHS.
712 *
713  END DO
714 *
715  RETURN
716  END
cla_lin_berr
subroutine cla_lin_berr(N, NZ, NRHS, RES, AYB, BERR)
CLA_LIN_BERR computes a component-wise relative backward error.
Definition: cla_lin_berr.f:103
cla_syamv
subroutine cla_syamv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CLA_SYAMV computes a matrix-vector product using a symmetric indefinite matrix to calculate error bou...
Definition: cla_syamv.f:181
csytrs
subroutine csytrs(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
CSYTRS
Definition: csytrs.f:122
cla_syrfsx_extended
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...
Definition: cla_syrfsx_extended.f:397
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
csymv
subroutine csymv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CSYMV computes a matrix-vector product for a complex symmetric matrix.
Definition: csymv.f:159
ccopy
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:83
caxpy
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:90
cla_wwaddw
subroutine cla_wwaddw(N, X, Y, W)
CLA_WWADDW adds a vector into a doubled-single vector.
Definition: cla_wwaddw.f:83