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