Pseudoinverse¶
General¶
Computes the pseudoinverse of a general matrix through computing its SVD, modifying the singular values with the function
where \(\epsilon\) is the relative machine precision, \(n\) is the height of \(A\), and \(\| A \|_2\) is the maximum singular value. If a nonzero value for tolerance was specified, it is used instead of \(\epsilon n \| A \|_2\).
C++ API¶
-
void
Pseudoinverse
(AbstractDistMatrix<F> &A, Base<F> tolerance = 0)¶
C API¶
Hermitian¶
Computes the pseudoinverse of a Hermitian matrix through a customized version
of HermitianFunction()
which used the eigenvalue mapping
function
where \(\epsilon\) is the relative machine precision, \(n\) is the height of \(A\), and \(\| A \|_2\) can be computed as the maximum absolute value of the eigenvalues of \(A\). If a nonzero value for tolerance is specified, it is used instead of \(\epsilon n \| A \|_2\).
C++ API¶
-
void
HermitianPseudoinverse
(UpperOrLower uplo, Matrix<F> &A, Base<F> tolerance = 0)¶
-
void
HermitianPseudoinverse
(UpperOrLower uplo, AbstractDistMatrix<F> &A, Base<F> tolerance = 0)¶
C API¶
-
ElError
ElHermitianPseudoinverse_s
(ElUpperOrLower uplo, ElMatrix_s A)¶
-
ElError
ElHermitianPseudoinverse_d
(ElUpperOrLower uplo, ElMatrix_d A)¶
-
ElError
ElHermitianPseudoinverse_c
(ElUpperOrLower uplo, ElMatrix_c A)¶
-
ElError
ElHermitianPseudoinverse_z
(ElUpperOrLower uplo, ElMatrix_z A)¶
-
ElError
ElHermitianPseudoinverseDist_s
(ElUpperOrLower uplo, ElDistMatrix_s A)¶
-
ElError
ElHermitianPseudoinverseDist_d
(ElUpperOrLower uplo, ElDistMatrix_d A)¶
-
ElError
ElHermitianPseudoinverseDist_c
(ElUpperOrLower uplo, ElDistMatrix_c A)¶
-
ElError
ElHermitianPseudoinverseDist_z
(ElUpperOrLower uplo, ElDistMatrix_z A)¶