spqr¶
- sksparse.spqr.spqr(A, *, mode='full', order=None, tol=None)[source]¶
Compute the QR factorization.
This function computes the QR factorization of a sparse matrix \(A\) such that
\[Q R = A E\]where \(Q\) is an orthogonal matrix and \(R\) is an upper-triangular matrix. \(E\) is a column permutation matrix that reduces fill-in during the factorization.
- Parameters:
A ((M, N) array_like or sparse array) – An array convertible to a sparse matrix.
mode ({‘full’, ‘r’, ‘economic’, ‘householder’}, optional) – The mode of the returned Q and R matrices. Options are:
full:Qis size(M, M),Ris size(M, N).economic:Qis size(M, K),Ris size(K, N), whereK = min(M, N).r: Only return the upper-triangular matrixR.householder: Return the Householder vectors and coefficients used to buildQ. This option is similar tomode='raw'inscipy.linalg.qr().
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 ofamd,colamd,metisand pick the best,cholmod: Same asbest,bestamd: tryamdandcolamdand pick the best.
tol (float, optional) – If the 2-norm of a column in
Ais less thantol, that column is considered to be a zero column. Iftol = 0, no columns are treated as zero. IfNone, the default tolerance is used. The default istol =\(20 \epsilon (M + N) \sqrt{\max{\mathrm{diag}(A^{\top} A)}}\), where \(\epsilon\) is the machine precision.
- Returns:
Q (csc_array) – The orthogonal matrix \(Q\). Shape (M, M) or (M, K) if
mode='economic'. Not returned ifmode='r'. Replaced bySPQRHouseholderifmode='householder'.R (csc_array) – The upper-triangular matrix \(R\). Shape (M, N) or (K, N) if
mode in ['economic', 'householder'], where K = min(M, N).P (ndarray of int) – The permutation vector of shape (N,).
See also
Notes
This function is part of an interface to the SuiteSparse SPQR library [1].
Added in version 0.5.0.
References