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

◆ sposvxx()

subroutine sposvxx ( character  fact,
character  uplo,
integer  n,
integer  nrhs,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( ldaf, * )  af,
integer  ldaf,
character  equed,
real, dimension( * )  s,
real, dimension( ldb, * )  b,
integer  ldb,
real, dimension( ldx, * )  x,
integer  ldx,
real  rcond,
real  rpvgrw,
real, dimension( * )  berr,
integer  n_err_bnds,
real, dimension( nrhs, * )  err_bnds_norm,
real, dimension( nrhs, * )  err_bnds_comp,
integer  nparams,
real, dimension( * )  params,
real, dimension( * )  work,
integer, dimension( * )  iwork,
integer  info 
)

SPOSVXX computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
    SPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
    to compute the solution to a real system of linear equations
    A * X = B, where A is an N-by-N symmetric positive definite matrix
    and X and B are N-by-NRHS matrices.

    If requested, both normwise and maximum componentwise error bounds
    are returned. SPOSVXX will return a solution with a tiny
    guaranteed error (O(eps) where eps is the working machine
    precision) unless the matrix is very ill-conditioned, in which
    case a warning is returned. Relevant condition numbers also are
    calculated and returned.

    SPOSVXX accepts user-provided factorizations and equilibration
    factors; see the definitions of the FACT and EQUED options.
    Solving with refinement and using a factorization from a previous
    SPOSVXX call will also produce a solution with either O(eps)
    errors or warnings, but we cannot make that claim for general
    user-provided factorizations and equilibration factors if they
    differ from what SPOSVXX would itself produce.
Description:
    The following steps are performed:

    1. If FACT = 'E', real scaling factors are computed to equilibrate
    the system:

      diag(S)*A*diag(S)     *inv(diag(S))*X = diag(S)*B

    Whether or not the system will be equilibrated depends on the
    scaling of the matrix A, but if equilibration is used, A is
    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.

    2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
    factor the matrix A (after equilibration if FACT = 'E') as
       A = U**T* U,  if UPLO = 'U', or
       A = L * L**T,  if UPLO = 'L',
    where U is an upper triangular matrix and L is a lower triangular
    matrix.

    3. If the leading principal minor of order i is not positive,
    then the routine returns with INFO = i. Otherwise, the factored
    form of A is used to estimate the condition number of the matrix
    A (see argument RCOND).  If the reciprocal of the condition number
    is less than machine precision, the routine still goes on to solve
    for X and compute error bounds as described below.

    4. The system of equations is solved for X using the factored form
    of A.

    5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
    the routine will use iterative refinement to try to get a small
    error and error bounds.  Refinement calculates the residual to at
    least twice the working precision.

    6. If equilibration was used, the matrix X is premultiplied by
    diag(S) so that it solves the original system before
    equilibration.
     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]FACT
          FACT is CHARACTER*1
     Specifies whether or not the factored form of the matrix A is
     supplied on entry, and if not, whether the matrix A should be
     equilibrated before it is factored.
       = 'F':  On entry, AF contains the factored form of A.
               If EQUED is not 'N', the matrix A has been
               equilibrated with scaling factors given by S.
               A and AF are not modified.
       = 'N':  The matrix A will be copied to AF and factored.
       = 'E':  The matrix A will be equilibrated if necessary, then
               copied to AF and factored.
[in]UPLO
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
     The number of linear equations, i.e., 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,out]A
          A is REAL array, dimension (LDA,N)
     On entry, the symmetric matrix A, except if FACT = 'F' and EQUED =
     'Y', then A must contain the equilibrated matrix
     diag(S)*A*diag(S).  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.  A is
     not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
     'N' on exit.

     On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
     diag(S)*A*diag(S).
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in,out]AF
          AF is REAL array, dimension (LDAF,N)
     If FACT = 'F', then AF is an input argument and on entry
     contains the triangular factor U or L from the Cholesky
     factorization A = U**T*U or A = L*L**T, in the same storage
     format as A.  If EQUED .ne. 'N', then AF is the factored
     form of the equilibrated matrix diag(S)*A*diag(S).

     If FACT = 'N', then AF is an output argument and on exit
     returns the triangular factor U or L from the Cholesky
     factorization A = U**T*U or A = L*L**T of the original
     matrix A.

     If FACT = 'E', then AF is an output argument and on exit
     returns the triangular factor U or L from the Cholesky
     factorization A = U**T*U or A = L*L**T of the equilibrated
     matrix A (see the description of A for the form of the
     equilibrated matrix).
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in,out]EQUED
          EQUED is CHARACTER*1
     Specifies the form of equilibration that was done.
       = 'N':  No equilibration (always true if FACT = 'N').
       = 'Y':  Both row and column equilibration, i.e., A has been
               replaced by diag(S) * A * diag(S).
     EQUED is an input argument if FACT = 'F'; otherwise, it is an
     output argument.
