Developer Reference for Intel® oneAPI Math Kernel Library for Fortran

ID 766686
Date 12/16/2022
Public

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

Document Table of Contents

p?posvx

Solves a symmetric or Hermitian positive definite system of linear equations.

Syntax

call psposvx(fact, uplo, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, equed, sr, sc, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, iwork, liwork, info)

call pdposvx(fact, uplo, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, equed, sr, sc, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, iwork, liwork, info)

call pcposvx(fact, uplo, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, equed, sr, sc, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, iwork, liwork, info)

call pzposvx(fact, uplo, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, equed, sr, sc, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, iwork, liwork, info)

Include Files

Description

The p?posvxroutine uses the Cholesky factorization A=UT*U or A=L*LT to compute the solution to a real or complex system of linear equations

A(ia:ia+n-1, ja:ja+n-1)*X = B(ib:ib+n-1, jb:jb+nrhs-1),

where A(ia:ia+n-1, ja:ja+n-1) is a n-by-n matrix and X and B(ib:ib+n-1,jb:jb+nrhs-1) are n-by-nrhs matrices.

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

In the following comments y denotes Y(iy:iy+m-1, jy:jy+k-1), an m-by-k matrix where y can be a, af, b and x.

The routine p?posvx performs the following steps:

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

    diag(sr)*A*diag(sc)*inv(diag(sc))*X = diag(sr)*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(sr)*A*diag(sc) and B by diag(sr)*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, if uplo = 'U', or

    A = L*LT, if uplo = 'L',

    where U is an upper triangular matrix and L is a lower triangular matrix.

  3. 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, steps 4-6 are skipped

  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(sr) so that it solves the original system before equilibration.

Input Parameters
fact

(global) CHARACTER. 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, 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.

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

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

uplo

(global) CHARACTER. Must be 'U' or 'L'.

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

n

(global) INTEGER. The order of the distributed matrix sub(A) (n 0).

nrhs

(global) INTEGER. The number of right-hand sides; the number of columns of the distributed submatrices B and X. (nrhs 0).

a

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

Pointer into the local memory to an array of local size (lld_a, LOCc(ja+n-1)). On entry, the symmetric/Hermitian matrix A, except if fact = 'F' and equed = 'Y', then A must contain the equilibrated matrix diag(sr)*A*diag(sc).

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.

ia, ja

(global) INTEGER. The row and column indices in the global matrix A indicating the first row and the first column of the submatrix A, respectively.

desca

(global and local) INTEGER array of size dlen_. The array descriptor for the distributed matrix A.

af

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

Pointer into the local memory to an array of local size (lld_af, LOCc(ja+n-1)).

If fact = 'F', then af is an input argument and on entry contains the triangular factor U or L from the Cholesky factorization A = UT*U or A = L*LT, in the same storage format as A. If equed'N', then af is the factored form of the equilibrated matrix diag(sr)*A*diag(sc).

iaf, jaf

(global) INTEGER. The row and column indices in the global matrix AF indicating the first row and the first column of the submatrix AF, respectively.

descaf

(global and local) INTEGER array of size dlen_. The array descriptor for the distributed matrix AF.

equed

(global) CHARACTER. 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 and A has been replaced by diag(sr)*A*diag(sc).

sr

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

Array of size lld_a.

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.

b

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

Pointer into the local memory to an array of local size (lld_b, LOCc(jb+nrhs-1)). On entry, the n-by-nrhs right-hand side matrix B.

ib, jb

(global) INTEGER. The row and column indices in the global matrix B indicating the first row and the first column of the submatrix B, respectively.

descb

(global and local) INTEGER. Array of size dlen_. The array descriptor for the distributed matrix B.

x

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

Pointer into the local memory to an array of local size (lld_x, LOCc(jx+nrhs-1)).

ix, jx

(global) INTEGER. The row and column indices in the global matrix X indicating the first row and the first column of the submatrix X, respectively.

descx

(global and local) INTEGER array of size dlen_. The array descriptor for the distributed matrix X.

work

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

Workspace array of size lwork.

lwork

(local or global) INTEGER.

The size of the array work. lwork is local input and must be at least lwork = max(p?pocon(lwork), p?porfs(lwork)) + LOCr(n_a).

lwork = 3*desca(lld_).

If lwork = -1, then lwork is global input and a workspace query is assumed; the routine only calculates the minimum and optimal size for all work arrays. Each of these values is returned in the first entry of the corresponding work array, and no error message is issued by pxerbla.

iwork

(local) INTEGER. Workspace array of size liwork.

liwork

(local or global)

INTEGER. The size of the array iwork. liwork is local input and must be at least liwork = desca(lld_)liwork = LOCr(n_a).

If liwork = -1, then liwork is global input and a workspace query is assumed; the routine only calculates the minimum and optimal size for all work arrays. Each of these values is returned in the first entry of the corresponding work array, and no error message is issued by pxerbla.

Output Parameters
a

On exit, if fact = 'E' and equed = 'Y', a is overwritten by diag(sr)*a*diag(sc).

af

If fact = 'N', then af 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 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 = UT*U or A = L*LT of the equilibrated matrix A (see the description of A for the form of the equilibrated matrix).

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

sr

This array is an output argument if fact'F'.

See the description of sr in Input Arguments section.

sc

This array is an output argument if fact'F'.

See the description of sc in Input Arguments section.

b

On exit, if equed = 'N', b is not modified; if trans = 'N' and equed = 'R' or 'B', b is overwritten by diag(r)*b; if trans = 'T' or 'C' and equed = 'C' or 'B', b is overwritten by diag(c)*b.

x

(local)

REAL for psposvx

DOUBLE PRECISION for pdposvx

COMPLEX for pcposvx

DOUBLE COMPLEX for pzposvx.

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'N', and the solution to the equilibrated system is

inv(diag(sc))*X if trans = 'N' and equed = 'C' or 'B', or

inv(diag(sr))*X if trans = 'T' or 'C' and equed = 'R' or 'B'.

rcond

(global)

REAL for single precision flavors.

DOUBLE PRECISION for double precision flavors.

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

REAL for single precision flavors.

DOUBLE PRECISION for double precision flavors.

Arrays of size at least max(LOC,n_b). The estimated forward error bounds for each solution vector X(j) (the j-th column of the solution matrix X). If xtrue is the true solution, ferr(j) bounds the magnitude of the largest entry in (X(j) - xtrue) divided by the magnitude of the largest entry in X(j). The quality of the error bound depends on the quality of the estimate of norm(inv(A)) computed in the code; if the estimate of norm(inv(A)) is accurate, the error bound is guaranteed.

berr

(local)

REAL for single precision flavors.

DOUBLE PRECISION for double precision flavors.

Arrays of size at least max(LOC,n_b). The componentwise relative backward error of each solution vector X(j) (the smallest relative change in any entry of A or B that makes X(j) an exact solution).

work(1)

(local) On exit, work(1) returns the minimal and optimal liwork.

info

(global) INTEGER.

If info=0, the execution is successful.

< 0: if info = -i, the i-th argument had an illegal value

> 0: if info = i, and i is n: 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 and error bounds could not be computed.

= n+1: rcond is less than machine precision. The factorization has been completed, but the matrix is singular to working precision, and the solution and error bounds have not been computed.

See Also