Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 11/07/2023
Public

A newer version of this document is available. Customers should click here to go to the newest version.

Document Table of Contents

?ppsvx

Uses the Cholesky factorization to compute the solution to the system of linear equations with a symmetric (Hermitian) positive definite packed coefficient matrix A, and provides error bounds on the solution.

Syntax

lapack_int LAPACKE_sppsvx( int matrix_layout, char fact, char uplo, lapack_int n, lapack_int nrhs, float* ap, float* afp, char* equed, float* s, float* b, lapack_int ldb, float* x, lapack_int ldx, float* rcond, float* ferr, float* berr );

lapack_int LAPACKE_dppsvx( int matrix_layout, char fact, char uplo, lapack_int n, lapack_int nrhs, double* ap, double* afp, char* equed, double* s, double* b, lapack_int ldb, double* x, lapack_int ldx, double* rcond, double* ferr, double* berr );

lapack_int LAPACKE_cppsvx( int matrix_layout, char fact, char uplo, lapack_int n, lapack_int nrhs, lapack_complex_float* ap, lapack_complex_float* afp, char* equed, float* s, lapack_complex_float* b, lapack_int ldb, lapack_complex_float* x, lapack_int ldx, float* rcond, float* ferr, float* berr );

lapack_int LAPACKE_zppsvx( int matrix_layout, char fact, char uplo, lapack_int n, lapack_int nrhs, lapack_complex_double* ap, lapack_complex_double* afp, char* equed, double* s, lapack_complex_double* b, lapack_int ldb, lapack_complex_double* x, lapack_int ldx, double* rcond, double* ferr, double* berr );

Include Files

  • mkl.h

Description

The routine uses the Cholesky factorization A=UT*U (real flavors) / A=UH*U (complex flavors) or A=L*LT (real flavors) / A=L*LH (complex flavors) to compute the solution to a real or complex system of linear equations A*X = B, where A is a n-by-n symmetric or Hermitian positive-definite matrix stored in packed format, the columns of matrix B are individual right-hand sides, and the columns of X are the corresponding solutions.

Error bounds on the solution and a condition estimate are also provided.

The routine ?ppsvx performs the following steps:

  1. If fact = 'E', real scaling factors s 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 = UT*U (real), A = UH*U (complex), if uplo = 'U',

    or A = L*LT (real), A = L*LH (complex), 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.

Input Parameters

matrix_layout

Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).

fact

Must be 'F', 'N', or 'E'.

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.

If fact = 'F': on entry, afp contains the factored form of A. If equed = 'Y', the matrix A has been equilibrated with scaling factors given by s.
ap and afp will not be modified.

If fact = 'N', the matrix A will be copied to afp and factored.

If fact = 'E', the matrix A will be equilibrated if necessary, then copied to afp and factored.

uplo

Must be 'U' or 'L'.

Indicates whether the upper or lower triangular part of A is stored:

If uplo = 'U', the upper triangle of A is stored.

If uplo = 'L', the lower triangle of A is stored.

n

The order of matrix A; n 0.

nrhs

The number of right-hand sides; the number of columns in B; nrhs 0.

ap, afp, b

Arrays: (size max(1,n*(n+1)/2), afp (size max(1,n*(n+1)/2), bof size max(1, ldb*nrhs) for column major layout and max(1, ldb*n) for row major layout.

The array ap contains the upper or lower triangle of the original symmetric/Hermitian matrix A in packed storage (see Matrix Storage Schemes). In case when fact = 'F' and equed = 'Y', ap must contain the equilibrated matrix diag(s)*A*diag(s).

The array afp is an input argument if fact = 'F' and contains the triangular factor U or L from the Cholesky factorization of A in the same storage format as A. If equed is not 'N', then afp is the factored form of the equilibrated matrix A.

The array b contains the matrix B whose columns are the right-hand sides for the systems of equations.

ldb

The leading dimension of b; ldb max(1, n) for column major layout and ldbnrhs for row major layout.

equed

Must be 'N' or 'Y'.

equed is an input argument if fact = 'F'. It specifies the form of equilibration that was done:

if equed = 'N', no equilibration was done (always true if fact = 'N');

if equed = 'Y', equilibration was done, that is, A has been replaced by diag(s)A*diag(s).

s

Array, size (n). The array s contains the scale factors for A. This array is an input argument if fact = 'F' only; otherwise it is an output argument.

If equed = 'N', s is not accessed.

If fact = 'F' and equed = 'Y', each element of s must be positive.

ldx

The leading dimension of the output array x; ldx max(1, n) for column major layout and ldxnrhs for row major layout.

Output Parameters

x

Array, size max(1, ldx*nrhs) for column major layout and max(1, ldx*n) for row major layout.

If info = 0 or info = n+1, the array x contains the 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.

ap

Array ap is not modified on exit if fact = 'F' or 'N', or if fact = 'E'and equed = 'N'.

If fact = 'E' and equed = 'Y', ap is overwritten by diag(s)*A*diag(s).

afp

If fact = 'N'or 'E', then afp is an output argument and on exit returns the triangular factor U or L from the Cholesky factorization A=UT*U or A=L*LT (real routines), A=UH*U or A=L*LH (complex routines) of the original matrix A (if fact = 'N'), or of the equilibrated matrix A (if fact = 'E'). See the description of ap for the form of the equilibrated matrix.

b

Overwritten by diag(s)*B, if equed = 'Y'; not changed if equed = 'N'.

s

This array is an output argument if fact'F'. See the description of s in Input Arguments section.

rcond

An 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

Array, size at least max(1, nrhs). Contains the estimated forward error bound for each solution vector xj(the j-th column of the solution matrix X). If xtrue is the true solution corresponding to xj, ferr[j-1] is an estimated upper bound for the magnitude of the largest element in (xj - xtrue) divided by the magnitude of the largest element in xj. The estimate is as reliable as the estimate for rcond, and is almost always a slight overestimate of the true error.

berr

Array, size at least max(1, nrhs). Contains the component-wise relative backward error for each solution vector xj, that is, the smallest relative change in any element of A or B that makes xj an exact solution.

equed

If fact'F', then equed is an output argument. It specifies the form of equilibration that was done (see the description of equed in Input Arguments section).

Return Values

This function returns a value info.

If info=0, the execution is successful.

If info = -i, parameter i had an illegal value.

If info = i, and in, the leading minor of order i (and therefore the matrix A itself) is not positive-definite, so the factorization could not be completed, and the solution and error bounds could not be computed; rcond = 0 is returned.

If info = i, and i = n + 1, then 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.