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