LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zherfsx()

subroutine zherfsx ( character uplo,
character equed,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldaf, * ) af,
integer ldaf,
integer, dimension( * ) ipiv,
double precision, dimension( * ) s,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( ldx, * ) x,
integer ldx,
double precision rcond,
double precision, dimension( * ) berr,
integer n_err_bnds,
double precision, dimension( nrhs, * ) err_bnds_norm,
double precision, dimension( nrhs, * ) err_bnds_comp,
integer nparams,
double precision, dimension( * ) params,
complex*16, dimension( * ) work,
double precision, dimension( * ) rwork,
integer info )

ZHERFSX

Download ZHERFSX + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!>    ZHERFSX improves the computed solution to a system of linear
!>    equations when the coefficient matrix is Hermitian indefinite, and
!>    provides error bounds and backward error estimates for the
!>    solution.  In addition to normwise error bound, the code provides
!>    maximum componentwise error bound if possible.  See comments for
!>    ERR_BNDS_NORM and ERR_BNDS_COMP for details of the error bounds.
!>
!>    The original system of linear equations may have been equilibrated
!>    before calling this routine, as described by arguments EQUED and S
!>    below. In this case, the solution and error bounds returned are
!>    for the original unequilibrated system.
!> 
!>     Some optional parameters are bundled in the PARAMS array.  These
!>     settings determine how refinement is performed, but often the
!>     defaults are acceptable.  If the defaults are acceptable, users
!>     can pass NPARAMS = 0 which prevents the source code from accessing
!>     the PARAMS argument.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>       = 'U':  Upper triangle of A is stored;
!>       = 'L':  Lower triangle of A is stored.
!> 
[in]EQUED
!>          EQUED is CHARACTER*1
!>     Specifies the form of equilibration that was done to A
!>     before calling this routine. This is needed to compute
!>     the solution and error bounds correctly.
!>       = 'N':  No equilibration
!>       = 'Y':  Both row and column equilibration, i.e., A has been
!>               replaced by diag(S) * A * diag(S).
!>               The right hand side B has been changed accordingly.
!> 
[in]N
!>          N is INTEGER
!>     The order of the matrix A.  N >= 0.
!> 
[in]NRHS
!>          NRHS is INTEGER
!>     The number of right hand sides, i.e., the number of columns
!>     of the matrices B and X.  NRHS >= 0.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension (LDA,N)
!>     The Hermitian matrix A.  If UPLO = 'U', the leading N-by-N
!>     upper triangular part of A contains the upper triangular
!>     part of the matrix A, and the strictly lower triangular
!>     part of A is not referenced.  If UPLO = 'L', the leading
!>     N-by-N lower triangular part of A contains the lower
!>     triangular part of the matrix A, and the strictly upper
!>     triangular part of A is not referenced.
!> 
[in]LDA
!>          LDA is INTEGER
!>     The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in]AF
!>          AF is COMPLEX*16 array, dimension (LDAF,N)
!>     The factored form of the matrix A.  AF contains the block
!>     diagonal matrix D and the multipliers used to obtain the
!>     factor U or L from the factorization A = U*D*U**H or A =
!>     L*D*L**H as computed by ZHETRF.
!> 
[in]LDAF
!>          LDAF is INTEGER
!>     The leading dimension of the array AF.  LDAF >= max(1,N).
!> 
[in]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>     Details of the interchanges and the block structure of D
!>     as determined by ZHETRF.
!> 
[in,out]S
!>          S is DOUBLE PRECISION array, dimension (N)
!>     The scale factors for A.  If EQUED = 'Y', A is multiplied on
!>     the left and right by diag(S).  S is an input argument if FACT =
!>     'F'; otherwise, S is an output argument.  If FACT = 'F' and EQUED
!>     = 'Y', each element of S must be positive.  If S is output, each
!>     element of S is a power of the radix. If S is input, each element
!>     of S should be a power of the radix to ensure a reliable solution
!>     and error estimates. Scaling by powers of the radix does not cause
!>     rounding errors unless the result underflows or overflows.
!>     Rounding errors during scaling lead to refining with a matrix that
!>     is not equivalent to the input matrix, producing error estimates
!>     that may not be reliable.
!> 
[in]B
!>          B is COMPLEX*16 array, dimension (LDB,NRHS)
!>     The right hand side matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>     The leading dimension of the array B.  LDB >= max(1,N).
!> 
[in,out]X
!>          X is COMPLEX*16 array, dimension (LDX,NRHS)
!>     On entry, the solution matrix X, as computed by ZHETRS.
!>     On exit, the improved solution matrix X.
!> 
[in]LDX
!>          LDX is INTEGER
!>     The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]RCOND
!>          RCOND is DOUBLE PRECISION
!>     Reciprocal scaled condition number.  This is an estimate of the
!>     reciprocal Skeel condition number of the matrix A after
!>     equilibration (if done).  If this is less than the machine
!>     precision (in particular, if it is zero), the matrix is singular
!>     to working precision.  Note that the error may still be small even
!>     if this number is very small and the matrix appears ill-
!>     conditioned.
!> 
[out]BERR
!>          BERR is DOUBLE PRECISION array, dimension (NRHS)
!>     Componentwise relative backward error.  This is the
!>     componentwise relative backward error of each solution vector X(j)
!>     (i.e., the smallest relative change in any element of A or B that
!>     makes X(j) an exact solution).
!> 
[in]N_ERR_BNDS
!>          N_ERR_BNDS is INTEGER
!>     Number of error bounds to return for each right hand side
!>     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
!>     ERR_BNDS_COMP below.
!> 
[out]ERR_BNDS_NORM
!>          ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
!>     For each right-hand side, this array contains information about
!>     various error bounds and condition numbers corresponding to the
!>     normwise relative error, which is defined as follows:
!>
!>     Normwise relative error in the ith solution vector:
!>             max_j (abs(XTRUE(j,i) - X(j,i)))
!>            ------------------------------
!>                  max_j abs(X(j,i))
!>
!>     The array is indexed by the type of error information as described
!>     below. There currently are up to three pieces of information
!>     returned.
!>
!>     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
!>     right-hand side.
!>
!>     The second index in ERR_BNDS_NORM(:,err) contains the following
!>     three fields:
!>     err = 1  boolean. Trust the answer if the
!>              reciprocal condition number is less than the threshold
!>              sqrt(n) * dlamch('Epsilon').
!>
!>     err = 2  error bound: The estimated forward error,
!>              almost certainly within a factor of 10 of the true error
!>              so long as the next entry is greater than the threshold
!>              sqrt(n) * dlamch('Epsilon'). This error bound should only
!>              be trusted if the previous boolean is true.
!>
!>     err = 3  Reciprocal condition number: Estimated normwise
!>              reciprocal condition number.  Compared with the threshold
!>              sqrt(n) * dlamch('Epsilon') to determine if the error
!>              estimate is . These reciprocal condition
!>              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
!>              appropriately scaled matrix Z.
!>              Let Z = S*A, where S scales each row by a power of the
!>              radix so all absolute row sums of Z are approximately 1.
!>
!>     See Lapack Working Note 165 for further details and extra
!>     cautions.
!> 
[out]ERR_BNDS_COMP
!>          ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
!>     For each right-hand side, this array contains information about
!>     various error bounds and condition numbers corresponding to the
!>     componentwise relative error, which is defined as follows:
!>
!>     Componentwise relative error in the ith solution vector:
!>                    abs(XTRUE(j,i) - X(j,i))
!>             max_j ----------------------
!>                         abs(X(j,i))
!>
!>     The array is indexed by the right-hand side i (on which the
!>     componentwise relative error depends), and the type of error
!>     information as described below. There currently are up to three
!>     pieces of information returned for each right-hand side. If
!>     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
!>     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS < 3, then at most
!>     the first (:,N_ERR_BNDS) entries are returned.
!>
!>     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
!>     right-hand side.
!>
!>     The second index in ERR_BNDS_COMP(:,err) contains the following
!>     three fields:
!>     err = 1  boolean. Trust the answer if the
!>              reciprocal condition number is less than the threshold
!>              sqrt(n) * dlamch('Epsilon').
!>
!>     err = 2  error bound: The estimated forward error,
!>              almost certainly within a factor of 10 of the true error
!>              so long as the next entry is greater than the threshold
!>              sqrt(n) * dlamch('Epsilon'). This error bound should only
!>              be trusted if the previous boolean is true.
!>
!>     err = 3  Reciprocal condition number: Estimated componentwise
!>              reciprocal condition number.  Compared with the threshold
!>              sqrt(n) * dlamch('Epsilon') to determine if the error
!>              estimate is . These reciprocal condition
!>              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
!>              appropriately scaled matrix Z.
!>              Let Z = S*(A*diag(x)), where x is the solution for the
!>              current right-hand side and S scales each row of
!>              A*diag(x) by a power of the radix so all absolute row
!>              sums of Z are approximately 1.
!>
!>     See Lapack Working Note 165 for further details and extra
!>     cautions.
!> 
[in]NPARAMS
!>          NPARAMS is INTEGER
!>     Specifies the number of parameters set in PARAMS.  If <= 0, the
!>     PARAMS array is never referenced and default values are used.
!> 
[in,out]PARAMS
!>          PARAMS is DOUBLE PRECISION array, dimension NPARAMS
!>     Specifies algorithm parameters.  If an entry is < 0.0, then
!>     that entry will be filled with default value used for that
!>     parameter.  Only positions up to NPARAMS are accessed; defaults
!>     are used for higher-numbered parameters.
!>
!>       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
!>            refinement or not.
!>         Default: 1.0D+0
!>            = 0.0:  No refinement is performed, and no error bounds are
!>                    computed.
!>            = 1.0:  Use the double-precision refinement algorithm,
!>                    possibly with doubled-single computations if the
!>                    compilation environment does not support DOUBLE
!>                    PRECISION.
!>              (other values are reserved for future use)
!>
!>       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
!>            computations allowed for refinement.
!>         Default: 10
!>         Aggressive: Set to 100 to permit convergence using approximate
!>                     factorizations or factorizations other than LU. If
!>                     the factorization uses a technique other than
!>                     Gaussian elimination, the guarantees in
!>                     err_bnds_norm and err_bnds_comp may no longer be
!>                     trustworthy.
!>
!>       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
!>            will attempt to find a solution with small componentwise
!>            relative error in the double-precision algorithm.  Positive
!>            is true, 0.0 is false.
!>         Default: 1.0 (attempt componentwise convergence)
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (2*N)
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (2*N)
!> 
[out]INFO
!>          INFO is INTEGER
!>       = 0:  Successful exit. The solution to every right-hand side is
!>         guaranteed.
!>       < 0:  If INFO = -i, the i-th argument had an illegal value
!>       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
!>         has been completed, but the factor U is exactly singular, so
!>         the solution and error bounds could not be computed. RCOND = 0
!>         is returned.
!>       = N+J: The solution corresponding to the Jth right-hand side is
!>         not guaranteed. The solutions corresponding to other right-
!>         hand sides K with K > J may not be guaranteed as well, but
!>         only the first such right-hand side is reported. If a small
!>         componentwise error is not requested (PARAMS(3) = 0.0) then
!>         the Jth right-hand side is the first with a normwise error
!>         bound that is not guaranteed (the smallest J such
!>         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
!>         the Jth right-hand side is the first with either a normwise or
!>         componentwise error bound that is not guaranteed (the smallest
!>         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
!>         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
!>         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
!>         about all of the right-hand sides check ERR_BNDS_NORM or
!>         ERR_BNDS_COMP.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 395 of file zherfsx.f.

400*
401* -- LAPACK computational routine --
402* -- LAPACK is a software package provided by Univ. of Tennessee, --
403* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
404*
405* .. Scalar Arguments ..
406 CHARACTER UPLO, EQUED
407 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
408 $ N_ERR_BNDS
409 DOUBLE PRECISION RCOND
410* ..
411* .. Array Arguments ..
412 INTEGER IPIV( * )
413 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
414 $ X( LDX, * ), WORK( * )
415 DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
416 $ ERR_BNDS_NORM( NRHS, * ),
417 $ ERR_BNDS_COMP( NRHS, * )
418*
419* ==================================================================
420*
421* .. Parameters ..
422 DOUBLE PRECISION ZERO, ONE
423 parameter( zero = 0.0d+0, one = 1.0d+0 )
424 DOUBLE PRECISION ITREF_DEFAULT, ITHRESH_DEFAULT
425 DOUBLE PRECISION COMPONENTWISE_DEFAULT, RTHRESH_DEFAULT
426 DOUBLE PRECISION DZTHRESH_DEFAULT
427 parameter( itref_default = 1.0d+0 )
428 parameter( ithresh_default = 10.0d+0 )
429 parameter( componentwise_default = 1.0d+0 )
430 parameter( rthresh_default = 0.5d+0 )
431 parameter( dzthresh_default = 0.25d+0 )
432 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
433 $ LA_LINRX_CWISE_I
434 parameter( la_linrx_itref_i = 1,
435 $ la_linrx_ithresh_i = 2 )
436 parameter( la_linrx_cwise_i = 3 )
437 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
438 $ LA_LINRX_RCOND_I
439 parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
440 parameter( la_linrx_rcond_i = 3 )
441* ..
442* .. Local Scalars ..
443 CHARACTER(1) NORM
444 LOGICAL RCEQU
445 INTEGER J, PREC_TYPE, REF_TYPE
446 INTEGER N_NORMS
447 DOUBLE PRECISION ANORM, RCOND_TMP
448 DOUBLE PRECISION ILLRCOND_THRESH, ERR_LBND, CWISE_WRONG
449 LOGICAL IGNORE_CWISE
450 INTEGER ITHRESH
451 DOUBLE PRECISION RTHRESH, UNSTABLE_THRESH
452* ..
453* .. External Subroutines ..
455* ..
456* .. Intrinsic Functions ..
457 INTRINSIC max, sqrt, transfer
458* ..
459* .. External Functions ..
460 EXTERNAL lsame, ilaprec
461 EXTERNAL dlamch, zlanhe, zla_hercond_x,
463 DOUBLE PRECISION DLAMCH, ZLANHE, ZLA_HERCOND_X, ZLA_HERCOND_C
464 LOGICAL LSAME
465 INTEGER ILAPREC
466* ..
467* .. Executable Statements ..
468*
469* Check the input parameters.
470*
471 info = 0
472 ref_type = int( itref_default )
473 IF ( nparams .GE. la_linrx_itref_i ) THEN
474 IF ( params( la_linrx_itref_i ) .LT. 0.0d+0 ) THEN
475 params( la_linrx_itref_i ) = itref_default
476 ELSE
477 ref_type = params( la_linrx_itref_i )
478 END IF
479 END IF
480*
481* Set default parameters.
482*
483 illrcond_thresh = dble( n ) * dlamch( 'Epsilon' )
484 ithresh = int( ithresh_default )
485 rthresh = rthresh_default
486 unstable_thresh = dzthresh_default
487 ignore_cwise = componentwise_default .EQ. 0.0d+0
488*
489 IF ( nparams.GE.la_linrx_ithresh_i ) THEN
490 IF ( params( la_linrx_ithresh_i ).LT.0.0d+0 ) THEN
491 params( la_linrx_ithresh_i ) = ithresh
492 ELSE
493 ithresh = int( params( la_linrx_ithresh_i ) )
494 END IF
495 END IF
496 IF ( nparams.GE.la_linrx_cwise_i ) THEN
497 IF ( params(la_linrx_cwise_i ).LT.0.0d+0 ) THEN
498 IF ( ignore_cwise ) THEN
499 params( la_linrx_cwise_i ) = 0.0d+0
500 ELSE
501 params( la_linrx_cwise_i ) = 1.0d+0
502 END IF
503 ELSE
504 ignore_cwise = params( la_linrx_cwise_i ) .EQ. 0.0d+0
505 END IF
506 END IF
507 IF ( ref_type .EQ. 0 .OR. n_err_bnds .EQ. 0 ) THEN
508 n_norms = 0
509 ELSE IF ( ignore_cwise ) THEN
510 n_norms = 1
511 ELSE
512 n_norms = 2
513 END IF
514*
515 rcequ = lsame( equed, 'Y' )
516*
517* Test input parameters.
518*
519 IF (.NOT.lsame( uplo, 'U' ) .AND.
520 $ .NOT.lsame( uplo, 'L' ) ) THEN
521 info = -1
522 ELSE IF( .NOT.rcequ .AND. .NOT.lsame( equed, 'N' ) ) THEN
523 info = -2
524 ELSE IF( n.LT.0 ) THEN
525 info = -3
526 ELSE IF( nrhs.LT.0 ) THEN
527 info = -4
528 ELSE IF( lda.LT.max( 1, n ) ) THEN
529 info = -6
530 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
531 info = -8
532 ELSE IF( ldb.LT.max( 1, n ) ) THEN
533 info = -12
534 ELSE IF( ldx.LT.max( 1, n ) ) THEN
535 info = -14
536 END IF
537 IF( info.NE.0 ) THEN
538 CALL xerbla( 'ZHERFSX', -info )
539 RETURN
540 END IF
541*
542* Quick return if possible.
543*
544 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
545 rcond = 1.0d+0
546 DO j = 1, nrhs
547 berr( j ) = 0.0d+0
548 IF ( n_err_bnds .GE. 1 ) THEN
549 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
550 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
551 END IF
552 IF ( n_err_bnds .GE. 2 ) THEN
553 err_bnds_norm( j, la_linrx_err_i ) = 0.0d+0
554 err_bnds_comp( j, la_linrx_err_i ) = 0.0d+0
555 END IF
556 IF ( n_err_bnds .GE. 3 ) THEN
557 err_bnds_norm( j, la_linrx_rcond_i ) = 1.0d+0
558 err_bnds_comp( j, la_linrx_rcond_i ) = 1.0d+0
559 END IF
560 END DO
561 RETURN
562 END IF
563*
564* Default to failure.
565*
566 rcond = 0.0d+0
567 DO j = 1, nrhs
568 berr( j ) = 1.0d+0
569 IF ( n_err_bnds .GE. 1 ) THEN
570 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
571 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
572 END IF
573 IF ( n_err_bnds .GE. 2 ) THEN
574 err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
575 err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
576 END IF
577 IF ( n_err_bnds .GE. 3 ) THEN
578 err_bnds_norm( j, la_linrx_rcond_i ) = 0.0d+0
579 err_bnds_comp( j, la_linrx_rcond_i ) = 0.0d+0
580 END IF
581 END DO
582*
583* Compute the norm of A and the reciprocal of the condition
584* number of A.
585*
586 norm = 'I'
587 anorm = zlanhe( norm, uplo, n, a, lda, rwork )
588 CALL zhecon( uplo, n, af, ldaf, ipiv, anorm, rcond, work,
589 $ info )
590*
591* Perform refinement on each right-hand side
592*
593 IF ( ref_type .NE. 0 ) THEN
594
595 prec_type = ilaprec( 'E' )
596
597 CALL zla_herfsx_extended( prec_type, uplo, n,
598 $ nrhs, a, lda, af, ldaf, ipiv, rcequ, s, b,
599 $ ldb, x, ldx, berr, n_norms, err_bnds_norm, err_bnds_comp,
600 $ work, rwork, work(n+1),
601 $ transfer(rwork(1:2*n), (/ (zero, zero) /), n), rcond,
602 $ ithresh, rthresh, unstable_thresh, ignore_cwise,
603 $ info )
604 END IF
605
606 err_lbnd = max( 10.0d+0,
607 $ sqrt( dble( n ) ) ) * dlamch( 'Epsilon' )
608 IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 1 ) THEN
609*
610* Compute scaled normwise condition number cond(A*C).
611*
612 IF ( rcequ ) THEN
613 rcond_tmp = zla_hercond_c( uplo, n, a, lda, af, ldaf,
614 $ ipiv,
615 $ s, .true., info, work, rwork )
616 ELSE
617 rcond_tmp = zla_hercond_c( uplo, n, a, lda, af, ldaf,
618 $ ipiv,
619 $ s, .false., info, work, rwork )
620 END IF
621 DO j = 1, nrhs
622*
623* Cap the error at 1.0.
624*
625 IF ( n_err_bnds .GE. la_linrx_err_i
626 $ .AND. err_bnds_norm( j, la_linrx_err_i ) .GT. 1.0d+0 )
627 $ err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
628*
629* Threshold the error (see LAWN).
630*
631 IF (rcond_tmp .LT. illrcond_thresh) THEN
632 err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
633 err_bnds_norm( j, la_linrx_trust_i ) = 0.0d+0
634 IF ( info .LE. n ) info = n + j
635 ELSE IF ( err_bnds_norm( j, la_linrx_err_i ) .LT. err_lbnd )
636 $ THEN
637 err_bnds_norm( j, la_linrx_err_i ) = err_lbnd
638 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
639 END IF
640*
641* Save the condition number.
642*
643 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
644 err_bnds_norm( j, la_linrx_rcond_i ) = rcond_tmp
645 END IF
646 END DO
647 END IF
648
649 IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 2 ) THEN
650*
651* Compute componentwise condition number cond(A*diag(Y(:,J))) for
652* each right-hand side using the current solution as an estimate of
653* the true solution. If the componentwise error estimate is too
654* large, then the solution is a lousy estimate of truth and the
655* estimated RCOND may be too optimistic. To avoid misleading users,
656* the inverse condition number is set to 0.0 when the estimated
657* cwise error is at least CWISE_WRONG.
658*
659 cwise_wrong = sqrt( dlamch( 'Epsilon' ) )
660 DO j = 1, nrhs
661 IF ( err_bnds_comp( j, la_linrx_err_i ) .LT. cwise_wrong )
662 $ THEN
663 rcond_tmp = zla_hercond_x( uplo, n, a, lda, af, ldaf,
664 $ ipiv, x( 1, j ), info, work, rwork )
665 ELSE
666 rcond_tmp = 0.0d+0
667 END IF
668*
669* Cap the error at 1.0.
670*
671 IF ( n_err_bnds .GE. la_linrx_err_i
672 $ .AND. err_bnds_comp( j, la_linrx_err_i ) .GT. 1.0d+0 )
673 $ err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
674*
675* Threshold the error (see LAWN).
676*
677 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
678 err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
679 err_bnds_comp( j, la_linrx_trust_i ) = 0.0d+0
680 IF ( .NOT. ignore_cwise
681 $ .AND. info.LT.n + j ) info = n + j
682 ELSE IF ( err_bnds_comp( j, la_linrx_err_i )
683 $ .LT. err_lbnd ) THEN
684 err_bnds_comp( j, la_linrx_err_i ) = err_lbnd
685 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
686 END IF
687*
688* Save the condition number.
689*
690 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
691 err_bnds_comp( j, la_linrx_rcond_i ) = rcond_tmp
692 END IF
693
694 END DO
695 END IF
696*
697 RETURN
698*
699* End of ZHERFSX
700*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhecon(uplo, n, a, lda, ipiv, anorm, rcond, work, info)
ZHECON
Definition zhecon.f:123
integer function ilaprec(prec)
ILAPREC
Definition ilaprec.f:56
double precision function zla_hercond_x(uplo, n, a, lda, af, ldaf, ipiv, x, info, work, rwork)
ZLA_HERCOND_X computes the infinity norm condition number of op(A)*diag(x) for Hermitian indefinite m...
double precision function zla_hercond_c(uplo, n, a, lda, af, ldaf, ipiv, c, capply, info, work, rwork)
ZLA_HERCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for Hermitian indefin...
subroutine zla_herfsx_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_HERFSX_EXTENDED improves the computed solution to a system of linear equations for Hermitian inde...
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlanhe(norm, uplo, n, a, lda, work)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhe.f:122
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: