Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 3/22/2024
Public

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

Document Table of Contents

?gbsvx

Computes the solution to the real or complex system of linear equations with a band coefficient matrix A and multiple right-hand sides, and provides error bounds on the solution.

Syntax

lapack_int LAPACKE_sgbsvx( int matrix_layout, char fact, char trans, lapack_int n, lapack_int kl, lapack_int ku, lapack_int nrhs, float* ab, lapack_int ldab, float* afb, lapack_int ldafb, lapack_int* ipiv, char* equed, float* r, float* c, float* b, lapack_int ldb, float* x, lapack_int ldx, float* rcond, float* ferr, float* berr, float* rpivot );

lapack_int LAPACKE_dgbsvx( int matrix_layout, char fact, char trans, lapack_int n, lapack_int kl, lapack_int ku, lapack_int nrhs, double* ab, lapack_int ldab, double* afb, lapack_int ldafb, lapack_int* ipiv, char* equed, double* r, double* c, double* b, lapack_int ldb, double* x, lapack_int ldx, double* rcond, double* ferr, double* berr, double* rpivot );

lapack_int LAPACKE_cgbsvx( int matrix_layout, char fact, char trans, lapack_int n, lapack_int kl, lapack_int ku, lapack_int nrhs, lapack_complex_float* ab, lapack_int ldab, lapack_complex_float* afb, lapack_int ldafb, lapack_int* ipiv, char* equed, float* r, float* c, lapack_complex_float* b, lapack_int ldb, lapack_complex_float* x, lapack_int ldx, float* rcond, float* ferr, float* berr, float* rpivot );

lapack_int LAPACKE_zgbsvx( int matrix_layout, char fact, char trans, lapack_int n, lapack_int kl, lapack_int ku, lapack_int nrhs, lapack_complex_double* ab, lapack_int ldab, lapack_complex_double* afb, lapack_int ldafb, lapack_int* ipiv, char* equed, double* r, double* c, lapack_complex_double* b, lapack_int ldb, lapack_complex_double* x, lapack_int ldx, double* rcond, double* ferr, double* berr, double* rpivot );

Include Files

  • mkl.h

Description

The routine uses the LU factorization to compute the solution to a real or complex system of linear equations A*X = B, AT*X = B, or AH*X = B, where A is a band matrix of order n with kl subdiagonals and ku superdiagonals, 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 ?gbsvx 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) *inv(diag(c))*X = diag(r)*B

    trans = 'T': (diag(r)*A*diag(c))T *inv(diag(r))*X = diag(c)*B

    trans = 'C': (diag(r)*A*diag(c))H *inv(diag(r))*X = diag(c)*B

    Whether 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 = L*U, where L is a product of permutation and unit lower triangular matrices with kl subdiagonals, and U is upper triangular with kl+ku superdiagonals.

  3. If some Ui,i = 0, so that U is exactly singular, 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(c) (if trans = 'N') or diag(r) (if trans = 'T' or 'C') 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 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, afb and ipiv contain the factored form of A. If equed is not 'N', the matrix A is equilibrated with scaling factors given by r and c.
ab, afb, and ipiv are not modified.

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

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

trans

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 (Transpose for real flavors, conjugate transpose for complex flavors).

n

The number of linear equations, the order of the matrix A; n 0.

kl

The number of subdiagonals within the band of A; kl 0.

ku

The number of superdiagonals within the band of A; ku 0.

nrhs

The number of right hand sides, the number of columns of the matrices B and X; nrhs 0.

ab, afb, b

Arrays: ab (max(ldab*n)), afb (max(ldafb*n)), b(max(1, ldb*nrhs) for column major layout and max(1, ldb*n) for row major layout).

The array ab contains the matrix A in band storage (see Matrix Storage Schemes). 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 afb is an input argument if fact = 'F'. It contains the factored form of the matrix A, that is, the factors L and U from the factorization A = P*L*U as computed by ?gbtrf. U is stored as an upper triangular band matrix with kl + ku superdiagonals.L is stored as lower triangular band matrix with kl subdiagonals. If equed is not 'N', then afb 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.

ldab

The leading dimension of ab; ldabkl+ku+1.

ldafb

The leading dimension of afb; ldafb 2*kl+ku+1.

ldb

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

ipiv

Array, size at least max(1, n). The array ipiv is an input argument if fact = 'F'. It contains the pivot indices from the factorization A = P*L*U as computed by ?gbtrf; row i of the matrix was interchanged with row ipiv[i-1].

equed

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, that is, A has been replaced by diag(r)*A*diag(c).

r, c

Arrays: r (size n), c (size n).

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.

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 A and B are modified on exit if equed'N', and the solution to the equilibrated system is: inv(diag(c))*X, if trans = 'N' and equed = 'C'or 'B'; inv(diag(r))*X, if trans = 'T' or 'C' and equed = 'R' or 'B'.

ab

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

afb

If fact = 'N' or 'E', then afb is an output argument and on exit returns details of the LU factorization of the original matrix A (if fact = 'N') or of the equilibrated matrix A (if fact = 'E'). See the description of ab 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' or 'C' 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

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.

ipiv

If fact = 'N' or 'E', then ipiv is an output argument and on exit contains the pivot indices from the factorization A = 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).

rpivot

On exit, rpivot contains the reciprocal pivot growth factor:

If rpivot 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, condition estimator rcond, and forward error bound ferr could be unreliable. If factorization fails with 0 < infon, then rpivot contains the reciprocal pivot growth factor for the leading info columns of A.

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, then Ui, 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; 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.