sposvx


NAME

sposvx - use 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,


SYNOPSIS

  SUBROUTINE SPOSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, 
 *      SCALE, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, WORK2, INFO)
  CHARACTER * 1 FACT, UPLO, EQUED
  INTEGER N, NRHS, LDA, LDAF, LDB, LDX, INFO
  INTEGER WORK2(*)
  REAL RCOND
  REAL A(LDA,*), AF(LDAF,*), SCALE(*), B(LDB,*), X(LDX,*), FERR(*), BERR(*), WORK(*)
 
  SUBROUTINE SPOSVX_64( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, 
 *      SCALE, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, WORK2, INFO)
  CHARACTER * 1 FACT, UPLO, EQUED
  INTEGER*8 N, NRHS, LDA, LDAF, LDB, LDX, INFO
  INTEGER*8 WORK2(*)
  REAL RCOND
  REAL A(LDA,*), AF(LDAF,*), SCALE(*), B(LDB,*), X(LDX,*), FERR(*), BERR(*), WORK(*)
 

F95 INTERFACE

  SUBROUTINE POSVX( FACT, UPLO, [N], [NRHS], A, [LDA], AF, [LDAF], 
 *       EQUED, SCALE, B, [LDB], X, [LDX], RCOND, FERR, BERR, [WORK], 
 *       [WORK2], [INFO])
  CHARACTER(LEN=1) :: FACT, UPLO, EQUED
  INTEGER :: N, NRHS, LDA, LDAF, LDB, LDX, INFO
  INTEGER, DIMENSION(:) :: WORK2
  REAL :: RCOND
  REAL, DIMENSION(:) :: SCALE, FERR, BERR, WORK
  REAL, DIMENSION(:,:) :: A, AF, B, X
 
  SUBROUTINE POSVX_64( FACT, UPLO, [N], [NRHS], A, [LDA], AF, [LDAF], 
 *       EQUED, SCALE, B, [LDB], X, [LDX], RCOND, FERR, BERR, [WORK], 
 *       [WORK2], [INFO])
  CHARACTER(LEN=1) :: FACT, UPLO, EQUED
  INTEGER(8) :: N, NRHS, LDA, LDAF, LDB, LDX, INFO
  INTEGER(8), DIMENSION(:) :: WORK2
  REAL :: RCOND
  REAL, DIMENSION(:) :: SCALE, FERR, BERR, WORK
  REAL, DIMENSION(:,:) :: A, AF, B, X
 

C INTERFACE

#include <sunperf.h>

void sposvx(char fact, char uplo, int n, int nrhs, float *a, int lda, float *af, int ldaf, char equed, float *scale, float *b, int ldb, float *x, int ldx, float *rcond, float *ferr, float *berr, int *info);

void sposvx_64(char fact, char uplo, long n, long nrhs, float *a, long lda, float *af, long ldaf, char equed, float *scale, float *b, long ldb, float *x, long ldx, float *rcond, float *ferr, float *berr, long *info);


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.

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.


ARGUMENTS

* FACT (input)
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.

* UPLO (input)
* N (input)
The number of linear equations, i.e., the order of the matrix A. N >= 0.

* NRHS (input)
The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >= 0.

* A (input/output)
On entry, the symmetric matrix A, except if FACT = 'F' and EQUED = 'Y', then A must contain the equilibrated matrix diag(SCALE)*A*diag(SCALE). 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(SCALE)*A*diag(SCALE).

* LDA (input)
The leading dimension of the array A. LDA >= max(1,N).

* AF (input/output)
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(SCALE)*A*diag(SCALE).

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).

* LDAF (input)
The leading dimension of the array AF. LDAF >= max(1,N).

* EQUED (input)
Specifies the form of equilibration that was done.

* SCALE (input/output)
The scale factors for A; not accessed if EQUED = 'N'. SCALE is an input argument if FACT = 'F'; otherwise, SCALE is an output argument. If FACT = 'F' and EQUED = 'Y', each element of SCALE must be positive.

* B (input/output)
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(SCALE) * B.

* LDB (input)
The leading dimension of the array B. LDB >= max(1,N).

* X (output)
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(SCALE))*X.

* LDX (input)
The leading dimension of the array X. LDX >= max(1,N).

* RCOND (output)
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.

* FERR (output)
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.

* BERR (output)
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).

* WORK (workspace)
dimension(2*N)

* WORK2 (workspace)
dimension(N)

* INFO (output)