Processing math: 100%
fortran-lapack
All Classes Namespaces Files Functions Variables Pages
la_cholesky::chol Interface Reference

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)
 

Detailed Description

Computes the Cholesky factorization A = L \cdot L^T , or A = U^T \cdot U .

Summary

Pure function interface for computing the Cholesky triangular factors.

Description

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.

Note
The solution is based on LAPACK's POTRF methods.
Parameters
[in]aThe 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.
Returns
c The output matrix of size [n, n] , containing the Cholesky factor.

Member Function/Subroutine Documentation

◆ la_c_cholesky_fun()

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 )
Parameters
[in]aInput matrix a[n,n]
Returns
Output matrix with Cholesky factors c[n,n]
Parameters
[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

◆ la_d_cholesky_fun()

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 )
Parameters
[in]aInput matrix a[n,n]
Returns
Output matrix with Cholesky factors c[n,n]
Parameters
[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

◆ la_q_cholesky_fun()

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 )
Parameters
[in]aInput matrix a[n,n]
Returns
Output matrix with Cholesky factors c[n,n]
Parameters
[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

◆ la_s_cholesky_fun()

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 )
Parameters
[in]aInput matrix a[n,n]
Returns
Output matrix with Cholesky factors c[n,n]
Parameters
[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

◆ la_w_cholesky_fun()

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 )
Parameters
[in]aInput matrix a[n,n]
Returns
Output matrix with Cholesky factors c[n,n]
Parameters
[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

◆ la_z_cholesky_fun()

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 )
Parameters
[in]aInput matrix a[n,n]
Returns
Output matrix with Cholesky factors c[n,n]
Parameters
[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

The documentation for this interface was generated from the following file: