spqr_factor

sksparse.spqr.spqr_factor(A, *, use_singletons=False, order=None, tol=None)[source]

Compute the SPQR factorization of a sparse matrix.

Compute the numeric factorization of the matrix and determine a fill-reducing ordering such that:

\[Q R = A E\]

where \(E\) is a column permutation matrix, \(Q\) is an orthogonal matrix, and \(R\) is an upper-triangular matrix.

This function returns a SPQRFactor object that contains the SPQR factorization of the input matrix. It is not currently possible to extract the individual factors \(Q\) and \(R\) explicitly, but the object provides methods to reuse the factorization to solve linear systems or multiply by \(Q\).

Parameters:
  • A ((M, N) array_like or sparse array) – An array convertible to a sparse matrix.

  • use_singletons (bool, optional) – If True, directly compute the numeric factorization to exploit singleton rows. Otherwise, only perform symbolic analysis. Default is False, so that the factor can be reused efficiently for multiple numeric factorizations.

  • order (str, optional) – The column ordering strategy to use. Let \(S\) be the matrix \(A\) with singleton rows/columns removed, the ordering options are:

    • default: COLAMD(S),

    • fixed: identity permutation (i.e. no singletons removed),

    • natural: singletons removed, but no fill-reducing ordering applied,

    • colamd: COLAMD(S),

    • amd: AMD(\(S^{\top} S\)),

    • metis: METIS(\(S^{\top} S\)),

    • best: try all of amd, colamd, metis and pick the best,

    • cholmod: Same as best,

    • bestamd: try amd and colamd and pick the best.

  • tol (float, optional) – If the 2-norm of a column in A is less than tol, that column is considered to be a zero column. If tol = 0, no columns are treated as zero. If None, the default tolerance is used. The default is tol = \(20 \epsilon (M + N) \sqrt{\max{\mathrm{diag}(A^{\top} A)}}\), where \(\epsilon\) is the machine precision.

Returns:

SPQRFactor – The SPQR factorization of the input matrix.

Notes

This function is part of an interface to the SuiteSparse SPQR library [1].

Added in version 0.5.0.

References