[in,out]S
          S is REAL array, dimension (N)
     The row 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,out]B
          B is REAL array, dimension (LDB,NRHS)
     On entry, the N-by-NRHS right hand side matrix B.
     On exit,
     if EQUED = 'N', B is not modified;
     if EQUED = 'Y', B is overwritten by diag(S)*B;
[in]LDB
          LDB is INTEGER
     The leading dimension of the array B.  LDB >= max(1,N).
[out]X
          X is REAL array, dimension (LDX,NRHS)
     If INFO = 0, the N-by-NRHS solution matrix X to the original
     system of equations.  Note that A and B are modified on exit if
     EQUED .ne. 'N', and the solution to the equilibrated system is
     inv(diag(S))*X.
[in]LDX
          LDX is INTEGER
     The leading dimension of the array X.  LDX >= max(1,N).
[out]RCOND
          RCOND is REAL
     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]RPVGRW
          RPVGRW is REAL
     Reciprocal pivot growth.  On exit, this contains the reciprocal
     pivot growth factor norm(A)/norm(U). The "max absolute element"
     norm is used.  If this is much less than 1, then the stability of
     the LU factorization of the (equilibrated) matrix A could be poor.
     This also means that the solution X, estimated condition numbers,
     and error bounds could be unreliable. If factorization fails with
     0<INFO<=N, then this contains the reciprocal pivot growth factor
     for the leading INFO columns of A.
[out]BERR
          BERR is REAL 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 REAL 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 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" 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) * slamch('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) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". 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 REAL 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 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" 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) * slamch('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) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". 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 REAL 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.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 REAL array, dimension (4*N)
