gbtf2#

Functions

void sgbtf2(
    const INT           m,
    const INT           n,
    const INT           kl,
    const INT           ku,
          f32* restrict AB,
    const INT           ldab,
          INT* restrict ipiv,
          INT*          info
);
void sgbtf2(const INT m, const INT n, const INT kl, const INT ku, f32 *restrict AB, const INT ldab, INT *restrict ipiv, INT *info)#

SGBTF2 computes an LU factorization of a real m-by-n band matrix A using partial pivoting with row interchanges.

This is the unblocked version of the algorithm, calling Level 2 BLAS.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

On exit, details of the factorization: U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 0 to kl+ku, and the multipliers used during the factorization are stored in rows kl+ku+1 to 2*kl+ku.

On entry (0-based indexing, rows 0 to 5, kl=2, ku=1, kv=kl+ku=3):

Row 0:  *    *    *   (fill-in storage)
Row 1:  *    *   a02  a13  a24  a35   (U superdiagonals)
Row 2:  *   a01  a12  a23  a34  a45   (U superdiagonals)
Row 3: a00  a11  a22  a33  a44  a55   (diagonal)
Row 4: a10  a21  a32  a43  a54   *    (L multipliers after factorization)
Row 5: a20  a31  a42  a53   *    *    (L multipliers after factorization)
On exit: Row 0: * * * u03 u14 u25 (U fill-in from pivoting) Row 1: * * u02 u13 u24 u35 (U superdiagonals) Row 2: * u01 u12 u23 u34 u45 (U superdiagonals) Row 3: u00 u11 u22 u33 u44 u55 (U diagonal) Row 4: m10 m21 m32 m43 m54 * (L multipliers) Row 5: m20 m31 m42 m53 * * (L multipliers)
Further Details:

The band storage scheme is illustrated by the following example, when m = n = 6, kl = 2, ku = 1:

Array elements marked * are not used by the routine.

Parameters

in
m

The number of rows of the matrix A. m >= 0.

in
n

The number of columns of the matrix A. n >= 0.

in
kl

The number of subdiagonals within the band of A. kl >= 0.

in
ku

The number of superdiagonals within the band of A. ku >= 0.

inout
AB

Double precision array, dimension (ldab, n). On entry, the matrix A in band storage, in rows kl to 2*kl+ku; rows 0 to kl-1 of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB[kl+ku+i-j + j*ldab] = A(i,j) for max(0,j-ku) <= i <= min(m-1,j+kl).

in
ldab

The leading dimension of the array AB. ldab >= 2*kl+ku+1.

out
ipiv

Integer array, dimension (min(m,n)). The pivot indices; for 0 <= i < min(m,n), row i of the matrix was interchanged with row ipiv[i]. 0-based indexing.

out
info

Exit status:

  • = 0: successful exit

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

  • > 0: if info = i, U(i-1,i-1) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

Functions

void dgbtf2(
    const INT           m,
    const INT           n,
    const INT           kl,
    const INT           ku,
          f64* restrict AB,
    const INT           ldab,
          INT* restrict ipiv,
          INT*          info
);
void dgbtf2(const INT m, const INT n, const INT kl, const INT ku, f64 *restrict AB, const INT ldab, INT *restrict ipiv, INT *info)#

DGBTF2 computes an LU factorization of a real m-by-n band matrix A using partial pivoting with row interchanges.

This is the unblocked version of the algorithm, calling Level 2 BLAS.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

On exit, details of the factorization: U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 0 to kl+ku, and the multipliers used during the factorization are stored in rows kl+ku+1 to 2*kl+ku.

On entry (0-based indexing, rows 0 to 5, kl=2, ku=1, kv=kl+ku=3):

Row 0:  *    *    *   (fill-in storage)
Row 1:  *    *   a02  a13  a24  a35   (U superdiagonals)
Row 2:  *   a01  a12  a23  a34  a45   (U superdiagonals)
Row 3: a00  a11  a22  a33  a44  a55   (diagonal)
Row 4: a10  a21  a32  a43  a54   *    (L multipliers after factorization)
Row 5: a20  a31  a42  a53   *    *    (L multipliers after factorization)
On exit: Row 0: * * * u03 u14 u25 (U fill-in from pivoting) Row 1: * * u02 u13 u24 u35 (U superdiagonals) Row 2: * u01 u12 u23 u34 u45 (U superdiagonals) Row 3: u00 u11 u22 u33 u44 u55 (U diagonal) Row 4: m10 m21 m32 m43 m54 * (L multipliers) Row 5: m20 m31 m42 m53 * * (L multipliers)
Further Details:

The band storage scheme is illustrated by the following example, when m = n = 6, kl = 2, ku = 1:

Array elements marked * are not used by the routine.

Parameters

in
m

The number of rows of the matrix A. m >= 0.

in
n

The number of columns of the matrix A. n >= 0.

in
kl

The number of subdiagonals within the band of A. kl >= 0.

in
ku

The number of superdiagonals within the band of A. ku >= 0.

inout
AB

Double precision array, dimension (ldab, n). On entry, the matrix A in band storage, in rows kl to 2*kl+ku; rows 0 to kl-1 of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB[kl+ku+i-j + j*ldab] = A(i,j) for max(0,j-ku) <= i <= min(m-1,j+kl).

in
ldab

The leading dimension of the array AB. ldab >= 2*kl+ku+1.

out
ipiv

Integer array, dimension (min(m,n)). The pivot indices; for 0 <= i < min(m,n), row i of the matrix was interchanged with row ipiv[i]. 0-based indexing.

out
info

Exit status:

  • = 0: successful exit

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

  • > 0: if info = i, U(i-1,i-1) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

Functions

void cgbtf2(
    const INT           m,
    const INT           n,
    const INT           kl,
    const INT           ku,
          c64* restrict AB,
    const INT           ldab,
          INT* restrict ipiv,
          INT*          info
);
void cgbtf2(const INT m, const INT n, const INT kl, const INT ku, c64 *restrict AB, const INT ldab, INT *restrict ipiv, INT *info)#

CGBTF2 computes an LU factorization of a complex m-by-n band matrix A using partial pivoting with row interchanges.

This is the unblocked version of the algorithm, calling Level 2 BLAS.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

On exit, details of the factorization: U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 0 to kl+ku, and the multipliers used during the factorization are stored in rows kl+ku+1 to 2*kl+ku.

On entry (0-based indexing, rows 0 to 5, kl=2, ku=1, kv=kl+ku=3):

Row 0:  *    *    *   (fill-in storage)
Row 1:  *    *   a02  a13  a24  a35   (U superdiagonals)
Row 2:  *   a01  a12  a23  a34  a45   (U superdiagonals)
Row 3: a00  a11  a22  a33  a44  a55   (diagonal)
Row 4: a10  a21  a32  a43  a54   *    (L multipliers after factorization)
Row 5: a20  a31  a42  a53   *    *    (L multipliers after factorization)
On exit: Row 0: * * * u03 u14 u25 (U fill-in from pivoting) Row 1: * * u02 u13 u24 u35 (U superdiagonals) Row 2: * u01 u12 u23 u34 u45 (U superdiagonals) Row 3: u00 u11 u22 u33 u44 u55 (U diagonal) Row 4: m10 m21 m32 m43 m54 * (L multipliers) Row 5: m20 m31 m42 m53 * * (L multipliers)
Further Details:

The band storage scheme is illustrated by the following example, when m = n = 6, kl = 2, ku = 1:

Array elements marked * are not used by the routine.

Parameters

in
m

The number of rows of the matrix A. m >= 0.

in
n

The number of columns of the matrix A. n >= 0.

in
kl

The number of subdiagonals within the band of A. kl >= 0.

in
ku

The number of superdiagonals within the band of A. ku >= 0.

inout
AB

Single complex array, dimension (ldab, n). On entry, the matrix A in band storage, in rows kl to 2*kl+ku; rows 0 to kl-1 of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB[kl+ku+i-j + j*ldab] = A(i,j) for max(0,j-ku) <= i <= min(m-1,j+kl).

in
ldab

The leading dimension of the array AB. ldab >= 2*kl+ku+1.

out
ipiv

Integer array, dimension (min(m,n)). The pivot indices; for 0 <= i < min(m,n), row i of the matrix was interchanged with row ipiv[i]. 0-based indexing.

out
info

Exit status:

  • = 0: successful exit

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

  • > 0: if info = i, U(i-1,i-1) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

Functions

void zgbtf2(
    const INT            m,
    const INT            n,
    const INT            kl,
    const INT            ku,
          c128* restrict AB,
    const INT            ldab,
          INT*  restrict ipiv,
          INT*           info
);
void zgbtf2(const INT m, const INT n, const INT kl, const INT ku, c128 *restrict AB, const INT ldab, INT *restrict ipiv, INT *info)#

ZGBTF2 computes an LU factorization of a complex m-by-n band matrix A using partial pivoting with row interchanges.

This is the unblocked version of the algorithm, calling Level 2 BLAS.

The factorization has the form A = P * L * U where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

On exit, details of the factorization: U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 0 to kl+ku, and the multipliers used during the factorization are stored in rows kl+ku+1 to 2*kl+ku.

On entry (0-based indexing, rows 0 to 5, kl=2, ku=1, kv=kl+ku=3):

Row 0:  *    *    *   (fill-in storage)
Row 1:  *    *   a02  a13  a24  a35   (U superdiagonals)
Row 2:  *   a01  a12  a23  a34  a45   (U superdiagonals)
Row 3: a00  a11  a22  a33  a44  a55   (diagonal)
Row 4: a10  a21  a32  a43  a54   *    (L multipliers after factorization)
Row 5: a20  a31  a42  a53   *    *    (L multipliers after factorization)
On exit: Row 0: * * * u03 u14 u25 (U fill-in from pivoting) Row 1: * * u02 u13 u24 u35 (U superdiagonals) Row 2: * u01 u12 u23 u34 u45 (U superdiagonals) Row 3: u00 u11 u22 u33 u44 u55 (U diagonal) Row 4: m10 m21 m32 m43 m54 * (L multipliers) Row 5: m20 m31 m42 m53 * * (L multipliers)
Further Details:

The band storage scheme is illustrated by the following example, when m = n = 6, kl = 2, ku = 1:

Array elements marked * are not used by the routine.

Parameters

in
m

The number of rows of the matrix A. m >= 0.

in
n

The number of columns of the matrix A. n >= 0.

in
kl

The number of subdiagonals within the band of A. kl >= 0.

in
ku

The number of superdiagonals within the band of A. ku >= 0.

inout
AB

Double complex array, dimension (ldab, n). On entry, the matrix A in band storage, in rows kl to 2*kl+ku; rows 0 to kl-1 of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB[kl+ku+i-j + j*ldab] = A(i,j) for max(0,j-ku) <= i <= min(m-1,j+kl).

in
ldab

The leading dimension of the array AB. ldab >= 2*kl+ku+1.

out
ipiv

Integer array, dimension (min(m,n)). The pivot indices; for 0 <= i < min(m,n), row i of the matrix was interchanged with row ipiv[i]. 0-based indexing.

out
info

Exit status:

  • = 0: successful exit

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

  • > 0: if info = i, U(i-1,i-1) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.