spilu#
- scipy.sparse.linalg.spilu(A, drop_tol=None, fill_factor=None, drop_rule=None, permc_spec=None, diag_pivot_thresh=None, relax=None, panel_size=None, options=None)[source]#
Compute an incomplete LU decomposition for a sparse, square matrix.
The resulting object is an approximation to the inverse of A.
- Parameters:
- A(N, N) array_like
Sparse array to factorize. Most efficient when provided in CSC format. Other formats will be converted to CSC before factorization.
- drop_tolfloat, optional
Drop tolerance (0 <= tol <= 1) for an incomplete LU decomposition. (default: 1e-4)
- fill_factorfloat, optional
Specifies the fill ratio upper bound (>= 1.0) for ILU. (default: 10)
- drop_rulestr, optional
Comma-separated string of drop rules to use. Available rules:
basic
,prows
,column
,area
,secondary
,dynamic
,interp
. (Default:basic,area
)See SuperLU documentation for details.
- Remaining other options
Same as for
splu
- Returns:
- invA_approxscipy.sparse.linalg.SuperLU
Object, which has a
solve
method.
See also
splu
complete 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.To improve the better approximation to the inverse, you may need to increase fill_factor AND decrease drop_tol.
This function uses the SuperLU library.
Examples
>>> import numpy as np >>> from scipy.sparse import csc_array >>> from scipy.sparse.linalg import spilu >>> A = csc_array([[1., 0., 0.], [5., 0., 2.], [0., -1., 0.]], dtype=float) >>> B = spilu(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.])