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

◆ zgbrfsx()

subroutine zgbrfsx ( character trans,
character equed,
integer n,
integer kl,
integer ku,
integer nrhs,
complex*16, dimension( ldab, * ) ab,
integer ldab,
complex*16, dimension( ldafb, * ) afb,
integer ldafb,
integer, dimension( * ) ipiv,
double precision, dimension( * ) r,
double precision, dimension( * ) c,
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 )

ZGBRFSX

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

Purpose:
!>
!>    ZGBRFSX improves the computed solution to a system of linear
!>    equations 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, R
!>    and C 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]TRANS
!>          TRANS is CHARACTER*1
!>     Specifies the form of the system of equations:
!>       = 'N':  A * X = B     (No transpose)
!>       = 'T':  A**T * X = B  (Transpose)
!>       = 'C':  A**H * X = B  (Conjugate transpose)
!> 
[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
!>       = 'R':  Row equilibration, i.e., A has been premultiplied by
!>               diag(R).
!>       = 'C':  Column equilibration, i.e., A has been postmultiplied
!>               by diag(C).
!>       = 'B':  Both row and column equilibration, i.e., A has been
!>               replaced by diag(R) * A * diag(C).
!>               The right hand side B has been changed accordingly.
!> 
[in]N
!>          N is INTEGER
!>     The order of the matrix A.  N >= 0.
!> 
[in]KL
!>          KL is INTEGER
!>     The number of subdiagonals within the band of A.  KL >= 0.
!> 
[in]KU
!>          KU is INTEGER
!>     The number of superdiagonals within the band of A.  KU >= 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]AB
!>          AB is COMPLEX*16 array, dimension (LDAB,N)
!>     The original band matrix A, stored in rows 1 to KL+KU+1.
!>     The j-th column of A is stored in the j-th column of the
!>     array AB as follows:
!>     AB(ku+1+i-j,j) = A(i,j) for max(1,j-ku)<=i<=min(n,j+kl).
!> 
[in]LDAB
!>          LDAB is INTEGER
!>     The leading dimension of the array AB.  LDAB >= KL+KU+1.
!> 
[in]AFB
!>          AFB is COMPLEX*16 array, dimension (LDAFB,N)
!>     Details of the LU factorization of the band matrix A, as
!>     computed by ZGBTRF.  U is stored as an upper triangular band
!>     matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and
!>     the multipliers used during the factorization are stored in
!>     rows KL+KU+2 to 2*KL+KU+1.
!> 
[in]LDAFB
!>          LDAFB is INTEGER
!>     The leading dimension of the array AFB.  LDAFB >= 2*KL*KU+1.
!> 
[in]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>     The pivot indices from ZGETRF; for 1<=i<=N, row i of the
!>     matrix was interchanged with row IPIV(i).
!> 
[in,out]R
!>          R is DOUBLE PRECISION array, dimension (N)
!>     The row scale factors for A.  If EQUED = 'R' or 'B', A is
!>     multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
!>     is not accessed.  R is an input argument if FACT = 'F';
!>     otherwise, R is an output argument.  If FACT = 'F' and
!>     EQUED = 'R' or 'B', each element of R must be positive.
!>     If R is output, each element of R is a power of the radix.
!>     If R is input, each element of R 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,out]C
!>          C is DOUBLE PRECISION array, dimension (N)
!>     The column scale factors for A.  If EQUED = 'C' or 'B', A is
!>     multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
!>     is not accessed.  C is an input argument if FACT = 'F';
!>     otherwise, C is an output argument.  If FACT = 'F' and
!>     EQUED = 'C' or 'B', each element of C must be positive.
!>     If C is output, each element of C is a power of the radix.
!>     If C is input, each element of C 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 ZGETRS.
!>     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 433 of file zgbrfsx.f.

