LAPACK  3.6.0
LAPACK: Linear Algebra PACKage
Collaboration diagram for real:

Functions

subroutine sposv (UPLO, N, NRHS, A, LDA, B, LDB, INFO)
  SPOSV computes the solution to system of linear equations A * X = B for PO matrices More...
 
subroutine sposvx (FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
  SPOSVX computes the solution to system of linear equations A * X = B for PO matrices More...
 
subroutine sposvxx (FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO)
  SPOSVXX computes the solution to system of linear equations A * X = B for PO matrices More...
 

Detailed Description

This is the group of real solve driver functions for PO matrices

Function Documentation

subroutine sposv ( character  UPLO,
integer  N,
integer  NRHS,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( ldb, * )  B,
integer  LDB,
integer  INFO 
)

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

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

Purpose:
 SPOSV computes 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.

 The Cholesky decomposition is used to factor A 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.  The factored form of A is then used to solve the system of
 equations A * X = B.
Parameters
[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 matrix B.  NRHS >= 0.
[in,out]A
          A is REAL array, dimension (LDA,N)
          On entry, the symmetric 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.

          On exit, if INFO = 0, the factor U or L from the Cholesky
          factorization A = U**T*U or A = L*L**T.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is REAL array, dimension (LDB,NRHS)
          On entry, the N-by-NRHS right hand side matrix B.
          On exit, if INFO = 0, the N-by-NRHS solution matrix X.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, the leading minor of order i of A is not
                positive definite, so the factorization could not be
                completed, and the solution has not been computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 132 of file sposv.f.

132 *
133 * -- LAPACK driver routine (version 3.4.0) --
134 * -- LAPACK is a software package provided by Univ. of Tennessee, --
135 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136 * November 2011
137 *
138 * .. Scalar Arguments ..
139  CHARACTER uplo
140  INTEGER info, lda, ldb, n, nrhs
141 * ..
142 * .. Array Arguments ..
143  REAL a( lda, * ), b( ldb, * )
144 * ..
145 *
146 * =====================================================================
147 *
148 * .. External Functions ..
149  LOGICAL lsame
150  EXTERNAL lsame
151 * ..
152 * .. External Subroutines ..
153  EXTERNAL spotrf, spotrs, xerbla
154 * ..
155 * .. Intrinsic Functions ..
156  INTRINSIC max
157 * ..
158 * .. Executable Statements ..
159 *
160 * Test the input parameters.
161 *
162  info = 0
163  IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
164  info = -1
165  ELSE IF( n.LT.0 ) THEN
166  info = -2
167  ELSE IF( nrhs.LT.0 ) THEN
168  info = -3
169  ELSE IF( lda.LT.max( 1, n ) ) THEN
170  info = -5
171  ELSE IF( ldb.LT.max( 1, n ) ) THEN
172  info = -7
173  END IF
174  IF( info.NE.0 ) THEN
175  CALL xerbla( 'SPOSV ', -info )
176  RETURN
177  END IF
178 *
179 * Compute the Cholesky factorization A = U**T*U or A = L*L**T.
180 *
181  CALL spotrf( uplo, n, a, lda, info )
182  IF( info.EQ.0 ) THEN
183 *
184 * Solve the system A*X = B, overwriting B with X.
185 *
186  CALL spotrs( uplo, n, nrhs, a, lda, b, ldb, info )
187 *
188  END IF
189  RETURN
190 *
191 * End of SPOSV
192 *
subroutine spotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
SPOTRS
Definition: spotrs.f:112
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:109
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine sposvx ( 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, dimension( * )  FERR,
real, dimension( * )  BERR,
real, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

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

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

Purpose:
 SPOSVX 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.

 Error bounds on the solution and a condition estimate are also
 provided.
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 i-by-i principal minor is not positive definite,
    then the routine returns with INFO = i. Otherwise, the factored
    form of A is used to estimate the condition number of the matrix
    A.  If the reciprocal of the condition number is less than machine
    precision, INFO = N+1 is returned as a warning, but 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. Iterative refinement is applied to improve the computed solution
    matrix and calculate error bounds and backward error estimates
    for it.

 6. If equilibration was used, the matrix X is premultiplied by
    diag(S) so that it solves the original system before
    equilibration.
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 = 'Y', the matrix A has been equilibrated
                  with scaling factors given by S.  A and AF will not
                  be 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':  Equilibration was done, 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 scale factors for A; not accessed if EQUED = 'N'.  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.
[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 or INFO = N+1, the N-by-NRHS solution matrix X to
          the original system of equations.  Note that if EQUED = 'Y',
          A and B are modified on exit, 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
          The estimate of the reciprocal condition number of the matrix
          A after equilibration (if done).  If RCOND is less than the
          machine precision (in particular, if RCOND = 0), the matrix
          is singular to working precision.  This condition is
          indicated by a return code of INFO > 0.
[out]FERR
          FERR is REAL array, dimension (NRHS)
          The estimated forward error bound for each solution vector
          X(j) (the j-th column of the solution matrix X).
          If XTRUE is the true solution corresponding to X(j), FERR(j)
          is an estimated upper bound for the magnitude of the largest
          element in (X(j) - XTRUE) divided by the magnitude of the
          largest element in X(j).  The estimate is as reliable as
          the estimate for RCOND, and is almost always a slight
          overestimate of the true error.
[out]BERR
          BERR is REAL array, dimension (NRHS)
          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).
[out]WORK
          WORK is REAL array, dimension (3*N)
[out]IWORK
          IWORK is INTEGER array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
          > 0: if INFO = i, and i is
                <= N:  the leading minor of order i of A is
                       not positive definite, so the factorization
                       could not be completed, and the solution has not
                       been computed. RCOND = 0 is returned.
                = N+1: U is nonsingular, but RCOND is less than machine
                       precision, meaning that the matrix is singular
                       to working precision.  Nevertheless, the
                       solution and error bounds are computed because
                       there are a number of situations where the
                       computed solution can be more accurate than the
                       value of RCOND would suggest.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
April 2012

Definition at line 309 of file sposvx.f.

309 *
310 * -- LAPACK driver routine (version 3.4.1) --
311 * -- LAPACK is a software package provided by Univ. of Tennessee, --
312 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
313 * April 2012
314 *
315 * .. Scalar Arguments ..
316  CHARACTER equed, fact, uplo
317  INTEGER info, lda, ldaf, ldb, ldx, n, nrhs
318  REAL rcond
319 * ..
320 * .. Array Arguments ..
321  INTEGER iwork( * )
322  REAL a( lda, * ), af( ldaf, * ), b( ldb, * ),
323  $ berr( * ), ferr( * ), s( * ), work( * ),
324  $ x( ldx, * )
325 * ..
326 *
327 * =====================================================================
328 *
329 * .. Parameters ..
330  REAL zero, one
331  parameter( zero = 0.0e+0, one = 1.0e+0 )
332 * ..
333 * .. Local Scalars ..
334  LOGICAL equil, nofact, rcequ
335  INTEGER i, infequ, j
336  REAL amax, anorm, bignum, scond, smax, smin, smlnum
337 * ..
338 * .. External Functions ..
339  LOGICAL lsame
340  REAL slamch, slansy
341  EXTERNAL lsame, slamch, slansy
342 * ..
343 * .. External Subroutines ..
344  EXTERNAL slacpy, slaqsy, spocon, spoequ, sporfs, spotrf,
345  $ spotrs, xerbla
346 * ..
347 * .. Intrinsic Functions ..
348  INTRINSIC max, min
349 * ..
350 * .. Executable Statements ..
351 *
352  info = 0
353  nofact = lsame( fact, 'N' )
354  equil = lsame( fact, 'E' )
355  IF( nofact .OR. equil ) THEN
356  equed = 'N'
357  rcequ = .false.
358  ELSE
359  rcequ = lsame( equed, 'Y' )
360  smlnum = slamch( 'Safe minimum' )
361  bignum = one / smlnum
362  END IF
363 *
364 * Test the input parameters.
365 *
366  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
367  $ THEN
368  info = -1
369  ELSE IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) )
370  $ THEN
371  info = -2
372  ELSE IF( n.LT.0 ) THEN
373  info = -3
374  ELSE IF( nrhs.LT.0 ) THEN
375  info = -4
376  ELSE IF( lda.LT.max( 1, n ) ) THEN
377  info = -6
378  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
379  info = -8
380  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
381  $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
382  info = -9
383  ELSE
384  IF( rcequ ) THEN
385  smin = bignum
386  smax = zero
387  DO 10 j = 1, n
388  smin = min( smin, s( j ) )
389  smax = max( smax, s( j ) )
390  10 CONTINUE
391  IF( smin.LE.zero ) THEN
392  info = -10
393  ELSE IF( n.GT.0 ) THEN
394  scond = max( smin, smlnum ) / min( smax, bignum )
395  ELSE
396  scond = one
397  END IF
398  END IF
399  IF( info.EQ.0 ) THEN
400  IF( ldb.LT.max( 1, n ) ) THEN
401  info = -12
402  ELSE IF( ldx.LT.max( 1, n ) ) THEN
403  info = -14
404  END IF
405  END IF
406  END IF
407 *
408  IF( info.NE.0 ) THEN
409  CALL xerbla( 'SPOSVX', -info )
410  RETURN
411  END IF
412 *
413  IF( equil ) THEN
414 *
415 * Compute row and column scalings to equilibrate the matrix A.
416 *
417  CALL spoequ( n, a, lda, s, scond, amax, infequ )
418  IF( infequ.EQ.0 ) THEN
419 *
420 * Equilibrate the matrix.
421 *
422  CALL slaqsy( uplo, n, a, lda, s, scond, amax, equed )
423  rcequ = lsame( equed, 'Y' )
424  END IF
425  END IF
426 *
427 * Scale the right hand side.
428 *
429  IF( rcequ ) THEN
430  DO 30 j = 1, nrhs
431  DO 20 i = 1, n
432  b( i, j ) = s( i )*b( i, j )
433  20 CONTINUE
434  30 CONTINUE
435  END IF
436 *
437  IF( nofact .OR. equil ) THEN
438 *
439 * Compute the Cholesky factorization A = U**T *U or A = L*L**T.
440 *
441  CALL slacpy( uplo, n, n, a, lda, af, ldaf )
442  CALL spotrf( uplo, n, af, ldaf, info )
443 *
444 * Return if INFO is non-zero.
445 *
446  IF( info.GT.0 )THEN
447  rcond = zero
448  RETURN
449  END IF
450  END IF
451 *
452 * Compute the norm of the matrix A.
453 *
454  anorm = slansy( '1', uplo, n, a, lda, work )
455 *
456 * Compute the reciprocal of the condition number of A.
457 *
458  CALL spocon( uplo, n, af, ldaf, anorm, rcond, work, iwork, info )
459 *
460 * Compute the solution matrix X.
461 *
462  CALL slacpy( 'Full', n, nrhs, b, ldb, x, ldx )
463  CALL spotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
464 *
465 * Use iterative refinement to improve the computed solution and
466 * compute error bounds and backward error estimates for it.
467 *
468  CALL sporfs( uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx,
469  $ ferr, berr, work, iwork, info )
470 *
471 * Transform the solution matrix X to a solution of the original
472 * system.
473 *
474  IF( rcequ ) THEN
475  DO 50 j = 1, nrhs
476  DO 40 i = 1, n
477  x( i, j ) = s( i )*x( i, j )
478  40 CONTINUE
479  50 CONTINUE
480  DO 60 j = 1, nrhs
481  ferr( j ) = ferr( j ) / scond
482  60 CONTINUE
483  END IF
484 *
485 * Set INFO = N+1 if the matrix is singular to working precision.
486 *
487  IF( rcond.LT.slamch( 'Epsilon' ) )
488  $ info = n + 1
489 *
490  RETURN
491 *
492 * End of SPOSVX
493 *
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine spotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
SPOTRS
Definition: spotrs.f:112
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:109
subroutine spoequ(N, A, LDA, S, SCOND, AMAX, INFO)
SPOEQU
Definition: spoequ.f:114
subroutine spocon(UPLO, N, A, LDA, ANORM, RCOND, WORK, IWORK, INFO)
SPOCON
Definition: spocon.f:123
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine sporfs(UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
SPORFS
Definition: sporfs.f:185
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
real function slansy(NORM, UPLO, N, A, LDA, WORK)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric matrix.
Definition: slansy.f:124
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:135

Here is the call graph for this function:

Here is the caller graph for this function:

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 i-by-i principal minor is not positive definite,
    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 .LT. 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 .LE. 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 .LT. 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.
Date
April 2012

Definition at line 499 of file sposvxx.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: