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?gesvx

Uses the LU factorization to compute the solution to the system of linear equations with a square matrix A and multiple right-hand sides, and provides error bounds on the solution.

Syntax

call psgesvx(fact, trans, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, ipiv, equed, r, c, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, iwork, liwork, info)

call pdgesvx(fact, trans, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, ipiv, equed, r, c, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, iwork, liwork, info)

call pcgesvx(fact, trans, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, ipiv, equed, r, c, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, rwork, lrwork, info)

call pzgesvx(fact, trans, n, nrhs, a, ia, ja, desca, af, iaf, jaf, descaf, ipiv, equed, r, c, b, ib, jb, descb, x, ix, jx, descx, rcond, ferr, berr, work, lwork, rwork, lrwork, info)

Include Files

Description

The p?gesvx routine uses the LU factorization to compute the solution to a real or complex system of linear equations AX = B, where A denotes the n-by-n submatrix A(ia:ia+n-1, ja:ja+n-1), B denotes the n-by-nrhs submatrix B(ib:ib+n-1, jb:jb+nrhs-1) and X denotes the n-by-nrhs submatrix X(ix:ix+n-1, jx:jx+nrhs-1).

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

In the following description, af stands for the subarray af(iaf:iaf+n-1, jaf:jaf+n-1).

The routine p?gesvx performs the following steps:

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

    trans = 'N': diag(R)*A*diag(C) *diag(C)-1*X = diag(R)*B

    trans = 'T': (diag(R)*A*diag(C))T *diag(R)-1*X = diag(C)*B

    trans = 'C': (diag(R)*A*diag(C))H *diag(R)-1*X = diag(C)*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(R)*A*diag(C) and B by diag(R)*B (if trans='N') or diag(c)*B (if trans = 'T' or 'C').

  2. If fact = 'N' or 'E', the LU decomposition is used to factor the matrix A (after equilibration if fact = 'E') as A = PLU, where P is a permutation matrix, L is a unit lower triangular matrix, and U is upper triangular.

  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 relative 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(C) (if trans = 'N') or diag(R) (if trans = 'T' or 'C') so that it solves the original system before equilibration.

Input Parameters
fact

(global) CHARACTER*1. 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' then, on entry, af and ipiv contain the factored form of A. If equed is not 'N', the matrix A has been equilibrated with scaling factors given by r and c. Arrays a, af, and ipiv are not modified.

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

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

trans

(global) CHARACTER*1. Must be 'N', 'T', or 'C'.

Specifies the form of the system of equations:

If trans = 'N', the system has the form A*X = B (No transpose);

If trans = 'T', the system has the form AT*X = B (Transpose);

If trans = 'C', the system has the form AH*X = B (Conjugate transpose);

n

(global) INTEGER. The number of linear equations; the order of the submatrix 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, af, b, work

(local)

REAL for psgesvx

DOUBLE PRECISION for pdgesvx

COMPLEX for pcgesvx

DOUBLE COMPLEX for pzgesvx.

Pointers into the local memory to arrays of local size a(lld_a,LOCc(ja+n-1)), af(lld_af,LOCc(ja+n-1)), b(lld_b,LOCc(jb+nrhs-1)), work(lwork).

The array a contains the matrix A. If fact = 'F' and equed is not 'N', then A must have been equilibrated by the scaling factors in r and/or c.

The array af is an input argument if fact = 'F'. In this case it contains on entry the factored form of the matrix A, that is, the factors L and U from the factorization A = P*L*U as computed by p?getrf. If equed is not 'N', then af is the factored form of the equilibrated matrix A.

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

work is a workspace array. The size of work is (lwork).

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(ia:ia+n-1, ja:ja+n-1), respectively.

desca

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

iaf, jaf

(global) INTEGER. The row and column indices in the global matrix AF indicating the first row and the first column of the subarray af(iaf:iaf+n-1, jaf:jaf+n-1), respectively.

descaf

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

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(ib:ib+n-1, jb:jb+nrhs-1), respectively.

descb

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

ipiv

(local) INTEGER Array of size LOCr(m_a)+mb_a.

The array ipiv is an input argument if fact = 'F' .

On entry, it contains the pivot indices from the factorization A = P*L*U as computed by p?getrf; (local) row i of the matrix was interchanged with the (global) row ipiv(i).

This array must be aligned with A(ia:ia+n-1, *).

equed

(global) CHARACTER*1. Must be 'N', 'R', 'C', or 'B'. 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 = 'R', row equilibration was done, that is, A has been premultiplied by diag(r);

If equed = 'C', column equilibration was done, that is, A has been postmultiplied by diag(c);

If equed = 'B', both row and column equilibration was done; A has been replaced by diag(r)*A*diag(c).

r, c

(local) REAL for single precision flavors;

DOUBLE PRECISION for double precision flavors.

Arrays of size LOCr(m_a) and LOCc(n_a), respectively.

The array r contains the row scale factors for A, and the array c contains the column scale factors for A. These arrays are input arguments if fact = 'F' only; otherwise they are output arguments. If equed = 'R' or 'B', A is multiplied on the left by diag(r); if equed = 'N' or 'C', r is not accessed.

If fact = 'F' and equed = 'R' or 'B', each element of r must be positive.

If equed = 'C' or 'B', A is multiplied on the right by diag(c); if equed = 'N' or 'R', c is not accessed.

If fact = 'F' and equed = 'C' or 'B', each element of c must be positive. Array r is replicated in every process column, and is aligned with the distributed matrix A. Array c is replicated in every process row, and is aligned with the distributed matrix A.

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(ix:ix+n-1, jx:jx+nrhs-1), respectively.

descx

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

lwork

(local or global) INTEGER. The size of the array work ; must be at least max(p?gecon(lwork), p?gerfs(lwork))+LOCr(n_a).

iwork

(local, psgesvx/pdgesvx only)INTEGER. Workspace array. The size of iwork is (liwork).

liwork

(local, psgesvx/pdgesvx only)INTEGER. The size of the array iwork , must be at least LOCr(n_a).

rwork

(local) REAL for pcgesvx

DOUBLE PRECISION for pzgesvx.

Workspace array, used in complex flavors only.

The size of rwork is (lrwork).

lrwork

(local or global, pcgesvx/pzgesvx only)INTEGER. The size of the array rwork;must be at least 2*LOCc(n_a) .

Output Parameters
x

(local)

REAL for psgesvx

DOUBLE PRECISION for pdgesvx

COMPLEX for pcgesvx

DOUBLE COMPLEX for pzgesvx.

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

If info = 0, the array x contains the 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:

diag(C)-1*X, if trans = 'N' and equed = 'C' or 'B'; and diag(R)-1*X, if trans = 'T' or 'C' and equed = 'R' or 'B'.

a

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

If equed'N', A is scaled on exit as follows:

equed = 'R': A = diag(R)*A

equed = 'C': A = A*diag(c)

equed = 'B': A = diag(R)*A*diag(c)

af

If fact = 'N' or 'E', then af is an output argument and on exit returns the factors L and U from the factorization A = P*L*U of the original matrix A (if fact = 'N') or of the equilibrated matrix A (if fact = 'E'). See the description of a for the form of the equilibrated matrix.

b

Overwritten by diag(R)*B if trans = 'N' and equed = 'R' or 'B';

overwritten by diag(c)*B if trans = 'T' and equed = 'C' or 'B'; not changed if equed = 'N'.

r, c

These arrays are output arguments if fact'F'.

See the description of r, c in Input Arguments section.

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). The routine sets rcond =0 if the estimate underflows; in this case the matrix is singular (to working precision). However, anytime rcond is small compared to 1.0, for the working precision, the matrix may be poorly conditioned or even singular.

ferr, berr

(local)REAL for single precision flavors

DOUBLE PRECISION for double precision flavors.

Arrays of size LOCc(n_b) each. Contain the component-wise forward and relative backward errors, respectively, for each solution vector.

Arrays ferr and berr are both replicated in every process row, and are aligned with the matrices B and X.

ipiv

If fact = 'N' or 'E', then ipiv is an output argument and on exit contains the pivot indices from the factorization A = P*L*U of the original matrix A (if fact = 'N') or of the equilibrated matrix A (if fact = 'E').

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

work(1)

If info=0, on exit work(1) returns the minimum value of lwork required for optimum performance.

iwork(1)

If info=0, on exit iwork(1) returns the minimum value of liwork required for optimum performance.

rwork(1)

If info=0, on exit rwork(1) returns the minimum value of lrwork required for optimum performance.

info

INTEGER. If info=0, the execution is successful.

info < 0: if the ith argument is an array and the jth entry had an illegal value, then info = -(i*100+j); if the ith argument is a scalar and had an illegal value, then info = -i. If info = i, and in, then U(i,i) 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. If info = i, and i = n +1, then U is nonsingular, but 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