Developer Reference for Intel® oneAPI Math Kernel Library for Fortran

ID 766686
Date 10/31/2024
Public
Document Table of Contents

?pbsvx

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

Syntax

call spbsvx( fact, uplo, n, kd, nrhs, ab, ldab, afb, ldafb, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info )

call dpbsvx( fact, uplo, n, kd, nrhs, ab, ldab, afb, ldafb, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info )

call cpbsvx( fact, uplo, n, kd, nrhs, ab, ldab, afb, ldafb, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info )

call zpbsvx( fact, uplo, n, kd, nrhs, ab, ldab, afb, ldafb, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info )

call pbsvx( ab, b, x [,uplo] [,afb] [,fact] [,equed] [,s] [,ferr] [,berr] [,rcond] [,info] )

Include Files

  • mkl.fi, mkl_lapack.f90

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 band matrix, 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 ?pbsvx 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 band matrix and L is a lower triangular band 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

fact

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': on entry, afb contains the factored form of A. If equed = 'Y', the matrix A has been equilibrated with scaling factors given by s.
ab and afb will not be 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.

uplo

CHARACTER*1. 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

INTEGER. The order of matrix A; n 0.

kd

INTEGER. The number of superdiagonals or subdiagonals in the matrix A; kd 0.

nrhs

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

ab, afb, b, work

REAL for spbsvx

DOUBLE PRECISION for dpbsvx

COMPLEX for cpbsvx

DOUBLE COMPLEX for zpbsvx.

Arrays: ab(size ldab by *), afb(size ldafb by *), b(size ldb by *), work(*).

The array ab contains the upper or lower triangle of the matrix A in band storage (see Matrix Storage Schemes).

If fact = 'F' and equed = 'Y', then ab must contain the equilibrated matrix diag(s)*A*diag(s). The second dimension of ab must be at least max(1, n).

The array afb is an input argument if fact = 'F'. It contains the triangular factor U or L from the Cholesky factorization of the band matrix A in the same storage format as A. If equed = 'Y', then afb is the factored form of the equilibrated matrix A. The second dimension of afb must be at least max(1, n).

The array b contains the matrix B whose columns are the right-hand sides for the systems of equations. The second dimension of b must be at least max(1, nrhs).

work(*) is a workspace array.

The dimension of work must be at least max(1,3*n) for real flavors, and at least max(1,2*n) for complex flavors.

ldab

INTEGER. The leading dimension of ab; ldabkd+1.

ldafb

INTEGER. The leading dimension of afb; ldafbkd+1.

ldb

INTEGER. The leading dimension of b; ldb max(1, n).

equed

CHARACTER*1. 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

REAL for single precision flavors

DOUBLE PRECISION for double precision flavors.

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

INTEGER. The leading dimension of the output array x; ldx max(1, n).

iwork

INTEGER. Workspace array, size at least max(1, n); used in real flavors only.

rwork

REAL for cpbsvx

DOUBLE PRECISION for zpbsvx.

Workspace array, size at least max(1, n); used in complex flavors only.

Output Parameters

x

REAL for spbsvx

DOUBLE PRECISION for dpbsvx

COMPLEX for cpbsvx

DOUBLE COMPLEX for zpbsvx.

Array, size ldx by *.

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. The second dimension of x must be at least max(1,nrhs).

ab

On exit, if fact = 'E'and equed = 'Y', A is overwritten by diag(s)*A*diag(s).

afb

If fact = 'N'or 'E', then afb 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 ab 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

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.

Array, size at least max(1, nrhs). Contains 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

REAL for single precision flavors

DOUBLE PRECISION for double precision flavors.

Array, size at least max(1, nrhs). Contains the component-wise relative backward error for each solution vector x(j), that is, the smallest relative change in any element of A or B that makes x(j) 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).

info

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

If info = -i, the i-th parameter 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.

LAPACK 95 Interface Notes

Routines in Fortran 95 interface have fewer arguments in the calling sequence than their FORTRAN 77 counterparts. For general conventions applied to skip redundant or reconstructible arguments, see LAPACK 95 Interface Conventions.

Specific details for the routine pbsvx interface are as follows:

ab

Holds the array A of size (kd+1,n).

b

Holds the matrix B of size (n,nrhs).

x

Holds the matrix X of size (n,nrhs).

afb

Holds the array AF of size (kd+1,n).

s

Holds the vector with the number of elements n. Default value for each element is s(i) = 1.0_WP.

ferr

Holds the vector with the number of elements nrhs.

berr

Holds the vector with the number of elements nrhs.

uplo

Must be 'U' or 'L'. The default value is 'U'.

fact

Must be 'N', 'E', or 'F'. The default value is 'N'. If fact = 'F', then af must be present; otherwise, an error is returned.

equed

Must be 'N' or 'Y'. The default value is 'N'.