[out]IWORK
          IWORK is INTEGER array, dimension (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 493 of file sposvxx.f.

497*
498* -- LAPACK driver routine --
499* -- LAPACK is a software package provided by Univ. of Tennessee, --
500* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
501*
502* .. Scalar Arguments ..
503 CHARACTER EQUED, FACT, UPLO
504 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
505 $ N_ERR_BNDS
506 REAL RCOND, RPVGRW
507* ..
508* .. Array Arguments ..
509 INTEGER IWORK( * )
510 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
511 $ X( LDX, * ), WORK( * )
512 REAL S( * ), PARAMS( * ), BERR( * ),
513 $ ERR_BNDS_NORM( NRHS, * ),
514 $ ERR_BNDS_COMP( NRHS, * )
515* ..
516*
517* ==================================================================
518*
519* .. Parameters ..
520 REAL ZERO, ONE
521 parameter( zero = 0.0e+0, one = 1.0e+0 )
522 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
523 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
524 INTEGER CMP_ERR_I, PIV_GROWTH_I
525 parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
526 $ berr_i = 3 )
527 parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
528 parameter( cmp_rcond_i = 7, cmp_err_i = 8,
529 $ piv_growth_i = 9 )
530* ..
531* .. Local Scalars ..
532 LOGICAL EQUIL, NOFACT, RCEQU
533 INTEGER INFEQU, J
534 REAL AMAX, BIGNUM, SMIN, SMAX,
535 $ SCOND, SMLNUM
536* ..
537* .. External Functions ..
538 EXTERNAL lsame, slamch, sla_porpvgrw
539 LOGICAL LSAME
540 REAL SLAMCH, SLA_PORPVGRW
541* ..
542* .. External Subroutines ..
543 EXTERNAL spoequb, spotrf, spotrs, slacpy, slaqsy,
545* ..
546* .. Intrinsic Functions ..
547 INTRINSIC max, min
548* ..
549* .. Executable Statements ..
550*
551 info = 0
552 nofact = lsame( fact, 'N' )
553 equil = lsame( fact, 'E' )
554 smlnum = slamch( 'Safe minimum' )
555 bignum = one / smlnum
556 IF( nofact .OR. equil ) THEN
557 equed = 'N'
558 rcequ = .false.
559 ELSE
560 rcequ = lsame( equed, 'Y' )
561 ENDIF
562*
563* Default is failure. If an input parameter is wrong or
564* factorization fails, make everything look horrible. Only the
565* pivot growth is set here, the rest is initialized in SPORFSX.
566*
567 rpvgrw = zero
568*
569* Test the input parameters. PARAMS is not tested until SPORFSX.
570*
571 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.
572 $ lsame( fact, 'F' ) ) THEN
573 info = -1
574 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
575 $ .NOT.lsame( uplo, 'L' ) ) THEN
576 info = -2
577 ELSE IF( n.LT.0 ) THEN
578 info = -3
579 ELSE IF( nrhs.LT.0 ) THEN
580 info = -4
581 ELSE IF( lda.LT.max( 1, n ) ) THEN
582 info = -6
583 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
584 info = -8
585 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
586 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
587 info = -9
588 ELSE
589 IF ( rcequ ) THEN
590 smin = bignum
591 smax = zero
592 DO 10 j = 1, n
593 smin = min( smin, s( j ) )
594 smax = max( smax, s( j ) )
595 10 CONTINUE
596 IF( smin.LE.zero ) THEN
597 info = -10
598 ELSE IF( n.GT.0 ) THEN
599 scond = max( smin, smlnum ) / min( smax, bignum )
600 ELSE
601 scond = one
602 END IF
603 END IF
604 IF( info.EQ.0 ) THEN
605 IF( ldb.LT.max( 1, n ) ) THEN
606 info = -12
607 ELSE IF( ldx.LT.max( 1, n ) ) THEN
608 info = -14
609 END IF
610 END IF
611 END IF
612*
613 IF( info.NE.0 ) THEN
614 CALL xerbla( 'SPOSVXX', -info )
615 RETURN
616 END IF
617*
618 IF( equil ) THEN
619*
620* Compute row and column scalings to equilibrate the matrix A.
621*
622 CALL spoequb( n, a, lda, s, scond, amax, infequ )
623 IF( infequ.EQ.0 ) THEN
624*
625* Equilibrate the matrix.
626*
627 CALL slaqsy( uplo, n, a, lda, s, scond, amax, equed )
628 rcequ = lsame( equed, 'Y' )
629 END IF
630 END IF
631*
632* Scale the right-hand side.
633*
634 IF( rcequ ) CALL slascl2( n, nrhs, s, b, ldb )
635*
636 IF( nofact .OR. equil ) THEN
637*
638* Compute the Cholesky factorization of A.
639*
640 CALL slacpy( uplo, n, n, a, lda, af, ldaf )
641 CALL spotrf( uplo, n, af, ldaf, info )
642*
643* Return if INFO is non-zero.
644*
645 IF( info.NE.0 ) THEN
646*
647* Pivot in column INFO is exactly 0
648* Compute the reciprocal pivot growth factor of the
649* leading rank-deficient INFO columns of A.
650*
651 rpvgrw = sla_porpvgrw( uplo, info, a, lda, af, ldaf, work )
652 RETURN
653 ENDIF
654 END IF
655*
656* Compute the reciprocal growth factor RPVGRW.
657*
658 rpvgrw = sla_porpvgrw( uplo, n, a, lda, af, ldaf, work )
659*
660* Compute the solution matrix X.
661*
662 CALL slacpy( 'Full', n, nrhs, b, ldb, x, ldx )
663 CALL spotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
664*
665* Use iterative refinement to improve the computed solution and
666* compute error bounds and backward error estimates for it.
667*
668 CALL sporfsx( uplo, equed, n, nrhs, a, lda, af, ldaf,
669 $ s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm,
670 $ err_bnds_comp, nparams, params, work, iwork, info )
671
672*
673* Scale solutions.
674*
675 IF ( rcequ ) THEN
676 CALL slascl2 ( n, nrhs, s, x, ldx )
677 END IF
678*
679 RETURN
680*
681* End of SPOSVXX
682*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
real function sla_porpvgrw(uplo, ncols, a, lda, af, ldaf, work)
SLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian...
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
subroutine slaqsy(uplo, n, a, lda, s, scond, amax, equed)
SLAQSY scales a symmetric/Hermitian matrix, using scaling factors computed by spoequ.
Definition slaqsy.f:133
subroutine slascl2(m, n, d, x, ldx)
SLASCL2 performs diagonal scaling on a matrix.
Definition slascl2.f:90
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine spoequb(n, a, lda, s, scond, amax, info)
SPOEQUB
Definition spoequb.f:118
subroutine sporfsx(uplo, equed, n, nrhs, a, lda, af, ldaf, s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
SPORFSX
Definition sporfsx.f:394
subroutine spotrf(uplo, n, a, lda, info)
SPOTRF
Definition spotrf.f:107
subroutine spotrs(uplo, n, nrhs, a, lda, b, ldb, info)
SPOTRS
Definition spotrs.f:110
Here is the call graph for this function:
Here is the caller graph for this function: