fortran-lapack
|
Compute the pseudo-inverse of a matrix. More...
Public Member Functions | |
real(sp) function, dimension(size(a, 2, kind=ilp), size(a, 1, kind=ilp)) | la_pseudoinverse_s (a, rtol, err) |
real(dp) function, dimension(size(a, 2, kind=ilp), size(a, 1, kind=ilp)) | la_pseudoinverse_d (a, rtol, err) |
real(qp) function, dimension(size(a, 2, kind=ilp), size(a, 1, kind=ilp)) | la_pseudoinverse_q (a, rtol, err) |
complex(sp) function, dimension(size(a, 2, kind=ilp), size(a, 1, kind=ilp)) | la_pseudoinverse_c (a, rtol, err) |
complex(dp) function, dimension(size(a, 2, kind=ilp), size(a, 1, kind=ilp)) | la_pseudoinverse_z (a, rtol, err) |
complex(qp) function, dimension(size(a, 2, kind=ilp), size(a, 1, kind=ilp)) | la_pseudoinverse_w (a, rtol, err) |
Compute the pseudo-inverse of a matrix.
This function computes the Moore-Penrose pseudo-inverse of a real or complex matrix A . The pseudo-inverse is computed using Singular Value Decomposition (SVD):
A^+ = V \Sigma^+ U^T
where U and V are unitary matrices, and \Sigma^+ is the pseudo-inverse of the singular values.
[in] | A | The input matrix of size [m, n] . |
[in] | rtol | (Optional) Relative tolerance for singular value truncation. If not provided, a default value is used. |
[out] | err | (Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution. |
rtol
is too large, important singular values may be discarded, leading to inaccurate results.Compute the pseudo-inverse of a matrix (subroutine version).
This subroutine computes the Moore-Penrose pseudo-inverse of a real or complex matrix A , storing the result in a pre-allocated output matrix A^+ . The computation is based on Singular Value Decomposition (SVD):
A^+ = V \Sigma^+ U^T
where U and V are unitary matrices, and \Sigma^+ is the pseudo-inverse of the singular values.
[in,out] | A | The input matrix of size [m, n] . Its contents may be modified. |
[out] | pinva | The output pseudo-inverse matrix of size [n, m] . |
[in] | rtol | (Optional) Relative tolerance for singular value truncation. |
[out] | err | (Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution. |
pinva
is already allocated and avoids memory allocation inside the routine. A
may be modified during computation.Compute the pseudo-inverse of a matrix using the .pinv.
operator.
This operator computes the Moore-Penrose pseudo-inverse of a real or complex matrix A using Singular Value Decomposition (SVD):
A^+ = V \Sigma^+ U^T
where U and V are unitary matrices, and \Sigma^+ is the pseudo-inverse of the singular values.
[in] | A | The input matrix of size [m, n] . |
pinv(A)
, allowing expressions such as: X = .pinv.A Compute the pseudo-inverse of a matrix.
This function computes the Moore-Penrose pseudo-inverse of a real or complex matrix A . The pseudo-inverse is computed using Singular Value Decomposition (SVD):
A^+ = V \Sigma^+ U^T
where U and V are unitary matrices, and \Sigma^+ is the pseudo-inverse of the singular values.
[in] | A | The input matrix of size [m, n] . |
[in] | rtol | (Optional) Relative tolerance for singular value truncation. If not provided, a default value is used. |
[out] | err | (Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution. |
*GESVD
or *GESDD
). rtol
is too large, important singular values may be discarded, leading to inaccurate results. complex(sp) function, dimension(size(a,2,kind=ilp),size(a,1,kind=ilp)) la_pseudoinverse::pinv::la_pseudoinverse_c | ( | complex(sp), dimension(:,:), intent(in), target | a, |
real(sp), intent(in), optional | rtol, | ||
type(la_state), intent(out), optional | err ) |
[in] | a | Input matrix a[m,n] |
[in] | rtol | [optional] .... |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
real(dp) function, dimension(size(a,2,kind=ilp),size(a,1,kind=ilp)) la_pseudoinverse::pinv::la_pseudoinverse_d | ( | real(dp), dimension(:,:), intent(in), target | a, |
real(dp), intent(in), optional | rtol, | ||
type(la_state), intent(out), optional | err ) |
[in] | a | Input matrix a[m,n] |
[in] | rtol | [optional] .... |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
real(qp) function, dimension(size(a,2,kind=ilp),size(a,1,kind=ilp)) la_pseudoinverse::pinv::la_pseudoinverse_q | ( | real(qp), dimension(:,:), intent(in), target | a, |
real(qp), intent(in), optional | rtol, | ||
type(la_state), intent(out), optional | err ) |
[in] | a | Input matrix a[m,n] |
[in] | rtol | [optional] .... |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
real(sp) function, dimension(size(a,2,kind=ilp),size(a,1,kind=ilp)) la_pseudoinverse::pinv::la_pseudoinverse_s | ( | real(sp), dimension(:,:), intent(in), target | a, |
real(sp), intent(in), optional | rtol, | ||
type(la_state), intent(out), optional | err ) |
[in] | a | Input matrix a[m,n] |
[in] | rtol | [optional] .... |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
complex(qp) function, dimension(size(a,2,kind=ilp),size(a,1,kind=ilp)) la_pseudoinverse::pinv::la_pseudoinverse_w | ( | complex(qp), dimension(:,:), intent(in), target | a, |
real(qp), intent(in), optional | rtol, | ||
type(la_state), intent(out), optional | err ) |
[in] | a | Input matrix a[m,n] |
[in] | rtol | [optional] .... |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
complex(dp) function, dimension(size(a,2,kind=ilp),size(a,1,kind=ilp)) la_pseudoinverse::pinv::la_pseudoinverse_z | ( | complex(dp), dimension(:,:), intent(in), target | a, |
real(dp), intent(in), optional | rtol, | ||
type(la_state), intent(out), optional | err ) |
[in] | a | Input matrix a[m,n] |
[in] | rtol | [optional] .... |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |