lsqr.py 21 KB

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