LAPACK  3.9.0
LAPACK: Linear Algebra PACKage
zla_syrfsx_extended.f
Go to the documentation of this file.
1 *> \brief \b ZLA_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 ZLA_SYRFSX_EXTENDED + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zla_syrfsx_extended.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zla_syrfsx_extended.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zla_syrfsx_extended.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZLA_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 * DOUBLE PRECISION RTHRESH, DZ_UB
35 * ..
36 * .. Array Arguments ..
37 * INTEGER IPIV( * )
38 * COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
39 * $ Y( LDY, * ), RES( * ), DY( * ), Y_TAIL( * )
40 * DOUBLE PRECISION 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 *> ZLA_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 ZSYRFSX 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*16 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*16 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 ZSYTRF.
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 ZSYTRF.
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 DOUBLE PRECISION 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*16 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*16 array, dimension (LDY,NRHS)
165 *> On entry, the solution matrix X, as computed by ZSYTRS.
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 DOUBLE PRECISION 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 ZLA_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 DOUBLE PRECISION 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 DOUBLE PRECISION 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*16 array, dimension (N)
295 *> Workspace to hold the intermediate residual.
296 *> \endverbatim
297 *>
298 *> \param[in] AYB
299 *> \verbatim
300 *> AYB is DOUBLE PRECISION array, dimension (N)
301 *> Workspace.
302 *> \endverbatim
303 *>
304 *> \param[in] DY
305 *> \verbatim
306 *> DY is COMPLEX*16 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*16 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 DOUBLE PRECISION
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 DOUBLE PRECISION
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 DOUBLE PRECISION
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 ZLA_HERFSX_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 complex16SYcomputational
388 *
389 * =====================================================================
390  SUBROUTINE zla_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  DOUBLE PRECISION RTHRESH, DZ_UB
409 * ..
410 * .. Array Arguments ..
411  INTEGER IPIV( * )
412  COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
413  $ y( ldy, * ), res( * ), dy( * ), y_tail( * )
414  DOUBLE PRECISION 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  DOUBLE PRECISION 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*16 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 zaxpy, zcopy, zsytrs, zsymv, blas_zsymv_x,
464  $ blas_zsymv2_x, zla_syamv, zla_wwaddw,
465  $ zla_lin_berr
466  DOUBLE PRECISION DLAMCH
467 * ..
468 * .. Intrinsic Functions ..
469  INTRINSIC abs, real, dimag, max, min
470 * ..
471 * .. Statement Functions ..
472  DOUBLE PRECISION CABS1
473 * ..
474 * .. Statement Function Definitions ..
475  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( 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( 'ZLA_HERFSX_EXTENDED', -info )
498  RETURN
499  END IF
500  eps = dlamch( 'Epsilon' )
501  hugeval = dlamch( 'Overflow' )
502 * Force HUGEVAL to Inf
503  hugeval = hugeval * hugeval
504 * Using HUGEVAL may lead to spurious underflows.
505  incr_thresh = dble( 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.0d+0
518  END DO
519  END IF
520 
521  dxrat = 0.0d+0
522  dxratmax = 0.0d+0
523  dzrat = 0.0d+0
524  dzratmax = 0.0d+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 zcopy( n, b( 1, j ), 1, res, 1 )
542  IF ( y_prec_state .EQ. base_residual ) THEN
543  CALL zsymv( uplo, n, dcmplx(-1.0d+0), a, lda, y(1,j), 1,
544  $ dcmplx(1.0d+0), res, 1 )
545  ELSE IF ( y_prec_state .EQ. extra_residual ) THEN
546  CALL blas_zsymv_x( uplo2, n, dcmplx(-1.0d+0), a, lda,
547  $ y( 1, j ), 1, dcmplx(1.0d+0), res, 1, prec_type )
548  ELSE
549  CALL blas_zsymv2_x(uplo2, n, dcmplx(-1.0d+0), a, lda,
550  $ y(1, j), y_tail, 1, dcmplx(1.0d+0), res, 1,
551  $ prec_type)
552  END IF
553 
554 ! XXX: RES is no longer needed.
555  CALL zcopy( n, res, 1, dy, 1 )
556  CALL zsytrs( uplo, n, 1, af, ldaf, ipiv, dy, n, info )
557 *
558 * Calculate relative changes DX_X, DZ_Z and ratios DXRAT, DZRAT.
559 *
560  normx = 0.0d+0
561  normy = 0.0d+0
562  normdx = 0.0d+0
563  dz_z = 0.0d+0
564  ymin = hugeval
565 
566  DO i = 1, n
567  yk = cabs1( y( i, j ) )
568  dyk = cabs1( dy( i ) )
569 
570  IF ( yk .NE. 0.0d+0 ) THEN
571  dz_z = max( dz_z, dyk / yk )
572  ELSE IF ( dyk .NE. 0.0d+0 ) THEN
573  dz_z = hugeval
574  END IF
575 
576  ymin = min( ymin, yk )
577 
578  normy = max( normy, yk )
579 
580  IF ( colequ ) THEN
581  normx = max( normx, yk * c( i ) )
582  normdx = max( normdx, dyk * c( i ) )
583  ELSE
584  normx = normy
585  normdx = max( normdx, dyk )
586  END IF
587  END DO
588 
589  IF ( normx .NE. 0.0d+0 ) THEN
590  dx_x = normdx / normx
591  ELSE IF ( normdx .EQ. 0.0d+0 ) THEN
592  dx_x = 0.0d+0
593  ELSE
594  dx_x = hugeval
595  END IF
596 
597  dxrat = normdx / prevnormdx
598  dzrat = dz_z / prev_dz_z
599 *
600 * Check termination criteria.
601 *
602  IF ( ymin*rcond .LT. incr_thresh*normy
603  $ .AND. y_prec_state .LT. extra_y )
604  $ incr_prec = .true.
605 
606  IF ( x_state .EQ. noprog_state .AND. dxrat .LE. rthresh )
607  $ x_state = working_state
608  IF ( x_state .EQ. working_state ) THEN
609  IF ( dx_x .LE. eps ) THEN
610  x_state = conv_state
611  ELSE IF ( dxrat .GT. rthresh ) THEN
612  IF ( y_prec_state .NE. extra_y ) THEN
613  incr_prec = .true.
614  ELSE
615  x_state = noprog_state
616  END IF
617  ELSE
618  IF (dxrat .GT. dxratmax) dxratmax = dxrat
619  END IF
620  IF ( x_state .GT. working_state ) final_dx_x = dx_x
621  END IF
622 
623  IF ( z_state .EQ. unstable_state .AND. dz_z .LE. dz_ub )
624  $ z_state = working_state
625  IF ( z_state .EQ. noprog_state .AND. dzrat .LE. rthresh )
626  $ z_state = working_state
627  IF ( z_state .EQ. working_state ) THEN
628  IF ( dz_z .LE. eps ) THEN
629  z_state = conv_state
630  ELSE IF ( dz_z .GT. dz_ub ) THEN
631  z_state = unstable_state
632  dzratmax = 0.0d+0
633  final_dz_z = hugeval
634  ELSE IF ( dzrat .GT. rthresh ) THEN
635  IF ( y_prec_state .NE. extra_y ) THEN
636  incr_prec = .true.
637  ELSE
638  z_state = noprog_state
639  END IF
640  ELSE
641  IF ( dzrat .GT. dzratmax ) dzratmax = dzrat
642  END IF
643  IF ( z_state .GT. working_state ) final_dz_z = dz_z
644  END IF
645 
646  IF ( x_state.NE.working_state.AND.
647  $ ( ignore_cwise.OR.z_state.NE.working_state ) )
648  $ GOTO 666
649 
650  IF ( incr_prec ) THEN
651  incr_prec = .false.
652  y_prec_state = y_prec_state + 1
653  DO i = 1, n
654  y_tail( i ) = 0.0d+0
655  END DO
656  END IF
657 
658  prevnormdx = normdx
659  prev_dz_z = dz_z
660 *
661 * Update soluton.
662 *
663  IF ( y_prec_state .LT. extra_y ) THEN
664  CALL zaxpy( n, dcmplx(1.0d+0), dy, 1, y(1,j), 1 )
665  ELSE
666  CALL zla_wwaddw( n, y(1,j), y_tail, dy )
667  END IF
668 
669  END DO
670 * Target of "IF (Z_STOP .AND. X_STOP)". Sun's f77 won't EXIT.
671  666 CONTINUE
672 *
673 * Set final_* when cnt hits ithresh.
674 *
675  IF ( x_state .EQ. working_state ) final_dx_x = dx_x
676  IF ( z_state .EQ. working_state ) final_dz_z = dz_z
677 *
678 * Compute error bounds.
679 *
680  IF ( n_norms .GE. 1 ) THEN
681  err_bnds_norm( j, la_linrx_err_i ) =
682  $ final_dx_x / (1 - dxratmax)
683  END IF
684  IF ( n_norms .GE. 2 ) THEN
685  err_bnds_comp( j, la_linrx_err_i ) =
686  $ final_dz_z / (1 - dzratmax)
687  END IF
688 *
689 * Compute componentwise relative backward error from formula
690 * max(i) ( abs(R(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
691 * where abs(Z) is the componentwise absolute value of the matrix
692 * or vector Z.
693 *
694 * Compute residual RES = B_s - op(A_s) * Y,
695 * op(A) = A, A**T, or A**H depending on TRANS (and type).
696 *
697  CALL zcopy( n, b( 1, j ), 1, res, 1 )
698  CALL zsymv( uplo, n, dcmplx(-1.0d+0), a, lda, y(1,j), 1,
699  $ dcmplx(1.0d+0), res, 1 )
700 
701  DO i = 1, n
702  ayb( i ) = cabs1( b( i, j ) )
703  END DO
704 *
705 * Compute abs(op(A_s))*abs(Y) + abs(B_s).
706 *
707  CALL zla_syamv ( uplo2, n, 1.0d+0,
708  $ a, lda, y(1, j), 1, 1.0d+0, ayb, 1 )
709 
710  CALL zla_lin_berr ( n, n, 1, res, ayb, berr_out( j ) )
711 *
712 * End of loop for each RHS.
713 *
714  END DO
715 *
716  RETURN
717  END
zla_lin_berr
subroutine zla_lin_berr(N, NZ, NRHS, RES, AYB, BERR)
ZLA_LIN_BERR computes a component-wise relative backward error.
Definition: zla_lin_berr.f:103
zaxpy
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:90
zla_wwaddw
subroutine zla_wwaddw(N, X, Y, W)
ZLA_WWADDW adds a vector into a doubled-single vector.
Definition: zla_wwaddw.f:83
zla_syrfsx_extended
subroutine zla_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)
ZLA_SYRFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric inde...
Definition: zla_syrfsx_extended.f:397
zsymv
subroutine zsymv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZSYMV computes a matrix-vector product for a complex symmetric matrix.
Definition: zsymv.f:159
zsytrs
subroutine zsytrs(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
ZSYTRS
Definition: zsytrs.f:122
zcopy
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:83
zla_syamv
subroutine zla_syamv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZLA_SYAMV computes a matrix-vector product using a symmetric indefinite matrix to calculate error bou...
Definition: zla_syamv.f:181
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62