Visible to Intel only — GUID: GUID-E040EE45-DBF5-4A09-A319-15DA1D1B8F4D
Visible to Intel only — GUID: GUID-E040EE45-DBF5-4A09-A319-15DA1D1B8F4D
?pbtrf
Computes the Cholesky factorization of a symmetric (Hermitian) positive-definite band matrix.
Syntax
lapack_int LAPACKE_spbtrf (int matrix_layout , char uplo , lapack_int n , lapack_int kd , float * ab , lapack_int ldab );
lapack_int LAPACKE_dpbtrf (int matrix_layout , char uplo , lapack_int n , lapack_int kd , double * ab , lapack_int ldab );
lapack_int LAPACKE_cpbtrf (int matrix_layout , char uplo , lapack_int n , lapack_int kd , lapack_complex_float * ab , lapack_int ldab );
lapack_int LAPACKE_zpbtrf (int matrix_layout , char uplo , lapack_int n , lapack_int kd , lapack_complex_double * ab , lapack_int ldab );
Include Files
- mkl.h
Description
The routine forms the Cholesky factorization of a symmetric positive-definite or, for complex data, Hermitian positive-definite band matrix A:
A = UT*U for real data, A = UH*U for complex data | if uplo='U' |
A = L*LT for real data, A = L*LH for complex data | if uplo='L' |
where L is a lower triangular matrix and U is upper triangular.
This routine supports the Progress Routine feature. See Progress Function for details.
Input Parameters
matrix_layout |
Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR). |
uplo |
Must be 'U' or 'L'. Indicates whether the upper or lower triangular part of A is stored in the array ab, and how A is factored: If uplo = 'U', the upper triangle of A is stored. If uplo = 'L', the lower triangle of A is stored. |
n |
The order of matrix A; n≥ 0. |
kd |
The number of superdiagonals or subdiagonals in the matrix A; kd≥ 0. |
ab |
Array, size max(1, ldab*n). The array ab contains either the upper or the lower triangular part of the matrix A (as specified by uplo) in band storage (see Matrix Storage Schemes). |
ldab |
The leading dimension of the array ab. (ldab≥kd + 1) |
Output Parameters
ab |
The upper or lower triangular part of A (in band storage) is overwritten by the Cholesky factor U or L, as specified by uplo. |
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, the leading minor of order i (and therefore the matrix A itself) is not positive-definite, and the factorization could not be completed. This may indicate an error in forming the matrix A.
Application Notes
If uplo = 'U', the computed factor U is the exact factor of a perturbed matrix A + E, where
c(n) is a modest linear function of n, and ε is the machine precision.
A similar estimate holds for uplo = 'L'.
The total number of floating-point operations for real flavors is approximately n(kd+1)2. The number of operations for complex flavors is 4 times greater. All these estimates assume that kd is much less than n.
After calling this routine, you can call the following routines:
to solve A*X = B |
|
to estimate the condition number of A. |