scipy.sparse.linalg.lsmr¶
-
scipy.sparse.linalg.
lsmr
(A, b, damp=0.0, atol=1e-06, btol=1e-06, conlim=100000000.0, maxiter=None, show=False)[source]¶ Iterative solver for least-squares problems.
lsmr solves the system of linear equations
Ax = b
. If the system is inconsistent, it solves the least-squares problemmin ||b - Ax||_2
. A is a rectangular matrix of dimension m-by-n, where all cases are allowed: m = n, m > n, or m < n. B is a vector of length m. The matrix A may be dense or sparse (usually sparse).Parameters: A : {matrix, sparse matrix, ndarray, LinearOperator}
Matrix A in the linear system.
b : array_like, shape (m,)
Vector b in the linear system.
damp : float
Damping factor for regularized least-squares.
lsmr
solves the regularized least-squares problem:min ||(b) - ( A )x|| ||(0) (damp*I) ||_2
where damp is a scalar. If damp is None or 0, the system is solved without regularization.
atol, btol : float, optional
Stopping tolerances.
lsmr
continues iterations until a certain backward error estimate is smaller than some quantity depending on atol and btol. Letr = b - Ax
be the residual vector for the current approximate solutionx
. IfAx = b
seems to be consistent,lsmr
terminates whennorm(r) <= atol * norm(A) * norm(x) + btol * norm(b)
. Otherwise, lsmr terminates whennorm(A^{T} r) <= atol * norm(A) * norm(r)
. If both tolerances are 1.0e-6 (say), the finalnorm(r)
should be accurate to about 6 digits. (The final x will usually have fewer correct digits, depending oncond(A)
and the size of LAMBDA.) If atol or btol is None, a default value of 1.0e-6 will be used. Ideally, they should be estimates of the relative error in the entries of A and B respectively. For example, if the entries of A have 7 correct digits, set atol = 1e-7. This prevents the algorithm from doing unnecessary work beyond the uncertainty of the input data.conlim : float, optional
lsmr
terminates if an estimate ofcond(A)
exceeds conlim. For compatible systemsAx = b
, conlim could be as large as 1.0e+12 (say). For least-squares problems, conlim should be less than 1.0e+8. If conlim is None, the default value is 1e+8. Maximum precision can be obtained by settingatol = btol = conlim = 0
, but the number of iterations may then be excessive.maxiter : int, optional
lsmr
terminates if the number of iterations reaches maxiter. The default ismaxiter = min(m, n)
. For ill-conditioned systems, a larger value of maxiter may be needed.show : bool, optional
Print iterations logs if
show=True
.Returns: x : ndarray of float
Least-square solution returned.
istop : int
istop gives the reason for stopping:
istop = 0 means x=0 is a solution. = 1 means x is an approximate solution to A*x = B, according to atol and btol. = 2 means x approximately solves the least-squares problem according to atol. = 3 means COND(A) seems to be greater than CONLIM. = 4 is the same as 1 with atol = btol = eps (machine precision) = 5 is the same as 2 with atol = eps. = 6 is the same as 3 with CONLIM = 1/eps. = 7 means ITN reached maxiter before the other stopping conditions were satisfied.
itn : int
Number of iterations used.
normr : float
norm(b-Ax)
normar : float
norm(A^T (b - Ax))
norma : float
norm(A)
conda : float
Condition number of A.
normx : float
norm(x)
Notes
New in version 0.11.0.
References
[R339] D. C.-L. Fong and M. A. Saunders, “LSMR: An iterative algorithm for sparse least-squares problems”, SIAM J. Sci. Comput., vol. 33, pp. 2950-2971, 2011. http://arxiv.org/abs/1006.0758 [R340] LSMR Software, http://web.stanford.edu/group/SOL/software/lsmr/