Computes the Cholesky factorization , or . (Specification)
Pure subroutine 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 factorization is computed in-place if only one matrix argument is present; or returned into
a second matrix argument, if present. The lower
logical
flag allows to select between upper or
lower factorization; the other_zeroed
optional logical
flag allows to choose whether the unused
part of the triangular matrix should be filled with zeroes.
Note
The solution is based on LAPACK's *POTRF
methods.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=sp), | intent(in) | :: | a(:,:) |
Input matrix a[n,n] |
||
complex(kind=sp), | intent(out) | :: | c(:,:) |
Output matrix with Cholesky factors c[n,n] |
||
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=sp), | intent(inout), | target | :: | a(:,:) |
Input matrix a[m,n] |
|
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | a(:,:) |
Input matrix a[n,n] |
||
real(kind=dp), | intent(out) | :: | c(:,:) |
Output matrix with Cholesky factors c[n,n] |
||
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(inout), | target | :: | a(:,:) |
Input matrix a[m,n] |
|
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=sp), | intent(in) | :: | a(:,:) |
Input matrix a[n,n] |
||
real(kind=sp), | intent(out) | :: | c(:,:) |
Output matrix with Cholesky factors c[n,n] |
||
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=sp), | intent(inout), | target | :: | a(:,:) |
Input matrix a[m,n] |
|
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=dp), | intent(in) | :: | a(:,:) |
Input matrix a[n,n] |
||
complex(kind=dp), | intent(out) | :: | c(:,:) |
Output matrix with Cholesky factors c[n,n] |
||
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex(kind=dp), | intent(inout), | target | :: | a(:,:) |
Input matrix a[m,n] |
|
logical(kind=lk), | intent(in), | optional | :: | lower |
[optional] is the lower or upper triangular factor required? Default = lower |
|
logical(kind=lk), | intent(in), | optional | :: | other_zeroed |
[optional] should the unused half of the return matrix be zeroed out? Default: yes |
|
type(linalg_state_type), | intent(out), | optional | :: | err |
[optional] state return flag. On error if not requested, the code will stop |