lsqr.py 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589
  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. __all__ = ['lsqr']
  39. import numpy as np
  40. from math import sqrt
  41. from scipy.sparse.linalg._interface import aslinearoperator
  42. from scipy.sparse._sputils import convert_pydata_sparse_to_scipy
  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-6, btol=1e-6, 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 ||Ax - b||^2`` or
  79. ``min ||Ax - b||^2 + d^2 ||x - x0||^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 Ax = b
  84. 2. Linear least squares -- solve Ax = b
  85. in the least-squares sense
  86. 3. Damped least squares -- solve ( A )*x = ( b )
  87. ( damp*I ) ( damp*x0 )
  88. in the least-squares sense
  89. Parameters
  90. ----------
  91. A : {sparse array, ndarray, LinearOperator}
  92. Representation of an m-by-n matrix.
  93. Alternatively, ``A`` can be a linear operator which can
  94. produce ``Ax`` and ``A^T x`` using, e.g.,
  95. ``scipy.sparse.linalg.LinearOperator``.
  96. b : array_like, shape (m,)
  97. Right-hand side vector ``b``.
  98. damp : float
  99. Damping coefficient. Default is 0.
  100. atol, btol : float, optional
  101. Stopping tolerances. `lsqr` continues iterations until a
  102. certain backward error estimate is smaller than some quantity
  103. depending on atol and btol. Let ``r = b - Ax`` be the
  104. residual vector for the current approximate solution ``x``.
  105. If ``Ax = b`` seems to be consistent, `lsqr` terminates
  106. when ``norm(r) <= atol * norm(A) * norm(x) + btol * norm(b)``.
  107. Otherwise, `lsqr` terminates when ``norm(A^H r) <=
  108. atol * norm(A) * norm(r)``. If both tolerances are 1.0e-6 (default),
  109. the final ``norm(r)`` should be accurate to about 6
  110. digits. (The final ``x`` will usually have fewer correct digits,
  111. depending on ``cond(A)`` and the size of LAMBDA.) If `atol`
  112. or `btol` is None, a default value of 1.0e-6 will be used.
  113. Ideally, they should be estimates of the relative error in the
  114. entries of ``A`` and ``b`` respectively. For example, if the entries
  115. of ``A`` have 7 correct digits, set ``atol = 1e-7``. This prevents
  116. the algorithm from doing unnecessary work beyond the
  117. uncertainty of the input data.
  118. conlim : float, optional
  119. Another stopping tolerance. lsqr terminates if an estimate of
  120. ``cond(A)`` exceeds `conlim`. For compatible systems ``Ax =
  121. b``, `conlim` could be as large as 1.0e+12 (say). For
  122. least-squares problems, conlim should be less than 1.0e+8.
  123. Maximum precision can be obtained by setting ``atol = btol =
  124. conlim = zero``, but the number of iterations may then be
  125. excessive. Default is 1e8.
  126. iter_lim : int, optional
  127. Explicit limitation on number of iterations (for safety).
  128. show : bool, optional
  129. Display an iteration log. Default is False.
  130. calc_var : bool, optional
  131. Whether to estimate diagonals of ``(A'A + damp^2*I)^{-1}``.
  132. x0 : array_like, shape (n,), optional
  133. Initial guess of x, if None zeros are used. Default is None.
  134. .. versionadded:: 1.0.0
  135. Returns
  136. -------
  137. x : ndarray of float
  138. The final solution.
  139. istop : int
  140. Gives the reason for termination.
  141. 1 means x is an approximate solution to Ax = b.
  142. 2 means x approximately solves the least-squares problem.
  143. itn : int
  144. Iteration number upon termination.
  145. r1norm : float
  146. ``norm(r)``, where ``r = b - Ax``.
  147. r2norm : float
  148. ``sqrt( norm(r)^2 + damp^2 * norm(x - x0)^2 )``. Equal to `r1norm`
  149. if ``damp == 0``.
  150. anorm : float
  151. Estimate of Frobenius norm of ``Abar = [[A]; [damp*I]]``.
  152. acond : float
  153. Estimate of ``cond(Abar)``.
  154. arnorm : float
  155. Estimate of ``norm(A'@r - damp^2*(x - x0))``.
  156. xnorm : float
  157. ``norm(x)``
  158. var : ndarray of float
  159. If ``calc_var`` is True, estimates all diagonals of
  160. ``(A'A)^{-1}`` (if ``damp == 0``) or more generally ``(A'A +
  161. damp^2*I)^{-1}``. This is well defined if A has full column
  162. rank or ``damp > 0``. (Not sure what var means if ``rank(A)
  163. < n`` and ``damp = 0.``)
  164. Notes
  165. -----
  166. LSQR uses an iterative method to approximate the solution. The
  167. number of iterations required to reach a certain accuracy depends
  168. strongly on the scaling of the problem. Poor scaling of the rows
  169. or columns of A should therefore be avoided where possible.
  170. For example, in problem 1 the solution is unaltered by
  171. row-scaling. If a row of A is very small or large compared to
  172. the other rows of A, the corresponding row of ( A b ) should be
  173. scaled up or down.
  174. In problems 1 and 2, the solution x is easily recovered
  175. following column-scaling. Unless better information is known,
  176. the nonzero columns of A should be scaled so that they all have
  177. the same Euclidean norm (e.g., 1.0).
  178. In problem 3, there is no freedom to re-scale if damp is
  179. nonzero. However, the value of damp should be assigned only
  180. after attention has been paid to the scaling of A.
  181. The parameter damp is intended to help regularize
  182. ill-conditioned systems, by preventing the true solution from
  183. being very large. Another aid to regularization is provided by
  184. the parameter acond, which may be used to terminate iterations
  185. before the computed solution becomes very large.
  186. If some initial estimate ``x0`` is known and if ``damp == 0``,
  187. one could proceed as follows:
  188. 1. Compute a residual vector ``r0 = b - A@x0``.
  189. 2. Use LSQR to solve the system ``A@dx = r0``.
  190. 3. Add the correction dx to obtain a final solution ``x = x0 + dx``.
  191. This requires that ``x0`` be available before and after the call
  192. to LSQR. To judge the benefits, suppose LSQR takes k1 iterations
  193. to solve A@x = b and k2 iterations to solve A@dx = r0.
  194. If x0 is "good", norm(r0) will be smaller than norm(b).
  195. If the same stopping tolerances atol and btol are used for each
  196. system, k1 and k2 will be similar, but the final solution x0 + dx
  197. should be more accurate. The only way to reduce the total work
  198. is to use a larger stopping tolerance for the second system.
  199. If some value btol is suitable for A@x = b, the larger value
  200. btol*norm(b)/norm(r0) should be suitable for A@dx = r0.
  201. Preconditioning is another way to reduce the number of iterations.
  202. If it is possible to solve a related system ``M@x = b``
  203. efficiently, where M approximates A in some helpful way (e.g. M -
  204. A has low rank or its elements are small relative to those of A),
  205. LSQR may converge more rapidly on the system ``A@M(inverse)@z =
  206. b``, after which x can be recovered by solving M@x = z.
  207. If A is symmetric, LSQR should not be used!
  208. Alternatives are the symmetric conjugate-gradient method (cg)
  209. and/or SYMMLQ. SYMMLQ is an implementation of symmetric cg that
  210. applies to any symmetric A and will converge more rapidly than
  211. LSQR. If A is positive definite, there are other implementations
  212. of symmetric cg that require slightly less work per iteration than
  213. SYMMLQ (but will take the same number of iterations).
  214. References
  215. ----------
  216. .. [1] C. C. Paige and M. A. Saunders (1982a).
  217. "LSQR: An algorithm for sparse linear equations and
  218. sparse least squares", ACM TOMS 8(1), 43-71.
  219. .. [2] C. C. Paige and M. A. Saunders (1982b).
  220. "Algorithm 583. LSQR: Sparse linear equations and least
  221. squares problems", ACM TOMS 8(2), 195-209.
  222. .. [3] M. A. Saunders (1995). "Solution of sparse rectangular
  223. systems using LSQR and CRAIG", BIT 35, 588-604.
  224. Examples
  225. --------
  226. >>> import numpy as np
  227. >>> from scipy.sparse import csc_array
  228. >>> from scipy.sparse.linalg import lsqr
  229. >>> A = csc_array([[1., 0.], [1., 1.], [0., 1.]], dtype=float)
  230. The first example has the trivial solution ``[0, 0]``
  231. >>> b = np.array([0., 0., 0.], dtype=float)
  232. >>> x, istop, itn, normr = lsqr(A, b)[:4]
  233. >>> istop
  234. 0
  235. >>> x
  236. array([ 0., 0.])
  237. The stopping code ``istop=0`` returned indicates that a vector of zeros was
  238. found as a solution. The returned solution `x` indeed contains
  239. ``[0., 0.]``. The next example has a non-trivial solution:
  240. >>> b = np.array([1., 0., -1.], dtype=float)
  241. >>> x, istop, itn, r1norm = lsqr(A, b)[:4]
  242. >>> istop
  243. 1
  244. >>> x
  245. array([ 1., -1.])
  246. >>> itn
  247. 1
  248. >>> r1norm
  249. 4.440892098500627e-16
  250. As indicated by ``istop=1``, `lsqr` found a solution obeying the tolerance
  251. limits. The given solution ``[1., -1.]`` obviously solves the equation. The
  252. remaining return values include information about the number of iterations
  253. (`itn=1`) and the remaining difference of left and right side of the solved
  254. equation.
  255. The final example demonstrates the behavior in the case where there is no
  256. solution for the equation:
  257. >>> b = np.array([1., 0.01, -1.], dtype=float)
  258. >>> x, istop, itn, r1norm = lsqr(A, b)[:4]
  259. >>> istop
  260. 2
  261. >>> x
  262. array([ 1.00333333, -0.99666667])
  263. >>> A.dot(x)-b
  264. array([ 0.00333333, -0.00333333, 0.00333333])
  265. >>> r1norm
  266. 0.005773502691896255
  267. `istop` indicates that the system is inconsistent and thus `x` is rather an
  268. approximate solution to the corresponding least-squares problem. `r1norm`
  269. contains the norm of the minimal residual that was found.
  270. """
  271. A = convert_pydata_sparse_to_scipy(A)
  272. A = aslinearoperator(A)
  273. b = np.atleast_1d(b)
  274. if b.ndim > 1:
  275. b = b.squeeze()
  276. m, n = A.shape
  277. if iter_lim is None:
  278. iter_lim = 2 * n
  279. var = np.zeros(n)
  280. msg = ('The exact solution is x = 0 ',
  281. 'Ax - b is small enough, given atol, btol ',
  282. 'The least-squares solution is good enough, given atol ',
  283. 'The estimate of cond(Abar) has exceeded conlim ',
  284. 'Ax - b is small enough for this machine ',
  285. 'The least-squares solution is good enough for this machine',
  286. 'Cond(Abar) seems to be too large for this machine ',
  287. 'The iteration limit has been reached ')
  288. if show:
  289. print(' ')
  290. print('LSQR Least-squares solution of Ax = b')
  291. str1 = f'The matrix A has {m} rows and {n} columns'
  292. str2 = f'damp = {damp:20.14e} calc_var = {calc_var:8g}'
  293. str3 = f'atol = {atol:8.2e} conlim = {conlim:8.2e}'
  294. str4 = f'btol = {btol:8.2e} iter_lim = {iter_lim:8g}'
  295. print(str1)
  296. print(str2)
  297. print(str3)
  298. print(str4)
  299. itn = 0
  300. istop = 0
  301. ctol = 0
  302. if conlim > 0:
  303. ctol = 1/conlim
  304. anorm = 0
  305. acond = 0
  306. dampsq = damp**2
  307. ddnorm = 0
  308. res2 = 0
  309. xnorm = 0
  310. xxnorm = 0
  311. z = 0
  312. cs2 = -1
  313. sn2 = 0
  314. # Set up the first vectors u and v for the bidiagonalization.
  315. # These satisfy beta*u = b - A@x, alfa*v = A'@u.
  316. u = b
  317. bnorm = np.linalg.norm(b)
  318. if x0 is None:
  319. x = np.zeros(n)
  320. beta = bnorm.copy()
  321. else:
  322. x = np.asarray(x0)
  323. u = u - A.matvec(x)
  324. beta = np.linalg.norm(u)
  325. if beta > 0:
  326. u = (1/beta) * u
  327. v = A.rmatvec(u)
  328. alfa = np.linalg.norm(v)
  329. else:
  330. v = x.copy()
  331. alfa = 0
  332. if alfa > 0:
  333. v = (1/alfa) * v
  334. w = v.copy()
  335. rhobar = alfa
  336. phibar = beta
  337. rnorm = beta
  338. r1norm = rnorm
  339. r2norm = rnorm
  340. # Reverse the order here from the original matlab code because
  341. # there was an error on return when arnorm==0
  342. arnorm = alfa * beta
  343. if arnorm == 0:
  344. if show:
  345. print(msg[0])
  346. return x, istop, itn, r1norm, r2norm, anorm, acond, arnorm, xnorm, var
  347. head1 = ' Itn x[0] r1norm r2norm '
  348. head2 = ' Compatible LS Norm A Cond A'
  349. if show:
  350. print(' ')
  351. print(head1, head2)
  352. test1 = 1
  353. test2 = alfa / beta
  354. str1 = f'{itn:6g} {x[0]:12.5e}'
  355. str2 = f' {r1norm:10.3e} {r2norm:10.3e}'
  356. str3 = f' {test1:8.1e} {test2:8.1e}'
  357. print(str1, str2, str3)
  358. # Main iteration loop.
  359. while itn < iter_lim:
  360. itn = itn + 1
  361. # Perform the next step of the bidiagonalization to obtain the
  362. # next beta, u, alfa, v. These satisfy the relations
  363. # beta*u = a@v - alfa*u,
  364. # alfa*v = A'@u - beta*v.
  365. u = A.matvec(v) - alfa * u
  366. beta = np.linalg.norm(u)
  367. if beta > 0:
  368. u = (1/beta) * u
  369. anorm = sqrt(anorm**2 + alfa**2 + beta**2 + dampsq)
  370. v = A.rmatvec(u) - beta * v
  371. alfa = np.linalg.norm(v)
  372. if alfa > 0:
  373. v = (1 / alfa) * v
  374. # Use a plane rotation to eliminate the damping parameter.
  375. # This alters the diagonal (rhobar) of the lower-bidiagonal matrix.
  376. if damp > 0:
  377. rhobar1 = sqrt(rhobar**2 + dampsq)
  378. cs1 = rhobar / rhobar1
  379. sn1 = damp / rhobar1
  380. psi = sn1 * phibar
  381. phibar = cs1 * phibar
  382. else:
  383. # cs1 = 1 and sn1 = 0
  384. rhobar1 = rhobar
  385. psi = 0.
  386. # Use a plane rotation to eliminate the subdiagonal element (beta)
  387. # of the lower-bidiagonal matrix, giving an upper-bidiagonal matrix.
  388. cs, sn, rho = _sym_ortho(rhobar1, beta)
  389. theta = sn * alfa
  390. rhobar = -cs * alfa
  391. phi = cs * phibar
  392. phibar = sn * phibar
  393. tau = sn * phi
  394. # Update x and w.
  395. t1 = phi / rho
  396. t2 = -theta / rho
  397. dk = (1 / rho) * w
  398. x = x + t1 * w
  399. w = v + t2 * w
  400. ddnorm = ddnorm + np.linalg.norm(dk)**2
  401. if calc_var:
  402. var = var + dk**2
  403. # Use a plane rotation on the right to eliminate the
  404. # super-diagonal element (theta) of the upper-bidiagonal matrix.
  405. # Then use the result to estimate norm(x).
  406. delta = sn2 * rho
  407. gambar = -cs2 * rho
  408. rhs = phi - delta * z
  409. zbar = rhs / gambar
  410. xnorm = sqrt(xxnorm + zbar**2)
  411. gamma = sqrt(gambar**2 + theta**2)
  412. cs2 = gambar / gamma
  413. sn2 = theta / gamma
  414. z = rhs / gamma
  415. xxnorm = xxnorm + z**2
  416. # Test for convergence.
  417. # First, estimate the condition of the matrix Abar,
  418. # and the norms of rbar and Abar'rbar.
  419. acond = anorm * sqrt(ddnorm)
  420. res1 = phibar**2
  421. res2 = res2 + psi**2
  422. rnorm = sqrt(res1 + res2)
  423. arnorm = alfa * abs(tau)
  424. # Distinguish between
  425. # r1norm = ||b - Ax|| and
  426. # r2norm = rnorm in current code
  427. # = sqrt(r1norm^2 + damp^2*||x - x0||^2).
  428. # Estimate r1norm from
  429. # r1norm = sqrt(r2norm^2 - damp^2*||x - x0||^2).
  430. # Although there is cancellation, it might be accurate enough.
  431. if damp > 0:
  432. r1sq = rnorm**2 - dampsq * xxnorm
  433. r1norm = sqrt(abs(r1sq))
  434. if r1sq < 0:
  435. r1norm = -r1norm
  436. else:
  437. r1norm = rnorm
  438. r2norm = rnorm
  439. # Now use these norms to estimate certain other quantities,
  440. # some of which will be small near a solution.
  441. test1 = rnorm / bnorm
  442. test2 = arnorm / (anorm * rnorm + eps)
  443. test3 = 1 / (acond + eps)
  444. t1 = test1 / (1 + anorm * xnorm / bnorm)
  445. rtol = btol + atol * anorm * xnorm / bnorm
  446. # The following tests guard against extremely small values of
  447. # atol, btol or ctol. (The user may have set any or all of
  448. # the parameters atol, btol, conlim to 0.)
  449. # The effect is equivalent to the normal tests using
  450. # atol = eps, btol = eps, conlim = 1/eps.
  451. if itn >= iter_lim:
  452. istop = 7
  453. if 1 + test3 <= 1:
  454. istop = 6
  455. if 1 + test2 <= 1:
  456. istop = 5
  457. if 1 + t1 <= 1:
  458. istop = 4
  459. # Allow for tolerances set by the user.
  460. if test3 <= ctol:
  461. istop = 3
  462. if test2 <= atol:
  463. istop = 2
  464. if test1 <= rtol:
  465. istop = 1
  466. if show:
  467. # See if it is time to print something.
  468. prnt = False
  469. if n <= 40:
  470. prnt = True
  471. if itn <= 10:
  472. prnt = True
  473. if itn >= iter_lim-10:
  474. prnt = True
  475. # if itn%10 == 0: prnt = True
  476. if test3 <= 2*ctol:
  477. prnt = True
  478. if test2 <= 10*atol:
  479. prnt = True
  480. if test1 <= 10*rtol:
  481. prnt = True
  482. if istop != 0:
  483. prnt = True
  484. if prnt:
  485. str1 = f'{itn:6g} {x[0]:12.5e}'
  486. str2 = f' {r1norm:10.3e} {r2norm:10.3e}'
  487. str3 = f' {test1:8.1e} {test2:8.1e}'
  488. str4 = f' {anorm:8.1e} {acond:8.1e}'
  489. print(str1, str2, str3, str4)
  490. if istop != 0:
  491. break
  492. # End of iteration loop.
  493. # Print the stopping condition.
  494. if show:
  495. print(' ')
  496. print('LSQR finished')
  497. print(msg[istop])
  498. print(' ')
  499. str1 = f'istop ={istop:8g} r1norm ={r1norm:8.1e}'
  500. str2 = f'anorm ={anorm:8.1e} arnorm ={arnorm:8.1e}'
  501. str3 = f'itn ={itn:8g} r2norm ={r2norm:8.1e}'
  502. str4 = f'acond ={acond:8.1e} xnorm ={xnorm:8.1e}'
  503. print(str1 + ' ' + str2)
  504. print(str3 + ' ' + str4)
  505. print(' ')
  506. return x, istop, itn, r1norm, r2norm, anorm, acond, arnorm, xnorm, var