lsqr.py 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568
  1. """Sparse Equations and Least Squares.
  2. The original Fortran code was written by C. C. Paige and M. A. Saunders as
  3. described in
  4. C. C. Paige and M. A. Saunders, LSQR: An algorithm for sparse linear
  5. equations and sparse least squares, TOMS 8(1), 43--71 (1982).
  6. C. C. Paige and M. A. Saunders, Algorithm 583; LSQR: Sparse linear
  7. equations and least-squares problems, TOMS 8(2), 195--209 (1982).
  8. It is licensed under the following BSD license:
  9. Copyright (c) 2006, Systems Optimization Laboratory
  10. All rights reserved.
  11. Redistribution and use in source and binary forms, with or without
  12. modification, are permitted provided that the following conditions are
  13. met:
  14. * Redistributions of source code must retain the above copyright
  15. notice, this list of conditions and the following disclaimer.
  16. * Redistributions in binary form must reproduce the above
  17. copyright notice, this list of conditions and the following
  18. disclaimer in the documentation and/or other materials provided
  19. with the distribution.
  20. * Neither the name of Stanford University nor the names of its
  21. contributors may be used to endorse or promote products derived
  22. from this software without specific prior written permission.
  23. THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
  24. "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
  25. LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
  26. A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
  27. OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
  28. SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
  29. LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
  30. DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
  31. THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
  32. (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
  33. OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
  34. The Fortran code was translated to Python for use in CVXOPT by Jeffery
  35. Kline with contributions by Mridul Aanjaneya and Bob Myhill.
  36. Adapted for SciPy by Stefan van der Walt.
  37. """
  38. from __future__ import division, print_function, absolute_import
  39. __all__ = ['lsqr']
  40. import numpy as np
  41. from math import sqrt
  42. from scipy.sparse.linalg.interface import aslinearoperator
  43. eps = np.finfo(np.float64).eps
  44. def _sym_ortho(a, b):
  45. """
  46. Stable implementation of Givens rotation.
  47. Notes
  48. -----
  49. The routine 'SymOrtho' was added for numerical stability. This is
  50. recommended by S.-C. Choi in [1]_. It removes the unpleasant potential of
  51. ``1/eps`` in some important places (see, for example text following
  52. "Compute the next plane rotation Qk" in minres.py).
  53. References
  54. ----------
  55. .. [1] S.-C. Choi, "Iterative Methods for Singular Linear Equations
  56. and Least-Squares Problems", Dissertation,
  57. http://www.stanford.edu/group/SOL/dissertations/sou-cheng-choi-thesis.pdf
  58. """
  59. if b == 0:
  60. return np.sign(a), 0, abs(a)
  61. elif a == 0:
  62. return 0, np.sign(b), abs(b)
  63. elif abs(b) > abs(a):
  64. tau = a / b
  65. s = np.sign(b) / sqrt(1 + tau * tau)
  66. c = s * tau
  67. r = b / s
  68. else:
  69. tau = b / a
  70. c = np.sign(a) / sqrt(1+tau*tau)
  71. s = c * tau
  72. r = a / c
  73. return c, s, r
  74. def lsqr(A, b, damp=0.0, atol=1e-8, btol=1e-8, conlim=1e8,
  75. iter_lim=None, show=False, calc_var=False, x0=None):
  76. """Find the least-squares solution to a large, sparse, linear system
  77. of equations.
  78. The function solves ``Ax = b`` or ``min ||b - Ax||^2`` or
  79. ``min ||Ax - b||^2 + d^2 ||x||^2``.
  80. The matrix A may be square or rectangular (over-determined or
  81. under-determined), and may have any rank.
  82. ::
  83. 1. Unsymmetric equations -- solve A*x = b
  84. 2. Linear least squares -- solve A*x = b
  85. in the least-squares sense
  86. 3. Damped least squares -- solve ( A )*x = ( b )
  87. ( damp*I ) ( 0 )
  88. in the least-squares sense
  89. Parameters
  90. ----------
  91. A : {sparse matrix, ndarray, LinearOperator}
  92. Representation of an m-by-n matrix. It is required that
  93. the linear operator can produce ``Ax`` and ``A^T x``.
  94. b : array_like, shape (m,)
  95. Right-hand side vector ``b``.
  96. damp : float
  97. Damping coefficient.
  98. atol, btol : float, optional
  99. Stopping tolerances. If both are 1.0e-9 (say), the final
  100. residual norm should be accurate to about 9 digits. (The
  101. final x will usually have fewer correct digits, depending on
  102. cond(A) and the size of damp.)
  103. conlim : float, optional
  104. Another stopping tolerance. lsqr terminates if an estimate of
  105. ``cond(A)`` exceeds `conlim`. For compatible systems ``Ax =
  106. b``, `conlim` could be as large as 1.0e+12 (say). For
  107. least-squares problems, conlim should be less than 1.0e+8.
  108. Maximum precision can be obtained by setting ``atol = btol =
  109. conlim = zero``, but the number of iterations may then be
  110. excessive.
  111. iter_lim : int, optional
  112. Explicit limitation on number of iterations (for safety).
  113. show : bool, optional
  114. Display an iteration log.
  115. calc_var : bool, optional
  116. Whether to estimate diagonals of ``(A'A + damp^2*I)^{-1}``.
  117. x0 : array_like, shape (n,), optional
  118. Initial guess of x, if None zeros are used.
  119. .. versionadded:: 1.0.0
  120. Returns
  121. -------
  122. x : ndarray of float
  123. The final solution.
  124. istop : int
  125. Gives the reason for termination.
  126. 1 means x is an approximate solution to Ax = b.
  127. 2 means x approximately solves the least-squares problem.
  128. itn : int
  129. Iteration number upon termination.
  130. r1norm : float
  131. ``norm(r)``, where ``r = b - Ax``.
  132. r2norm : float
  133. ``sqrt( norm(r)^2 + damp^2 * norm(x)^2 )``. Equal to `r1norm` if
  134. ``damp == 0``.
  135. anorm : float
  136. Estimate of Frobenius norm of ``Abar = [[A]; [damp*I]]``.
  137. acond : float
  138. Estimate of ``cond(Abar)``.
  139. arnorm : float
  140. Estimate of ``norm(A'*r - damp^2*x)``.
  141. xnorm : float
  142. ``norm(x)``
  143. var : ndarray of float
  144. If ``calc_var`` is True, estimates all diagonals of
  145. ``(A'A)^{-1}`` (if ``damp == 0``) or more generally ``(A'A +
  146. damp^2*I)^{-1}``. This is well defined if A has full column
  147. rank or ``damp > 0``. (Not sure what var means if ``rank(A)
  148. < n`` and ``damp = 0.``)
  149. Notes
  150. -----
  151. LSQR uses an iterative method to approximate the solution. The
  152. number of iterations required to reach a certain accuracy depends
  153. strongly on the scaling of the problem. Poor scaling of the rows
  154. or columns of A should therefore be avoided where possible.
  155. For example, in problem 1 the solution is unaltered by
  156. row-scaling. If a row of A is very small or large compared to
  157. the other rows of A, the corresponding row of ( A b ) should be
  158. scaled up or down.
  159. In problems 1 and 2, the solution x is easily recovered
  160. following column-scaling. Unless better information is known,
  161. the nonzero columns of A should be scaled so that they all have
  162. the same Euclidean norm (e.g., 1.0).
  163. In problem 3, there is no freedom to re-scale if damp is
  164. nonzero. However, the value of damp should be assigned only
  165. after attention has been paid to the scaling of A.
  166. The parameter damp is intended to help regularize
  167. ill-conditioned systems, by preventing the true solution from
  168. being very large. Another aid to regularization is provided by
  169. the parameter acond, which may be used to terminate iterations
  170. before the computed solution becomes very large.
  171. If some initial estimate ``x0`` is known and if ``damp == 0``,
  172. one could proceed as follows:
  173. 1. Compute a residual vector ``r0 = b - A*x0``.
  174. 2. Use LSQR to solve the system ``A*dx = r0``.
  175. 3. Add the correction dx to obtain a final solution ``x = x0 + dx``.
  176. This requires that ``x0`` be available before and after the call
  177. to LSQR. To judge the benefits, suppose LSQR takes k1 iterations
  178. to solve A*x = b and k2 iterations to solve A*dx = r0.
  179. If x0 is "good", norm(r0) will be smaller than norm(b).
  180. If the same stopping tolerances atol and btol are used for each
  181. system, k1 and k2 will be similar, but the final solution x0 + dx
  182. should be more accurate. The only way to reduce the total work
  183. is to use a larger stopping tolerance for the second system.
  184. If some value btol is suitable for A*x = b, the larger value
  185. btol*norm(b)/norm(r0) should be suitable for A*dx = r0.
  186. Preconditioning is another way to reduce the number of iterations.
  187. If it is possible to solve a related system ``M*x = b``
  188. efficiently, where M approximates A in some helpful way (e.g. M -
  189. A has low rank or its elements are small relative to those of A),
  190. LSQR may converge more rapidly on the system ``A*M(inverse)*z =
  191. b``, after which x can be recovered by solving M*x = z.
  192. If A is symmetric, LSQR should not be used!
  193. Alternatives are the symmetric conjugate-gradient method (cg)
  194. and/or SYMMLQ. SYMMLQ is an implementation of symmetric cg that
  195. applies to any symmetric A and will converge more rapidly than
  196. LSQR. If A is positive definite, there are other implementations
  197. of symmetric cg that require slightly less work per iteration than
  198. SYMMLQ (but will take the same number of iterations).
  199. References
  200. ----------
  201. .. [1] C. C. Paige and M. A. Saunders (1982a).
  202. "LSQR: An algorithm for sparse linear equations and
  203. sparse least squares", ACM TOMS 8(1), 43-71.
  204. .. [2] C. C. Paige and M. A. Saunders (1982b).
  205. "Algorithm 583. LSQR: Sparse linear equations and least
  206. squares problems", ACM TOMS 8(2), 195-209.
  207. .. [3] M. A. Saunders (1995). "Solution of sparse rectangular
  208. systems using LSQR and CRAIG", BIT 35, 588-604.
  209. Examples
  210. --------
  211. >>> from scipy.sparse import csc_matrix
  212. >>> from scipy.sparse.linalg import lsqr
  213. >>> A = csc_matrix([[1., 0.], [1., 1.], [0., 1.]], dtype=float)
  214. The first example has the trivial solution `[0, 0]`
  215. >>> b = np.array([0., 0., 0.], dtype=float)
  216. >>> x, istop, itn, normr = lsqr(A, b)[:4]
  217. The exact solution is x = 0
  218. >>> istop
  219. 0
  220. >>> x
  221. array([ 0., 0.])
  222. The stopping code `istop=0` returned indicates that a vector of zeros was
  223. found as a solution. The returned solution `x` indeed contains `[0., 0.]`.
  224. The next example has a non-trivial solution:
  225. >>> b = np.array([1., 0., -1.], dtype=float)
  226. >>> x, istop, itn, r1norm = lsqr(A, b)[:4]
  227. >>> istop
  228. 1
  229. >>> x
  230. array([ 1., -1.])
  231. >>> itn
  232. 1
  233. >>> r1norm
  234. 4.440892098500627e-16
  235. As indicated by `istop=1`, `lsqr` found a solution obeying the tolerance
  236. limits. The given solution `[1., -1.]` obviously solves the equation. The
  237. remaining return values include information about the number of iterations
  238. (`itn=1`) and the remaining difference of left and right side of the solved
  239. equation.
  240. The final example demonstrates the behavior in the case where there is no
  241. solution for the equation:
  242. >>> b = np.array([1., 0.01, -1.], dtype=float)
  243. >>> x, istop, itn, r1norm = lsqr(A, b)[:4]
  244. >>> istop
  245. 2
  246. >>> x
  247. array([ 1.00333333, -0.99666667])
  248. >>> A.dot(x)-b
  249. array([ 0.00333333, -0.00333333, 0.00333333])
  250. >>> r1norm
  251. 0.005773502691896255
  252. `istop` indicates that the system is inconsistent and thus `x` is rather an
  253. approximate solution to the corresponding least-squares problem. `r1norm`
  254. contains the norm of the minimal residual that was found.
  255. """
  256. A = aslinearoperator(A)
  257. b = np.atleast_1d(b)
  258. if b.ndim > 1:
  259. b = b.squeeze()
  260. m, n = A.shape
  261. if iter_lim is None:
  262. iter_lim = 2 * n
  263. var = np.zeros(n)
  264. msg = ('The exact solution is x = 0 ',
  265. 'Ax - b is small enough, given atol, btol ',
  266. 'The least-squares solution is good enough, given atol ',
  267. 'The estimate of cond(Abar) has exceeded conlim ',
  268. 'Ax - b is small enough for this machine ',
  269. 'The least-squares solution is good enough for this machine',
  270. 'Cond(Abar) seems to be too large for this machine ',
  271. 'The iteration limit has been reached ')
  272. if show:
  273. print(' ')
  274. print('LSQR Least-squares solution of Ax = b')
  275. str1 = 'The matrix A has %8g rows and %8g cols' % (m, n)
  276. str2 = 'damp = %20.14e calc_var = %8g' % (damp, calc_var)
  277. str3 = 'atol = %8.2e conlim = %8.2e' % (atol, conlim)
  278. str4 = 'btol = %8.2e iter_lim = %8g' % (btol, iter_lim)
  279. print(str1)
  280. print(str2)
  281. print(str3)
  282. print(str4)
  283. itn = 0
  284. istop = 0
  285. ctol = 0
  286. if conlim > 0:
  287. ctol = 1/conlim
  288. anorm = 0
  289. acond = 0
  290. dampsq = damp**2
  291. ddnorm = 0
  292. res2 = 0
  293. xnorm = 0
  294. xxnorm = 0
  295. z = 0
  296. cs2 = -1
  297. sn2 = 0
  298. """
  299. Set up the first vectors u and v for the bidiagonalization.
  300. These satisfy beta*u = b - A*x, alfa*v = A'*u.
  301. """
  302. u = b
  303. bnorm = np.linalg.norm(b)
  304. if x0 is None:
  305. x = np.zeros(n)
  306. beta = bnorm.copy()
  307. else:
  308. x = np.asarray(x0)
  309. u = u - A.matvec(x)
  310. beta = np.linalg.norm(u)
  311. if beta > 0:
  312. u = (1/beta) * u
  313. v = A.rmatvec(u)
  314. alfa = np.linalg.norm(v)
  315. else:
  316. v = x.copy()
  317. alfa = 0
  318. if alfa > 0:
  319. v = (1/alfa) * v
  320. w = v.copy()
  321. rhobar = alfa
  322. phibar = beta
  323. rnorm = beta
  324. r1norm = rnorm
  325. r2norm = rnorm
  326. # Reverse the order here from the original matlab code because
  327. # there was an error on return when arnorm==0
  328. arnorm = alfa * beta
  329. if arnorm == 0:
  330. print(msg[0])
  331. return x, istop, itn, r1norm, r2norm, anorm, acond, arnorm, xnorm, var
  332. head1 = ' Itn x[0] r1norm r2norm '
  333. head2 = ' Compatible LS Norm A Cond A'
  334. if show:
  335. print(' ')
  336. print(head1, head2)
  337. test1 = 1
  338. test2 = alfa / beta
  339. str1 = '%6g %12.5e' % (itn, x[0])
  340. str2 = ' %10.3e %10.3e' % (r1norm, r2norm)
  341. str3 = ' %8.1e %8.1e' % (test1, test2)
  342. print(str1, str2, str3)
  343. # Main iteration loop.
  344. while itn < iter_lim:
  345. itn = itn + 1
  346. """
  347. % Perform the next step of the bidiagonalization to obtain the
  348. % next beta, u, alfa, v. These satisfy the relations
  349. % beta*u = a*v - alfa*u,
  350. % alfa*v = A'*u - beta*v.
  351. """
  352. u = A.matvec(v) - alfa * u
  353. beta = np.linalg.norm(u)
  354. if beta > 0:
  355. u = (1/beta) * u
  356. anorm = sqrt(anorm**2 + alfa**2 + beta**2 + damp**2)
  357. v = A.rmatvec(u) - beta * v
  358. alfa = np.linalg.norm(v)
  359. if alfa > 0:
  360. v = (1 / alfa) * v
  361. # Use a plane rotation to eliminate the damping parameter.
  362. # This alters the diagonal (rhobar) of the lower-bidiagonal matrix.
  363. rhobar1 = sqrt(rhobar**2 + damp**2)
  364. cs1 = rhobar / rhobar1
  365. sn1 = damp / rhobar1
  366. psi = sn1 * phibar
  367. phibar = cs1 * phibar
  368. # Use a plane rotation to eliminate the subdiagonal element (beta)
  369. # of the lower-bidiagonal matrix, giving an upper-bidiagonal matrix.
  370. cs, sn, rho = _sym_ortho(rhobar1, beta)
  371. theta = sn * alfa
  372. rhobar = -cs * alfa
  373. phi = cs * phibar
  374. phibar = sn * phibar
  375. tau = sn * phi
  376. # Update x and w.
  377. t1 = phi / rho
  378. t2 = -theta / rho
  379. dk = (1 / rho) * w
  380. x = x + t1 * w
  381. w = v + t2 * w
  382. ddnorm = ddnorm + np.linalg.norm(dk)**2
  383. if calc_var:
  384. var = var + dk**2
  385. # Use a plane rotation on the right to eliminate the
  386. # super-diagonal element (theta) of the upper-bidiagonal matrix.
  387. # Then use the result to estimate norm(x).
  388. delta = sn2 * rho
  389. gambar = -cs2 * rho
  390. rhs = phi - delta * z
  391. zbar = rhs / gambar
  392. xnorm = sqrt(xxnorm + zbar**2)
  393. gamma = sqrt(gambar**2 + theta**2)
  394. cs2 = gambar / gamma
  395. sn2 = theta / gamma
  396. z = rhs / gamma
  397. xxnorm = xxnorm + z**2
  398. # Test for convergence.
  399. # First, estimate the condition of the matrix Abar,
  400. # and the norms of rbar and Abar'rbar.
  401. acond = anorm * sqrt(ddnorm)
  402. res1 = phibar**2
  403. res2 = res2 + psi**2
  404. rnorm = sqrt(res1 + res2)
  405. arnorm = alfa * abs(tau)
  406. # Distinguish between
  407. # r1norm = ||b - Ax|| and
  408. # r2norm = rnorm in current code
  409. # = sqrt(r1norm^2 + damp^2*||x||^2).
  410. # Estimate r1norm from
  411. # r1norm = sqrt(r2norm^2 - damp^2*||x||^2).
  412. # Although there is cancellation, it might be accurate enough.
  413. r1sq = rnorm**2 - dampsq * xxnorm
  414. r1norm = sqrt(abs(r1sq))
  415. if r1sq < 0:
  416. r1norm = -r1norm
  417. r2norm = rnorm
  418. # Now use these norms to estimate certain other quantities,
  419. # some of which will be small near a solution.
  420. test1 = rnorm / bnorm
  421. test2 = arnorm / (anorm * rnorm + eps)
  422. test3 = 1 / (acond + eps)
  423. t1 = test1 / (1 + anorm * xnorm / bnorm)
  424. rtol = btol + atol * anorm * xnorm / bnorm
  425. # The following tests guard against extremely small values of
  426. # atol, btol or ctol. (The user may have set any or all of
  427. # the parameters atol, btol, conlim to 0.)
  428. # The effect is equivalent to the normal tests using
  429. # atol = eps, btol = eps, conlim = 1/eps.
  430. if itn >= iter_lim:
  431. istop = 7
  432. if 1 + test3 <= 1:
  433. istop = 6
  434. if 1 + test2 <= 1:
  435. istop = 5
  436. if 1 + t1 <= 1:
  437. istop = 4
  438. # Allow for tolerances set by the user.
  439. if test3 <= ctol:
  440. istop = 3
  441. if test2 <= atol:
  442. istop = 2
  443. if test1 <= rtol:
  444. istop = 1
  445. # See if it is time to print something.
  446. prnt = False
  447. if n <= 40:
  448. prnt = True
  449. if itn <= 10:
  450. prnt = True
  451. if itn >= iter_lim-10:
  452. prnt = True
  453. # if itn%10 == 0: prnt = True
  454. if test3 <= 2*ctol:
  455. prnt = True
  456. if test2 <= 10*atol:
  457. prnt = True
  458. if test1 <= 10*rtol:
  459. prnt = True
  460. if istop != 0:
  461. prnt = True
  462. if prnt:
  463. if show:
  464. str1 = '%6g %12.5e' % (itn, x[0])
  465. str2 = ' %10.3e %10.3e' % (r1norm, r2norm)
  466. str3 = ' %8.1e %8.1e' % (test1, test2)
  467. str4 = ' %8.1e %8.1e' % (anorm, acond)
  468. print(str1, str2, str3, str4)
  469. if istop != 0:
  470. break
  471. # End of iteration loop.
  472. # Print the stopping condition.
  473. if show:
  474. print(' ')
  475. print('LSQR finished')
  476. print(msg[istop])
  477. print(' ')
  478. str1 = 'istop =%8g r1norm =%8.1e' % (istop, r1norm)
  479. str2 = 'anorm =%8.1e arnorm =%8.1e' % (anorm, arnorm)
  480. str3 = 'itn =%8g r2norm =%8.1e' % (itn, r2norm)
  481. str4 = 'acond =%8.1e xnorm =%8.1e' % (acond, xnorm)
  482. print(str1 + ' ' + str2)
  483. print(str3 + ' ' + str4)
  484. print(' ')
  485. return x, istop, itn, r1norm, r2norm, anorm, acond, arnorm, xnorm, var