LAPACK  3.9.0
LAPACK: Linear Algebra PACKage
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 P
68 *> = 'S': Single
69 *> = 'D': Double
70 *> = 'I': Indigenous
71 *> = 'X' or '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 T
79 *> = '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 (NRHS, N_ERR_BNDS)
199 *> For each right-hand side, this array contains information about
200 *> various error bounds and condition numbers corresponding to the
201 *> normwise relative error, which is defined as follows:
202 *>
203 *> Normwise relative error in the ith solution vector:
204 *> max_j (abs(XTRUE(j,i) - X(j,i)))
205 *> ------------------------------
206 *> max_j abs(X(j,i))
207 *>
208 *> The array is indexed by the type of error information as described
209 *> below. There currently are up to three pieces of information
210 *> returned.
211 *>
212 *> The first index in ERRS_N(i,:) corresponds to the ith
213 *> right-hand side.
214 *>
215 *> The second index in ERRS_N(:,err) contains the following
216 *> three fields:
217 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
218 *> reciprocal condition number is less than the threshold
219 *> sqrt(n) * slamch('Epsilon').
220 *>
221 *> err = 2 "Guaranteed" error bound: The estimated forward error,
222 *> almost certainly within a factor of 10 of the true error
223 *> so long as the next entry is greater than the threshold
224 *> sqrt(n) * slamch('Epsilon'). This error bound should only
225 *> be trusted if the previous boolean is true.
226 *>
227 *> err = 3 Reciprocal condition number: Estimated normwise
228 *> reciprocal condition number. Compared with the threshold
229 *> sqrt(n) * slamch('Epsilon') to determine if the error
230 *> estimate is "guaranteed". These reciprocal condition
231 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
232 *> appropriately scaled matrix Z.
233 *> Let Z = S*A, where S scales each row by a power of the
234 *> radix so all absolute row sums of Z are approximately 1.
235 *>
236 *> This subroutine is only responsible for setting the second field
237 *> above.
238 *> See Lapack Working Note 165 for further details and extra
239 *> cautions.
240 *> \endverbatim
241 *>
242 *> \param[in,out] ERRS_C
243 *> \verbatim
244 *> ERRS_C is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
245 *> For each right-hand side, this array contains information about
246 *> various error bounds and condition numbers corresponding to the
247 *> componentwise relative error, which is defined as follows:
248 *>
249 *> Componentwise relative error in the ith solution vector:
250 *> abs(XTRUE(j,i) - X(j,i))
251 *> max_j ----------------------
252 *> abs(X(j,i))
253 *>
254 *> The array is indexed by the right-hand side i (on which the
255 *> componentwise relative error depends), and the type of error
256 *> information as described below. There currently are up to three
257 *> pieces of information returned for each right-hand side. If
258 *> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
259 *> ERRS_C is not accessed. If N_ERR_BNDS < 3, then at most
260 *> the first (:,N_ERR_BNDS) entries are returned.
261 *>
262 *> The first index in ERRS_C(i,:) corresponds to the ith
263 *> right-hand side.
264 *>
265 *> The second index in ERRS_C(:,err) contains the following
266 *> three fields:
267 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
268 *> reciprocal condition number is less than the threshold
269 *> sqrt(n) * slamch('Epsilon').
270 *>
271 *> err = 2 "Guaranteed" error bound: The estimated forward error,
272 *> almost certainly within a factor of 10 of the true error
273 *> so long as the next entry is greater than the threshold
274 *> sqrt(n) * slamch('Epsilon'). This error bound should only
275 *> be trusted if the previous boolean is true.
276 *>
277 *> err = 3 Reciprocal condition number: Estimated componentwise
278 *> reciprocal condition number. Compared with the threshold
279 *> sqrt(n) * slamch('Epsilon') to determine if the error
280 *> estimate is "guaranteed". These reciprocal condition
281 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
282 *> appropriately scaled matrix Z.
283 *> Let Z = S*(A*diag(x)), where x is the solution for the
284 *> current right-hand side and S scales each row of
285 *> A*diag(x) by a power of the radix so all absolute row
286 *> sums of Z are approximately 1.
287 *>
288 *> This subroutine is only responsible for setting the second field
289 *> above.
290 *> See Lapack Working Note 165 for further details and extra
291 *> cautions.
292 *> \endverbatim
293 *>
294 *> \param[in] RES
295 *> \verbatim
296 *> RES is COMPLEX*16 array, dimension (N)
297 *> Workspace to hold the intermediate residual.
298 *> \endverbatim
299 *>
300 *> \param[in] AYB
301 *> \verbatim
302 *> AYB is DOUBLE PRECISION array, dimension (N)
303 *> Workspace.
304 *> \endverbatim
305 *>
306 *> \param[in] DY
307 *> \verbatim
308 *> DY is COMPLEX*16 array, dimension (N)
309 *> Workspace to hold the intermediate solution.
310 *> \endverbatim
311 *>
312 *> \param[in] Y_TAIL
313 *> \verbatim
314 *> Y_TAIL is COMPLEX*16 array, dimension (N)
315 *> Workspace to hold the trailing bits of the intermediate solution.
316 *> \endverbatim
317 *>
318 *> \param[in] RCOND
319 *> \verbatim
320 *> RCOND is DOUBLE PRECISION
321 *> Reciprocal scaled condition number. This is an estimate of the
322 *> reciprocal Skeel condition number of the matrix A after
323 *> equilibration (if done). If this is less than the machine
324 *> precision (in particular, if it is zero), the matrix is singular
325 *> to working precision. Note that the error may still be small even
326 *> if this number is very small and the matrix appears ill-
327 *> conditioned.
328 *> \endverbatim
329 *>
330 *> \param[in] ITHRESH
331 *> \verbatim
332 *> ITHRESH is INTEGER
333 *> The maximum number of residual computations allowed for
334 *> refinement. The default is 10. For 'aggressive' set to 100 to
335 *> permit convergence using approximate factorizations or
336 *> factorizations other than LU. If the factorization uses a
337 *> technique other than Gaussian elimination, the guarantees in
338 *> ERRS_N and ERRS_C may no longer be trustworthy.
339 *> \endverbatim
340 *>
341 *> \param[in] RTHRESH
342 *> \verbatim
343 *> RTHRESH is DOUBLE PRECISION
344 *> Determines when to stop refinement if the error estimate stops
345 *> decreasing. Refinement will stop when the next solution no longer
346 *> satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is
347 *> the infinity norm of Z. RTHRESH satisfies 0 < RTHRESH <= 1. The
348 *> default value is 0.5. For 'aggressive' set to 0.9 to permit
349 *> convergence on extremely ill-conditioned matrices. See LAWN 165
350 *> for more details.
351 *> \endverbatim
352 *>
353 *> \param[in] DZ_UB
354 *> \verbatim
355 *> DZ_UB is DOUBLE PRECISION
356 *> Determines when to start considering componentwise convergence.
357 *> Componentwise convergence is only considered after each component
358 *> of the solution Y is stable, which we definte as the relative
359 *> change in each component being less than DZ_UB. The default value
360 *> is 0.25, requiring the first bit to be stable. See LAWN 165 for
361 *> more details.
362 *> \endverbatim
363 *>
364 *> \param[in] IGNORE_CWISE
365 *> \verbatim
366 *> IGNORE_CWISE is LOGICAL
367 *> If .TRUE. then ignore componentwise convergence. Default value
368 *> is .FALSE..
369 *> \endverbatim
370 *>
371 *> \param[out] INFO
372 *> \verbatim
373 *> INFO is INTEGER
374 *> = 0: Successful exit.
375 *> < 0: if INFO = -i, the ith argument to ZGETRS had an illegal
376 *> value
377 *> \endverbatim
378 *
379 * Authors:
380 * ========
381 *
382 *> \author Univ. of Tennessee
383 *> \author Univ. of California Berkeley
384 *> \author Univ. of Colorado Denver
385 *> \author NAG Ltd.
386 *
387 *> \date June 2017
388 *
389 *> \ingroup complex16GEcomputational
390 *
391 * =====================================================================
392  SUBROUTINE zla_gerfsx_extended( PREC_TYPE, TRANS_TYPE, N, NRHS, A,
393  $ LDA, AF, LDAF, IPIV, COLEQU, C, B,
394  $ LDB, Y, LDY, BERR_OUT, N_NORMS,
395  $ ERRS_N, ERRS_C, RES, AYB, DY,
396  $ Y_TAIL, RCOND, ITHRESH, RTHRESH,
397  $ DZ_UB, IGNORE_CWISE, INFO )
398 *
399 * -- LAPACK computational routine (version 3.7.1) --
400 * -- LAPACK is a software package provided by Univ. of Tennessee, --
401 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
402 * June 2017
403 *
404 * .. Scalar Arguments ..
405  INTEGER INFO, LDA, LDAF, LDB, LDY, N, NRHS, PREC_TYPE,
406  $ TRANS_TYPE, N_NORMS
407  LOGICAL COLEQU, IGNORE_CWISE
408  INTEGER ITHRESH
409  DOUBLE PRECISION RTHRESH, DZ_UB
410 * ..
411 * .. Array Arguments
412  INTEGER IPIV( * )
413  COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
414  $ y( ldy, * ), res( * ), dy( * ), y_tail( * )
415  DOUBLE PRECISION C( * ), AYB( * ), RCOND, BERR_OUT( * ),
416  $ ERRS_N( NRHS, * ), ERRS_C( NRHS, * )
417 * ..
418 *
419 * =====================================================================
420 *
421 * .. Local Scalars ..
422  CHARACTER TRANS
423  INTEGER CNT, I, J, X_STATE, Z_STATE, 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
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,
437  $ noprog_state = 3 )
438  parameter( base_residual = 0, extra_residual = 1,
439  $ extra_y = 2 )
440  INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
441  INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
442  INTEGER CMP_ERR_I, PIV_GROWTH_I
443  parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
444  $ berr_i = 3 )
445  parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
446  parameter( cmp_rcond_i = 7, cmp_err_i = 8,
447  $ piv_growth_i = 9 )
448  INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
449  $ LA_LINRX_CWISE_I
450  parameter( la_linrx_itref_i = 1,
451  $ la_linrx_ithresh_i = 2 )
452  parameter( la_linrx_cwise_i = 3 )
453  INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
454  $ LA_LINRX_RCOND_I
455  parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
456  parameter( la_linrx_rcond_i = 3 )
457 * ..
458 * .. External Subroutines ..
459  EXTERNAL zaxpy, zcopy, zgetrs, zgemv, blas_zgemv_x,
460  $ blas_zgemv2_x, zla_geamv, zla_wwaddw, dlamch,
462  DOUBLE PRECISION DLAMCH
463  CHARACTER CHLA_TRANSTYPE
464 * ..
465 * .. Intrinsic Functions ..
466  INTRINSIC abs, max, min
467 * ..
468 * .. Statement Functions ..
469  DOUBLE PRECISION CABS1
470 * ..
471 * .. Statement Function Definitions ..
472  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
473 * ..
474 * .. Executable Statements ..
475 *
476  IF ( info.NE.0 ) RETURN
477  trans = chla_transtype(trans_type)
478  eps = dlamch( 'Epsilon' )
479  hugeval = dlamch( 'Overflow' )
480 * Force HUGEVAL to Inf
481  hugeval = hugeval * hugeval
482 * Using HUGEVAL may lead to spurious underflows.
483  incr_thresh = dble( n ) * eps
484 *
485  DO j = 1, nrhs
486  y_prec_state = extra_residual
487  IF ( y_prec_state .EQ. extra_y ) THEN
488  DO i = 1, n
489  y_tail( i ) = 0.0d+0
490  END DO
491  END IF
492 
493  dxrat = 0.0d+0
494  dxratmax = 0.0d+0
495  dzrat = 0.0d+0
496  dzratmax = 0.0d+0
497  final_dx_x = hugeval
498  final_dz_z = hugeval
499  prevnormdx = hugeval
500  prev_dz_z = hugeval
501  dz_z = hugeval
502  dx_x = hugeval
503 
504  x_state = working_state
505  z_state = unstable_state
506  incr_prec = .false.
507 
508  DO cnt = 1, ithresh
509 *
510 * Compute residual RES = B_s - op(A_s) * Y,
511 * op(A) = A, A**T, or A**H depending on TRANS (and type).
512 *
513  CALL zcopy( n, b( 1, j ), 1, res, 1 )
514  IF ( y_prec_state .EQ. base_residual ) THEN
515  CALL zgemv( trans, n, n, (-1.0d+0,0.0d+0), a, lda,
516  $ y( 1, j ), 1, (1.0d+0,0.0d+0), res, 1)
517  ELSE IF (y_prec_state .EQ. extra_residual) THEN
518  CALL blas_zgemv_x( trans_type, n, n, (-1.0d+0,0.0d+0), a,
519  $ lda, y( 1, j ), 1, (1.0d+0,0.0d+0),
520  $ res, 1, prec_type )
521  ELSE
522  CALL blas_zgemv2_x( trans_type, n, n, (-1.0d+0,0.0d+0),
523  $ a, lda, y(1, j), y_tail, 1, (1.0d+0,0.0d+0), res, 1,
524  $ prec_type)
525  END IF
526 
527 ! XXX: RES is no longer needed.
528  CALL zcopy( n, res, 1, dy, 1 )
529  CALL zgetrs( trans, n, 1, af, ldaf, ipiv, dy, n, info )
530 *
531 * Calculate relative changes DX_X, DZ_Z and ratios DXRAT, DZRAT.
532 *
533  normx = 0.0d+0
534  normy = 0.0d+0
535  normdx = 0.0d+0
536  dz_z = 0.0d+0
537  ymin = hugeval
538 *
539  DO i = 1, n
540  yk = cabs1( y( i, j ) )
541  dyk = cabs1( dy( i ) )
542 
543  IF ( yk .NE. 0.0d+0 ) THEN
544  dz_z = max( dz_z, dyk / yk )
545  ELSE IF ( dyk .NE. 0.0d+0 ) THEN
546  dz_z = hugeval
547  END IF
548 
549  ymin = min( ymin, yk )
550 
551  normy = max( normy, yk )
552 
553  IF ( colequ ) THEN
554  normx = max( normx, yk * c( i ) )
555  normdx = max( normdx, dyk * c( i ) )
556  ELSE
557  normx = normy
558  normdx = max(normdx, dyk)
559  END IF
560  END DO
561 
562  IF ( normx .NE. 0.0d+0 ) THEN
563  dx_x = normdx / normx
564  ELSE IF ( normdx .EQ. 0.0d+0 ) THEN
565  dx_x = 0.0d+0
566  ELSE
567  dx_x = hugeval
568  END IF
569 
570  dxrat = normdx / prevnormdx
571  dzrat = dz_z / prev_dz_z
572 *
573 * Check termination criteria
574 *
575  IF (.NOT.ignore_cwise
576  $ .AND. ymin*rcond .LT. incr_thresh*normy
577  $ .AND. y_prec_state .LT. extra_y )
578  $ incr_prec = .true.
579 
580  IF ( x_state .EQ. noprog_state .AND. dxrat .LE. rthresh )
581  $ x_state = working_state
582  IF ( x_state .EQ. working_state ) THEN
583  IF (dx_x .LE. eps) THEN
584  x_state = conv_state
585  ELSE IF ( dxrat .GT. rthresh ) THEN
586  IF ( y_prec_state .NE. extra_y ) THEN
587  incr_prec = .true.
588  ELSE
589  x_state = noprog_state
590  END IF
591  ELSE
592  IF ( dxrat .GT. dxratmax ) dxratmax = dxrat
593  END IF
594  IF ( x_state .GT. working_state ) final_dx_x = dx_x
595  END IF
596 
597  IF ( z_state .EQ. unstable_state .AND. dz_z .LE. dz_ub )
598  $ z_state = working_state
599  IF ( z_state .EQ. noprog_state .AND. dzrat .LE. rthresh )
600  $ z_state = working_state
601  IF ( z_state .EQ. working_state ) THEN
602  IF ( dz_z .LE. eps ) THEN
603  z_state = conv_state
604  ELSE IF ( dz_z .GT. dz_ub ) THEN
605  z_state = unstable_state
606  dzratmax = 0.0d+0
607  final_dz_z = hugeval
608  ELSE IF ( dzrat .GT. rthresh ) THEN
609  IF ( y_prec_state .NE. extra_y ) THEN
610  incr_prec = .true.
611  ELSE
612  z_state = noprog_state
613  END IF
614  ELSE
615  IF ( dzrat .GT. dzratmax ) dzratmax = dzrat
616  END IF
617  IF ( z_state .GT. working_state ) final_dz_z = dz_z
618  END IF
619 *
620 * Exit if both normwise and componentwise stopped working,
621 * but if componentwise is unstable, let it go at least two
622 * iterations.
623 *
624  IF ( x_state.NE.working_state ) THEN
625  IF ( ignore_cwise ) GOTO 666
626  IF ( z_state.EQ.noprog_state .OR. z_state.EQ.conv_state )
627  $ GOTO 666
628  IF ( z_state.EQ.unstable_state .AND. cnt.GT.1 ) GOTO 666
629  END IF
630 
631  IF ( incr_prec ) THEN
632  incr_prec = .false.
633  y_prec_state = y_prec_state + 1
634  DO i = 1, n
635  y_tail( i ) = 0.0d+0
636  END DO
637  END IF
638 
639  prevnormdx = normdx
640  prev_dz_z = dz_z
641 *
642 * Update soluton.
643 *
644  IF ( y_prec_state .LT. extra_y ) THEN
645  CALL zaxpy( n, (1.0d+0,0.0d+0), dy, 1, y(1,j), 1 )
646  ELSE
647  CALL zla_wwaddw( n, y( 1, j ), y_tail, dy )
648  END IF
649 
650  END DO
651 * Target of "IF (Z_STOP .AND. X_STOP)". Sun's f77 won't EXIT.
652  666 CONTINUE
653 *
654 * Set final_* when cnt hits ithresh
655 *
656  IF ( x_state .EQ. working_state ) final_dx_x = dx_x
657  IF ( z_state .EQ. working_state ) final_dz_z = dz_z
658 *
659 * Compute error bounds
660 *
661  IF (n_norms .GE. 1) THEN
662  errs_n( j, la_linrx_err_i ) = final_dx_x / (1 - dxratmax)
663 
664  END IF
665  IF ( n_norms .GE. 2 ) THEN
666  errs_c( j, la_linrx_err_i ) = final_dz_z / (1 - dzratmax)
667  END IF
668 *
669 * Compute componentwise relative backward error from formula
670 * max(i) ( abs(R(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
671 * where abs(Z) is the componentwise absolute value of the matrix
672 * or vector Z.
673 *
674 * Compute residual RES = B_s - op(A_s) * Y,
675 * op(A) = A, A**T, or A**H depending on TRANS (and type).
676 *
677  CALL zcopy( n, b( 1, j ), 1, res, 1 )
678  CALL zgemv( trans, n, n, (-1.0d+0,0.0d+0), a, lda, y(1,j), 1,
679  $ (1.0d+0,0.0d+0), res, 1 )
680 
681  DO i = 1, n
682  ayb( i ) = cabs1( b( i, j ) )
683  END DO
684 *
685 * Compute abs(op(A_s))*abs(Y) + abs(B_s).
686 *
687  CALL zla_geamv ( trans_type, n, n, 1.0d+0,
688  $ a, lda, y(1, j), 1, 1.0d+0, ayb, 1 )
689 
690  CALL zla_lin_berr ( n, n, 1, res, ayb, berr_out( j ) )
691 *
692 * End of loop for each RHS.
693 *
694  END DO
695 *
696  RETURN
697  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
zgemv
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
zla_wwaddw
subroutine zla_wwaddw(N, X, Y, W)
ZLA_WWADDW adds a vector into a doubled-single vector.
Definition: zla_wwaddw.f:83
zcopy
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:83
zla_geamv
subroutine zla_geamv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZLA_GEAMV computes a matrix-vector product using a general matrix to calculate error bounds.
Definition: zla_geamv.f:177
zla_gerfsx_extended
subroutine zla_gerfsx_extended(PREC_TYPE, TRANS_TYPE, N, NRHS, A, LDA, AF, LDAF, IPIV, COLEQU, C, B, LDB, Y, LDY, BERR_OUT, N_NORMS, ERRS_N, ERRS_C, RES, AYB, DY, Y_TAIL, RCOND, ITHRESH, RTHRESH, DZ_UB, IGNORE_CWISE, INFO)
ZLA_GERFSX_EXTENDED
Definition: zla_gerfsx_extended.f:398
zgetrs
subroutine zgetrs(TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
ZGETRS
Definition: zgetrs.f:123
dlamch
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:70
chla_transtype
character *1 function chla_transtype(TRANS)
CHLA_TRANSTYPE
Definition: chla_transtype.f:60