fortran-lapack
|
Computes the Cholesky factorization A = L \cdot L^T , or A = U^T \cdot U . More...
Public Member Functions | |
pure real(sp) function, dimension(size(a, 1), size(a, 2)) | la_s_cholesky_fun (a, lower, other_zeroed) |
pure real(dp) function, dimension(size(a, 1), size(a, 2)) | la_d_cholesky_fun (a, lower, other_zeroed) |
pure real(qp) function, dimension(size(a, 1), size(a, 2)) | la_q_cholesky_fun (a, lower, other_zeroed) |
pure complex(sp) function, dimension(size(a, 1), size(a, 2)) | la_c_cholesky_fun (a, lower, other_zeroed) |
pure complex(dp) function, dimension(size(a, 1), size(a, 2)) | la_z_cholesky_fun (a, lower, other_zeroed) |
pure complex(qp) function, dimension(size(a, 1), size(a, 2)) | la_w_cholesky_fun (a, lower, other_zeroed) |
Computes the Cholesky factorization A = L \cdot L^T , or A = U^T \cdot U .
Pure function interface for computing the Cholesky triangular factors.
This interface provides methods for computing the lower- or upper-triangular matrix from the Cholesky factorization of a real
symmetric or complex
Hermitian matrix. Supported data types include real
and complex
.
The function returns the Cholesky factor as a separate matrix. The lower
logical
flag allows selection between upper or lower factorization, and the other_zeroed
optional logical
flag determines whether the unused part of the triangular matrix should be zeroed.
[in] | a | The input matrix of size [n, n] . |
[in] | lower | (Optional) Logical flag indicating whether the lower or upper triangular factor is required. Default: lower triangular factor. |
[in] | other_zeroed | (Optional) Logical flag indicating whether the unused half of the matrix should be zeroed. Default: true. |
pure complex(sp) function, dimension(size(a,1),size(a,2)) la_cholesky::chol::la_c_cholesky_fun | ( | complex(sp), dimension(:,:), intent(in) | a, |
logical(lk), intent(in), optional | lower, | ||
logical(lk), intent(in), optional | other_zeroed ) |
[in] | a | Input matrix a[n,n] |
[in] | lower | [optional] is the lower or upper triangular factor required? Default = lower |
[in] | other_zeroed | [optional] should the unused half of the return matrix be zeroed out? Default: yes |
pure real(dp) function, dimension(size(a,1),size(a,2)) la_cholesky::chol::la_d_cholesky_fun | ( | real(dp), dimension(:,:), intent(in) | a, |
logical(lk), intent(in), optional | lower, | ||
logical(lk), intent(in), optional | other_zeroed ) |
[in] | a | Input matrix a[n,n] |
[in] | lower | [optional] is the lower or upper triangular factor required? Default = lower |
[in] | other_zeroed | [optional] should the unused half of the return matrix be zeroed out? Default: yes |
pure real(qp) function, dimension(size(a,1),size(a,2)) la_cholesky::chol::la_q_cholesky_fun | ( | real(qp), dimension(:,:), intent(in) | a, |
logical(lk), intent(in), optional | lower, | ||
logical(lk), intent(in), optional | other_zeroed ) |
[in] | a | Input matrix a[n,n] |
[in] | lower | [optional] is the lower or upper triangular factor required? Default = lower |
[in] | other_zeroed | [optional] should the unused half of the return matrix be zeroed out? Default: yes |
pure real(sp) function, dimension(size(a,1),size(a,2)) la_cholesky::chol::la_s_cholesky_fun | ( | real(sp), dimension(:,:), intent(in) | a, |
logical(lk), intent(in), optional | lower, | ||
logical(lk), intent(in), optional | other_zeroed ) |
[in] | a | Input matrix a[n,n] |
[in] | lower | [optional] is the lower or upper triangular factor required? Default = lower |
[in] | other_zeroed | [optional] should the unused half of the return matrix be zeroed out? Default: yes |
pure complex(qp) function, dimension(size(a,1),size(a,2)) la_cholesky::chol::la_w_cholesky_fun | ( | complex(qp), dimension(:,:), intent(in) | a, |
logical(lk), intent(in), optional | lower, | ||
logical(lk), intent(in), optional | other_zeroed ) |
[in] | a | Input matrix a[n,n] |
[in] | lower | [optional] is the lower or upper triangular factor required? Default = lower |
[in] | other_zeroed | [optional] should the unused half of the return matrix be zeroed out? Default: yes |
pure complex(dp) function, dimension(size(a,1),size(a,2)) la_cholesky::chol::la_z_cholesky_fun | ( | complex(dp), dimension(:,:), intent(in) | a, |
logical(lk), intent(in), optional | lower, | ||
logical(lk), intent(in), optional | other_zeroed ) |
[in] | a | Input matrix a[n,n] |
[in] | lower | [optional] is the lower or upper triangular factor required? Default = lower |
[in] | other_zeroed | [optional] should the unused half of the return matrix be zeroed out? Default: yes |