Determinant

Main header file

Subroutines

Though there are many different possible definitions of the determinant of a matrix \(A \in \mathbb{F}^{n \times n}\), the simplest one is in terms of the product of the eigenvalues (including multiplicity):

\[\mbox{det}(A) = \prod_{i=0}^{n-1} \lambda_i.\]

General

Since \(\mbox{det}(AB)=\mbox{det}(A)\mbox{det}(B)\), we can compute the determinant of an arbitrary matrix in \(\mathcal{O}(n^3)\) work by computing its LU decomposition (with partial pivoting), \(PA=LU\), recognizing that \(\mbox{det}(P)=\pm 1\) (the signature of the permutation), and computing

\[\mbox{det}(A) = \mbox{det}(P)\mbox{det}(L)\mbox{det}(U) = \mbox{det}(P) \prod_{i=0}^{n-1} \upsilon_{i,i} = \pm \prod_{i=0}^{n-1} \upsilon_{i,i},\]

where \(\upsilon_{i,i}\) is the i’th diagonal entry of \(U\).

C++ API

F Determinant(const Matrix<F> &A)
F Determinant(const AbstractDistMatrix<F> &A)
F Determinant(Matrix<F> &A, bool canOverwrite = false)
F Determinant(AbstractDistMatrix<F> &A, bool canOverwrite = false)

Return the determinant of the (fully populated) square matrix A. Some of the variants allow for overwriting the input matrix in order to avoid forming another temporary matrix.

SafeProduct<F> SafeDeterminant(const Matrix<F> &A)
SafeProduct<F> SafeDeterminant(const AbstractDistMatrix<F> &A)
SafeProduct<F> SafeDeterminant(Matrix<F> &A, bool canOverwrite = false)
SafeProduct<F> SafeDeterminant(AbstractDistMatrix<F> &A, bool canOverwrite = false)

Return the determinant of the square matrix A in an expanded form which is less likely to over/under-flow.

template<>
type SafeProduct<F>

Represents the product of n values as \(\rho \exp(\kappa n)\), where \(|\rho| \le 1\) and \(\kappa \in \mathbb{R}\).

F rho

For nonzero values, rho is the modulus and lies on the unit circle; in order to represent a value that is precisely zero, rho is set to zero.

Base<F> kappa

kappa can be an arbitrary real number.

Int n

The number of values in the product.

C API

ElError ElDeterminant_s(ElConstMatrix_s A, float* det)
ElError ElDeterminant_d(ElConstMatrix_d A, double* det)
ElError ElDeterminant_c(ElConstMatrix_c A, complex_float* det)
ElError ElDeterminant_z(ElConstMatrix_z A, complex_double* det)
ElError ElDeterminantDist_s(ElConstDistMatrix_s A, float* det)
ElError ElDeterminantDist_d(ElConstDistMatrix_d A, double* det)
ElError ElDeterminantDist_c(ElConstDistMatrix_c A, complex_float* det)
ElError ElDeterminantDist_z(ElConstDistMatrix_z A, complex_double* det)

Return the determinant of the (fully populated) square matrix A.

ElError ElSafeDeterminant_s(ElConstMatrix_s A, ElSafeProduct_s* det)
ElError ElSafeDeterminant_d(ElConstMatrix_d A, ElSafeProduct_d* det)
ElError ElSafeDeterminant_c(ElConstMatrix_c A, ElSafeProduct_c* det)
ElError ElSafeDeterminant_z(ElConstMatrix_z A, ElSafeProduct_z* det)
ElError ElSafeDeterminantDist_s(ElConstDistMatrix_s A, ElSafeProduct_s* det)
ElError ElSafeDeterminantDist_d(ElConstDistMatrix_d A, ElSafeProduct_d* det)
ElError ElSafeDeterminantDist_c(ElConstDistMatrix_c A, ElSafeProduct_c* det)
ElError ElSafeDeterminantDist_z(ElConstDistMatrix_z A, ElSafeProduct_z* det)

Return the determinant of the (fully populated) square matrix A in an expanded form which helps prevent under/overflow.

HPD

A version of the above determinant specialized for Hermitian positive-definite matrices (which will therefore have all positive eigenvalues and a positive determinant).

C++ API

Base<F> HPDDeterminant(UpperOrLower uplo, const Matrix<F> &A)
Base<F> HPDDeterminant(UpperOrLower uplo, const AbstractDistMatrix<F> &A)
Base<F> HPDDeterminant(UpperOrLower uplo, Matrix<F> &A, bool canOverwrite = false)
Base<F> HPDDeterminant(UpperOrLower uplo, AbstractDistMatrix<F> &A, bool canOverwrite = false)

Return the determinant of the (fully populated) Hermitian positive-definite matrix A. Some of the variants allow for overwriting the input matrix in order to avoid forming another temporary matrix.

SafeProduct<F> SafeHPDDeterminant(UpperOrLower uplo, const Matrix<F> &A)
SafeProduct<F> SafeHPDDeterminant(UpperOrLower uplo, const AbstractDistMatrix<F> &A)
SafeProduct<F> SafeHPDDeterminant(UpperOrLower uplo, Matrix<F> &A, bool canOverwrite = false)
SafeProduct<F> SafeHPDDeterminant(UpperOrLower uplo, AbstractDistMatrix<F> &A, bool canOverwrite = false)

Return the determinant of the Hermitian positive-definite matrix A in an expanded form which is less likely to over/under-flow.

C API

ElError ElHPDDeterminant_s(ElUpperOrLower uplo, ElConstMatrix_s A, float* det)
ElError ElHPDDeterminant_d(ElUpperOrLower uplo, ElConstMatrix_d A, double* det)
ElError ElHPDDeterminant_c(ElUpperOrLower uplo, ElConstMatrix_c A, float* det)
ElError ElHPDDeterminant_z(ElUpperOrLower uplo, ElConstMatrix_z A, double* det)
ElError ElHPDDeterminantDist_s(ElUpperOrLower uplo, ElConstDistMatrix_s A, float* det)
ElError ElHPDDeterminantDist_d(ElUpperOrLower uplo, ElConstDistMatrix_d A, double* det)
ElError ElHPDDeterminantDist_c(ElUpperOrLower uplo, ElConstDistMatrix_c A, float* det)
ElError ElHPDDeterminantDist_z(ElUpperOrLower uplo, ElConstDistMatrix_z A, double* det)

Return the determinant of the (fully populated) Hermitian positive-definite matrix A.

ElError ElHPDSafeDeterminant_s(ElUpperOrLower uplo, ElConstMatrix_s A, ElSafeProduct_s* det)
ElError ElHPDSafeDeterminant_d(ElUpperOrLower uplo, ElConstMatrix_d A, ElSafeProduct_d* det)
ElError ElHPDSafeDeterminant_c(ElUpperOrLower uplo, ElConstMatrix_c A, ElSafeProduct_s* det)
ElError ElHPDSafeDeterminant_z(ElUpperOrLower uplo, ElConstMatrix_z A, ElSafeProduct_d* det)
ElError ElHPDSafeDeterminantDist_s(ElUpperOrLower uplo, ElConstDistMatrix_s A, ElSafeProduct_s* det)
ElError ElHPDSafeDeterminantDist_d(ElUpperOrLower uplo, ElConstDistMatrix_d A, ElSafeProduct_d* det)
ElError ElHPDSafeDeterminantDist_c(ElUpperOrLower uplo, ElConstDistMatrix_c A, ElSafeProduct_s* det)
ElError ElHPDSafeDeterminantDist_z(ElUpperOrLower uplo, ElConstDistMatrix_z A, ElSafeProduct_d* det)

Return the determinant of the Hermitian positive-definite matrix A in an expanded form which is less likely to over/under-flow.