Visible to Intel only — GUID: GUID-07B01968-68F5-4540-B5C4-E3F9FC5A3E76
Visible to Intel only — GUID: GUID-07B01968-68F5-4540-B5C4-E3F9FC5A3E76
?pftrf
Computes the Cholesky factorization of a symmetric (Hermitian) positive-definite matrix using the Rectangular Full Packed (RFP) format .
Syntax
call spftrf( transr, uplo, n, a, info )
call dpftrf( transr, uplo, n, a, info )
call cpftrf( transr, uplo, n, a, info )
call zpftrf( transr, uplo, n, a, info )
Include Files
- mkl.fi, lapack.f90
Description
The routine forms the Cholesky factorization of a symmetric positive-definite or, for complex data, a Hermitian positive-definite 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.
The matrix A is in the Rectangular Full Packed (RFP) format. For the description of the RFP format, see Matrix Storage Schemes.
This is the block version of the algorithm, calling Level 3 BLAS.
Input Parameters
transr |
CHARACTER*1. Must be 'N', 'T' (for real data) or 'C' (for complex data). If transr = 'N', the Normal transr of RFP A is stored. If transr = 'T', the Transpose transr of RFP A is stored. If transr = 'C', the Conjugate-Transpose transr of RFP A is stored. |
uplo |
CHARACTER*1. Must be 'U' or 'L'. Indicates whether the upper or lower triangular part of A is stored: If uplo = 'U', the array a stores the upper triangular part of the matrix A. If uplo = 'L', the array a stores the lower triangular part of the matrix A. |
n |
INTEGER. The order of the matrix A; n≥ 0. |
a |
REAL for spftrf DOUBLE PRECISION for dpftrf COMPLEX for cpftrf DOUBLE COMPLEX for zpftrf. Array, size (n*(n+1)/2). The array a contains the matrix A in the RFP format. |
Output Parameters
a |
a is overwritten by the Cholesky factor U or L, as specified by uplo and trans. |
info |
INTEGER. If info=0, the execution is successful. If info = -i, the i-th parameter 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. |