scipy.sparse.linalg.

splu#

scipy.sparse.linalg.splu(A, permc_spec=None, diag_pivot_thresh=None, relax=None, panel_size=None, options=None)[source]#

Compute the LU decomposition of a sparse, square matrix.

Parameters:
Asparse array or matrix

Sparse array to factorize. Most efficient when provided in CSC format. Other formats will be converted to CSC before factorization.

permc_specstr, optional

How to permute the columns of the matrix for sparsity preservation. (default: ‘COLAMD’)

  • NATURAL: natural ordering.

  • MMD_ATA: minimum degree ordering on the structure of A^T A.

  • MMD_AT_PLUS_A: minimum degree ordering on the structure of A^T+A.

  • COLAMD: approximate minimum degree column ordering

diag_pivot_threshfloat, optional

Threshold used for a diagonal entry to be an acceptable pivot. See SuperLU user’s guide for details [1]

relaxint, optional

Expert option for customizing the degree of relaxing supernodes. See SuperLU user’s guide for details [1]

panel_sizeint, optional

Expert option for customizing the panel size. See SuperLU user’s guide for details [1]

optionsdict, optional

Dictionary containing additional expert options to SuperLU. See SuperLU user guide [1] (section 2.4 on the ‘Options’ argument) for more details. For example, you can specify options=dict(Equil=False, IterRefine='SINGLE')) to turn equilibration off and perform a single iterative refinement.

Returns:
invAscipy.sparse.linalg.SuperLU

Object, which has a solve method.

See also

spilu

incomplete LU decomposition

Notes

When a real array is factorized and the returned SuperLU object’s solve() method is used with complex arguments an error is generated. Instead, cast the initial array to complex and then factorize.

This function uses the SuperLU library.

References

Examples

>>> import numpy as np
>>> from scipy.sparse import csc_array
>>> from scipy.sparse.linalg import splu
>>> A = csc_array([[1., 0., 0.], [5., 0., 2.], [0., -1., 0.]], dtype=float)
>>> B = splu(A)
>>> x = np.array([1., 2., 3.], dtype=float)
>>> B.solve(x)
array([ 1. , -3. , -1.5])
>>> A.dot(B.solve(x))
array([ 1.,  2.,  3.])
>>> B.solve(A.dot(x))
array([ 1.,  2.,  3.])