439*
440* -- LAPACK computational routine --
441* -- LAPACK is a software package provided by Univ. of Tennessee, --
442* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
443*
444* .. Scalar Arguments ..
445 CHARACTER TRANS, EQUED
446 INTEGER INFO, LDAB, LDAFB, LDB, LDX, N, KL, KU, NRHS,
447 $ NPARAMS, N_ERR_BNDS
448 DOUBLE PRECISION RCOND
449* ..
450* .. Array Arguments ..
451 INTEGER IPIV( * )
452 COMPLEX*16 AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
453 $ X( LDX , * ),WORK( * )
454 DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
455 $ ERR_BNDS_NORM( NRHS, * ),
456 $ ERR_BNDS_COMP( NRHS, * ), RWORK( * )
457* ..
458*
459* ==================================================================
460*
461* .. Parameters ..
462 DOUBLE PRECISION ZERO, ONE
463 parameter( zero = 0.0d+0, one = 1.0d+0 )
464 DOUBLE PRECISION ITREF_DEFAULT, ITHRESH_DEFAULT
465 DOUBLE PRECISION COMPONENTWISE_DEFAULT, RTHRESH_DEFAULT
466 DOUBLE PRECISION DZTHRESH_DEFAULT
467 parameter( itref_default = 1.0d+0 )
468 parameter( ithresh_default = 10.0d+0 )
469 parameter( componentwise_default = 1.0d+0 )
470 parameter( rthresh_default = 0.5d+0 )
471 parameter( dzthresh_default = 0.25d+0 )
472 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
473 $ LA_LINRX_CWISE_I
474 parameter( la_linrx_itref_i = 1,
475 $ la_linrx_ithresh_i = 2 )
476 parameter( la_linrx_cwise_i = 3 )
477 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
478 $ LA_LINRX_RCOND_I
479 parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
480 parameter( la_linrx_rcond_i = 3 )
481* ..
482* .. Local Scalars ..
483 CHARACTER(1) NORM
484 LOGICAL ROWEQU, COLEQU, NOTRAN, IGNORE_CWISE
485 INTEGER J, TRANS_TYPE, PREC_TYPE, REF_TYPE, N_NORMS,
486 $ ITHRESH
487 DOUBLE PRECISION ANORM, RCOND_TMP, ILLRCOND_THRESH, ERR_LBND,
488 $ CWISE_WRONG, RTHRESH, UNSTABLE_THRESH
489* ..
490* .. External Subroutines ..
492* ..
493* .. Intrinsic Functions ..
494 INTRINSIC max, sqrt, transfer
495* ..
496* .. External Functions ..
497 EXTERNAL lsame, ilaprec
498 EXTERNAL dlamch, zlangb, zla_gbrcond_x,
500 DOUBLE PRECISION DLAMCH, ZLANGB, ZLA_GBRCOND_X, ZLA_GBRCOND_C
501 LOGICAL LSAME
502 INTEGER ILATRANS, ILAPREC
503* ..
504* .. Executable Statements ..
505*
506* Check the input parameters.
507*
508 info = 0
509 trans_type = ilatrans( trans )
510 ref_type = int( itref_default )
511 IF ( nparams .GE. la_linrx_itref_i ) THEN
512 IF ( params( la_linrx_itref_i ) .LT. 0.0d+0 ) THEN
513 params( la_linrx_itref_i ) = itref_default
514 ELSE
515 ref_type = params( la_linrx_itref_i )
516 END IF
517 END IF
518*
519* Set default parameters.
520*
521 illrcond_thresh = dble( n ) * dlamch( 'Epsilon' )
522 ithresh = int( ithresh_default )
523 rthresh = rthresh_default
524 unstable_thresh = dzthresh_default
525 ignore_cwise = componentwise_default .EQ. 0.0d+0
526*
527 IF ( nparams.GE.la_linrx_ithresh_i ) THEN
528 IF ( params( la_linrx_ithresh_i ).LT.0.0d+0 ) THEN
529 params( la_linrx_ithresh_i ) = ithresh
530 ELSE
531 ithresh = int( params( la_linrx_ithresh_i ) )
532 END IF
533 END IF
534 IF ( nparams.GE.la_linrx_cwise_i ) THEN
535 IF ( params( la_linrx_cwise_i ).LT.0.0d+0 ) THEN
536 IF ( ignore_cwise ) THEN
537 params( la_linrx_cwise_i ) = 0.0d+0
538 ELSE
539 params( la_linrx_cwise_i ) = 1.0d+0
540 END IF
541 ELSE
542 ignore_cwise = params( la_linrx_cwise_i ) .EQ. 0.0d+0
543 END IF
544 END IF
545 IF ( ref_type .EQ. 0 .OR. n_err_bnds .EQ. 0 ) THEN
546 n_norms = 0
547 ELSE IF ( ignore_cwise ) THEN
548 n_norms = 1
549 ELSE
550 n_norms = 2
551 END IF
552*
553 notran = lsame( trans, 'N' )
554 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
555 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
556*
557* Test input parameters.
558*
559 IF( trans_type.EQ.-1 ) THEN
560 info = -1
561 ELSE IF( .NOT.rowequ .AND. .NOT.colequ .AND.
562 $ .NOT.lsame( equed, 'N' ) ) THEN
563 info = -2
564 ELSE IF( n.LT.0 ) THEN
565 info = -3
566 ELSE IF( kl.LT.0 ) THEN
567 info = -4
568 ELSE IF( ku.LT.0 ) THEN
569 info = -5
570 ELSE IF( nrhs.LT.0 ) THEN
571 info = -6
572 ELSE IF( ldab.LT.kl+ku+1 ) THEN
573 info = -8
574 ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
575 info = -10
576 ELSE IF( ldb.LT.max( 1, n ) ) THEN
577 info = -13
578 ELSE IF( ldx.LT.max( 1, n ) ) THEN
579 info = -15
580 END IF
581 IF( info.NE.0 ) THEN
582 CALL xerbla( 'ZGBRFSX', -info )
583 RETURN
584 END IF
585*
586* Quick return if possible.
587*
588 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
589 rcond = 1.0d+0
590 DO j = 1, nrhs
591 berr( j ) = 0.0d+0
592 IF ( n_err_bnds .GE. 1 ) THEN
593 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
594 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
595 END IF
596 IF ( n_err_bnds .GE. 2 ) THEN
597 err_bnds_norm( j, la_linrx_err_i ) = 0.0d+0
598 err_bnds_comp( j, la_linrx_err_i ) = 0.0d+0
599 END IF
600 IF ( n_err_bnds .GE. 3 ) THEN
601 err_bnds_norm( j, la_linrx_rcond_i ) = 1.0d+0
602 err_bnds_comp( j, la_linrx_rcond_i ) = 1.0d+0
603 END IF
604 END DO
605 RETURN
606 END IF
607*
608* Default to failure.
609*
610 rcond = 0.0d+0
611 DO j = 1, nrhs
612 berr( j ) = 1.0d+0
613 IF ( n_err_bnds .GE. 1 ) THEN
614 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
615 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
616 END IF
617 IF ( n_err_bnds .GE. 2 ) THEN
618 err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
619 err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
620 END IF
621 IF ( n_err_bnds .GE. 3 ) THEN
622 err_bnds_norm( j, la_linrx_rcond_i ) = 0.0d+0
623 err_bnds_comp( j, la_linrx_rcond_i ) = 0.0d+0
624 END IF
625 END DO
626*
627* Compute the norm of A and the reciprocal of the condition
628* number of A.
629*
630 IF( notran ) THEN
631 norm = 'I'
632 ELSE
633 norm = '1'
634 END IF
635 anorm = zlangb( norm, n, kl, ku, ab, ldab, rwork )
636 CALL zgbcon( norm, n, kl, ku, afb, ldafb, ipiv, anorm, rcond,
637 $ work, rwork, info )
638*
639* Perform refinement on each right-hand side
640*
641 IF ( ref_type .NE. 0 .AND. info .EQ. 0 ) THEN
642
643 prec_type = ilaprec( 'E' )
644
645 IF ( notran ) THEN
646 CALL zla_gbrfsx_extended( prec_type, trans_type, n, kl,
647 $ ku,
648 $ nrhs, ab, ldab, afb, ldafb, ipiv, colequ, c, b,
649 $ ldb, x, ldx, berr, n_norms, err_bnds_norm,
650 $ err_bnds_comp, work, rwork, work(n+1),
651 $ transfer(rwork(1:2*n), (/ (zero, zero) /), n),
652 $ rcond, ithresh, rthresh, unstable_thresh, ignore_cwise,
653 $ info )
654 ELSE
655 CALL zla_gbrfsx_extended( prec_type, trans_type, n, kl,
656 $ ku,
657 $ nrhs, ab, ldab, afb, ldafb, ipiv, rowequ, r, b,
658 $ ldb, x, ldx, berr, n_norms, err_bnds_norm,
659 $ err_bnds_comp, work, rwork, work(n+1),
660 $ transfer(rwork(1:2*n), (/ (zero, zero) /), n),
661 $ rcond, ithresh, rthresh, unstable_thresh, ignore_cwise,
662 $ info )
663 END IF
664 END IF
665
666 err_lbnd = max( 10.0d+0,
667 $ sqrt( dble( n ) ) ) * dlamch( 'Epsilon' )
668 IF (n_err_bnds .GE. 1 .AND. n_norms .GE. 1) THEN
669*
670* Compute scaled normwise condition number cond(A*C).
671*
672 IF ( colequ .AND. notran ) THEN
673 rcond_tmp = zla_gbrcond_c( trans, n, kl, ku, ab, ldab,
674 $ afb,
675 $ ldafb, ipiv, c, .true., info, work, rwork )
676 ELSE IF ( rowequ .AND. .NOT. notran ) THEN
677 rcond_tmp = zla_gbrcond_c( trans, n, kl, ku, ab, ldab,
678 $ afb,
679 $ ldafb, ipiv, r, .true., info, work, rwork )
680 ELSE
681 rcond_tmp = zla_gbrcond_c( trans, n, kl, ku, ab, ldab,
682 $ afb,
683 $ ldafb, ipiv, c, .false., info, work, rwork )
684 END IF
685 DO j = 1, nrhs
686*
687* Cap the error at 1.0.
688*
689 IF ( n_err_bnds .GE. la_linrx_err_i
690 $ .AND. err_bnds_norm( j, la_linrx_err_i ) .GT. 1.0d+0)
691 $ err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
692*
693* Threshold the error (see LAWN).
694*
695 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
696 err_bnds_norm( j, la_linrx_err_i ) = 1.0d+0
697 err_bnds_norm( j, la_linrx_trust_i ) = 0.0d+0
698 IF ( info .LE. n ) info = n + j
699 ELSE IF ( err_bnds_norm( j, la_linrx_err_i ) .LT. err_lbnd )
700 $ THEN
701 err_bnds_norm( j, la_linrx_err_i ) = err_lbnd
702 err_bnds_norm( j, la_linrx_trust_i ) = 1.0d+0
703 END IF
704*
705* Save the condition number.
706*
707 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
708 err_bnds_norm( j, la_linrx_rcond_i ) = rcond_tmp
709 END IF
710
711 END DO
712 END IF
713
714 IF (n_err_bnds .GE. 1 .AND. n_norms .GE. 2) THEN
715*
716* Compute componentwise condition number cond(A*diag(Y(:,J))) for
717* each right-hand side using the current solution as an estimate of
718* the true solution. If the componentwise error estimate is too
719* large, then the solution is a lousy estimate of truth and the
720* estimated RCOND may be too optimistic. To avoid misleading users,
721* the inverse condition number is set to 0.0 when the estimated
722* cwise error is at least CWISE_WRONG.
723*
724 cwise_wrong = sqrt( dlamch( 'Epsilon' ) )
725 DO j = 1, nrhs
726 IF (err_bnds_comp( j, la_linrx_err_i ) .LT. cwise_wrong )
727 $ THEN
728 rcond_tmp = zla_gbrcond_x( trans, n, kl, ku, ab, ldab,
729 $ afb, ldafb, ipiv, x( 1, j ), info, work, rwork )
730 ELSE
731 rcond_tmp = 0.0d+0
732 END IF
733*
734* Cap the error at 1.0.
735*
736 IF ( n_err_bnds .GE. la_linrx_err_i
737 $ .AND. err_bnds_comp( j, la_linrx_err_i ) .GT. 1.0d+0 )
738 $ err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
739*
740* Threshold the error (see LAWN).
741*
742 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
743 err_bnds_comp( j, la_linrx_err_i ) = 1.0d+0
744 err_bnds_comp( j, la_linrx_trust_i ) = 0.0d+0
745 IF ( params( la_linrx_cwise_i ) .EQ. 1.0d+0
746 $ .AND. info.LT.n + j ) info = n + j
747 ELSE IF ( err_bnds_comp( j, la_linrx_err_i )
748 $ .LT. err_lbnd ) THEN
749 err_bnds_comp( j, la_linrx_err_i ) = err_lbnd
750 err_bnds_comp( j, la_linrx_trust_i ) = 1.0d+0
751 END IF
752*
753* Save the condition number.
754*
755 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
756 err_bnds_comp( j, la_linrx_rcond_i ) = rcond_tmp
757 END IF
758
759 END DO
760 END IF
761*
762 RETURN
763*
764* End of ZGBRFSX
765*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgbcon(norm, n, kl, ku, ab, ldab, ipiv, anorm, rcond, work, rwork, info)
ZGBCON
Definition zgbcon.f:146
integer function ilaprec(prec)
ILAPREC
Definition ilaprec.f:56
integer function ilatrans(trans)
ILATRANS
Definition ilatrans.f:56
double precision function zla_gbrcond_x(trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, x, info, work, rwork)
ZLA_GBRCOND_X computes the infinity norm condition number of op(A)*diag(x) for general banded matrice...
double precision function zla_gbrcond_c(trans, n, kl, ku, ab, ldab, afb, ldafb, ipiv, c, capply, info, work, rwork)
ZLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general banded ma...
subroutine zla_gbrfsx_extended(prec_type, trans_type, n, kl, ku, nrhs, ab, ldab, afb, ldafb, 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_GBRFSX_EXTENDED improves the computed solution to a system of linear equations for general banded...
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlangb(norm, n, kl, ku, ab, ldab, work)
ZLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlangb.f:123
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: