common.py 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735
  1. """Functions used by least-squares algorithms."""
  2. from __future__ import division, print_function, absolute_import
  3. from math import copysign
  4. import numpy as np
  5. from numpy.linalg import norm
  6. from scipy.linalg import cho_factor, cho_solve, LinAlgError
  7. from scipy.sparse import issparse
  8. from scipy.sparse.linalg import LinearOperator, aslinearoperator
  9. EPS = np.finfo(float).eps
  10. # Functions related to a trust-region problem.
  11. def intersect_trust_region(x, s, Delta):
  12. """Find the intersection of a line with the boundary of a trust region.
  13. This function solves the quadratic equation with respect to t
  14. ||(x + s*t)||**2 = Delta**2.
  15. Returns
  16. -------
  17. t_neg, t_pos : tuple of float
  18. Negative and positive roots.
  19. Raises
  20. ------
  21. ValueError
  22. If `s` is zero or `x` is not within the trust region.
  23. """
  24. a = np.dot(s, s)
  25. if a == 0:
  26. raise ValueError("`s` is zero.")
  27. b = np.dot(x, s)
  28. c = np.dot(x, x) - Delta**2
  29. if c > 0:
  30. raise ValueError("`x` is not within the trust region.")
  31. d = np.sqrt(b*b - a*c) # Root from one fourth of the discriminant.
  32. # Computations below avoid loss of significance, see "Numerical Recipes".
  33. q = -(b + copysign(d, b))
  34. t1 = q / a
  35. t2 = c / q
  36. if t1 < t2:
  37. return t1, t2
  38. else:
  39. return t2, t1
  40. def solve_lsq_trust_region(n, m, uf, s, V, Delta, initial_alpha=None,
  41. rtol=0.01, max_iter=10):
  42. """Solve a trust-region problem arising in least-squares minimization.
  43. This function implements a method described by J. J. More [1]_ and used
  44. in MINPACK, but it relies on a single SVD of Jacobian instead of series
  45. of Cholesky decompositions. Before running this function, compute:
  46. ``U, s, VT = svd(J, full_matrices=False)``.
  47. Parameters
  48. ----------
  49. n : int
  50. Number of variables.
  51. m : int
  52. Number of residuals.
  53. uf : ndarray
  54. Computed as U.T.dot(f).
  55. s : ndarray
  56. Singular values of J.
  57. V : ndarray
  58. Transpose of VT.
  59. Delta : float
  60. Radius of a trust region.
  61. initial_alpha : float, optional
  62. Initial guess for alpha, which might be available from a previous
  63. iteration. If None, determined automatically.
  64. rtol : float, optional
  65. Stopping tolerance for the root-finding procedure. Namely, the
  66. solution ``p`` will satisfy ``abs(norm(p) - Delta) < rtol * Delta``.
  67. max_iter : int, optional
  68. Maximum allowed number of iterations for the root-finding procedure.
  69. Returns
  70. -------
  71. p : ndarray, shape (n,)
  72. Found solution of a trust-region problem.
  73. alpha : float
  74. Positive value such that (J.T*J + alpha*I)*p = -J.T*f.
  75. Sometimes called Levenberg-Marquardt parameter.
  76. n_iter : int
  77. Number of iterations made by root-finding procedure. Zero means
  78. that Gauss-Newton step was selected as the solution.
  79. References
  80. ----------
  81. .. [1] More, J. J., "The Levenberg-Marquardt Algorithm: Implementation
  82. and Theory," Numerical Analysis, ed. G. A. Watson, Lecture Notes
  83. in Mathematics 630, Springer Verlag, pp. 105-116, 1977.
  84. """
  85. def phi_and_derivative(alpha, suf, s, Delta):
  86. """Function of which to find zero.
  87. It is defined as "norm of regularized (by alpha) least-squares
  88. solution minus `Delta`". Refer to [1]_.
  89. """
  90. denom = s**2 + alpha
  91. p_norm = norm(suf / denom)
  92. phi = p_norm - Delta
  93. phi_prime = -np.sum(suf ** 2 / denom**3) / p_norm
  94. return phi, phi_prime
  95. suf = s * uf
  96. # Check if J has full rank and try Gauss-Newton step.
  97. if m >= n:
  98. threshold = EPS * m * s[0]
  99. full_rank = s[-1] > threshold
  100. else:
  101. full_rank = False
  102. if full_rank:
  103. p = -V.dot(uf / s)
  104. if norm(p) <= Delta:
  105. return p, 0.0, 0
  106. alpha_upper = norm(suf) / Delta
  107. if full_rank:
  108. phi, phi_prime = phi_and_derivative(0.0, suf, s, Delta)
  109. alpha_lower = -phi / phi_prime
  110. else:
  111. alpha_lower = 0.0
  112. if initial_alpha is None or not full_rank and initial_alpha == 0:
  113. alpha = max(0.001 * alpha_upper, (alpha_lower * alpha_upper)**0.5)
  114. else:
  115. alpha = initial_alpha
  116. for it in range(max_iter):
  117. if alpha < alpha_lower or alpha > alpha_upper:
  118. alpha = max(0.001 * alpha_upper, (alpha_lower * alpha_upper)**0.5)
  119. phi, phi_prime = phi_and_derivative(alpha, suf, s, Delta)
  120. if phi < 0:
  121. alpha_upper = alpha
  122. ratio = phi / phi_prime
  123. alpha_lower = max(alpha_lower, alpha - ratio)
  124. alpha -= (phi + Delta) * ratio / Delta
  125. if np.abs(phi) < rtol * Delta:
  126. break
  127. p = -V.dot(suf / (s**2 + alpha))
  128. # Make the norm of p equal to Delta, p is changed only slightly during
  129. # this. It is done to prevent p lie outside the trust region (which can
  130. # cause problems later).
  131. p *= Delta / norm(p)
  132. return p, alpha, it + 1
  133. def solve_trust_region_2d(B, g, Delta):
  134. """Solve a general trust-region problem in 2 dimensions.
  135. The problem is reformulated as a 4-th order algebraic equation,
  136. the solution of which is found by numpy.roots.
  137. Parameters
  138. ----------
  139. B : ndarray, shape (2, 2)
  140. Symmetric matrix, defines a quadratic term of the function.
  141. g : ndarray, shape (2,)
  142. Defines a linear term of the function.
  143. Delta : float
  144. Radius of a trust region.
  145. Returns
  146. -------
  147. p : ndarray, shape (2,)
  148. Found solution.
  149. newton_step : bool
  150. Whether the returned solution is the Newton step which lies within
  151. the trust region.
  152. """
  153. try:
  154. R, lower = cho_factor(B)
  155. p = -cho_solve((R, lower), g)
  156. if np.dot(p, p) <= Delta**2:
  157. return p, True
  158. except LinAlgError:
  159. pass
  160. a = B[0, 0] * Delta**2
  161. b = B[0, 1] * Delta**2
  162. c = B[1, 1] * Delta**2
  163. d = g[0] * Delta
  164. f = g[1] * Delta
  165. coeffs = np.array(
  166. [-b + d, 2 * (a - c + f), 6 * b, 2 * (-a + c + f), -b - d])
  167. t = np.roots(coeffs) # Can handle leading zeros.
  168. t = np.real(t[np.isreal(t)])
  169. p = Delta * np.vstack((2 * t / (1 + t**2), (1 - t**2) / (1 + t**2)))
  170. value = 0.5 * np.sum(p * B.dot(p), axis=0) + np.dot(g, p)
  171. i = np.argmin(value)
  172. p = p[:, i]
  173. return p, False
  174. def update_tr_radius(Delta, actual_reduction, predicted_reduction,
  175. step_norm, bound_hit):
  176. """Update the radius of a trust region based on the cost reduction.
  177. Returns
  178. -------
  179. Delta : float
  180. New radius.
  181. ratio : float
  182. Ratio between actual and predicted reductions. Zero if predicted
  183. reduction is zero.
  184. """
  185. if predicted_reduction > 0:
  186. ratio = actual_reduction / predicted_reduction
  187. else:
  188. ratio = 0
  189. if ratio < 0.25:
  190. Delta = 0.25 * step_norm
  191. elif ratio > 0.75 and bound_hit:
  192. Delta *= 2.0
  193. return Delta, ratio
  194. # Construction and minimization of quadratic functions.
  195. def build_quadratic_1d(J, g, s, diag=None, s0=None):
  196. """Parameterize a multivariate quadratic function along a line.
  197. The resulting univariate quadratic function is given as follows:
  198. ::
  199. f(t) = 0.5 * (s0 + s*t).T * (J.T*J + diag) * (s0 + s*t) +
  200. g.T * (s0 + s*t)
  201. Parameters
  202. ----------
  203. J : ndarray, sparse matrix or LinearOperator shape (m, n)
  204. Jacobian matrix, affects the quadratic term.
  205. g : ndarray, shape (n,)
  206. Gradient, defines the linear term.
  207. s : ndarray, shape (n,)
  208. Direction vector of a line.
  209. diag : None or ndarray with shape (n,), optional
  210. Addition diagonal part, affects the quadratic term.
  211. If None, assumed to be 0.
  212. s0 : None or ndarray with shape (n,), optional
  213. Initial point. If None, assumed to be 0.
  214. Returns
  215. -------
  216. a : float
  217. Coefficient for t**2.
  218. b : float
  219. Coefficient for t.
  220. c : float
  221. Free term. Returned only if `s0` is provided.
  222. """
  223. v = J.dot(s)
  224. a = np.dot(v, v)
  225. if diag is not None:
  226. a += np.dot(s * diag, s)
  227. a *= 0.5
  228. b = np.dot(g, s)
  229. if s0 is not None:
  230. u = J.dot(s0)
  231. b += np.dot(u, v)
  232. c = 0.5 * np.dot(u, u) + np.dot(g, s0)
  233. if diag is not None:
  234. b += np.dot(s0 * diag, s)
  235. c += 0.5 * np.dot(s0 * diag, s0)
  236. return a, b, c
  237. else:
  238. return a, b
  239. def minimize_quadratic_1d(a, b, lb, ub, c=0):
  240. """Minimize a 1-d quadratic function subject to bounds.
  241. The free term `c` is 0 by default. Bounds must be finite.
  242. Returns
  243. -------
  244. t : float
  245. Minimum point.
  246. y : float
  247. Minimum value.
  248. """
  249. t = [lb, ub]
  250. if a != 0:
  251. extremum = -0.5 * b / a
  252. if lb < extremum < ub:
  253. t.append(extremum)
  254. t = np.asarray(t)
  255. y = a * t**2 + b * t + c
  256. min_index = np.argmin(y)
  257. return t[min_index], y[min_index]
  258. def evaluate_quadratic(J, g, s, diag=None):
  259. """Compute values of a quadratic function arising in least squares.
  260. The function is 0.5 * s.T * (J.T * J + diag) * s + g.T * s.
  261. Parameters
  262. ----------
  263. J : ndarray, sparse matrix or LinearOperator, shape (m, n)
  264. Jacobian matrix, affects the quadratic term.
  265. g : ndarray, shape (n,)
  266. Gradient, defines the linear term.
  267. s : ndarray, shape (k, n) or (n,)
  268. Array containing steps as rows.
  269. diag : ndarray, shape (n,), optional
  270. Addition diagonal part, affects the quadratic term.
  271. If None, assumed to be 0.
  272. Returns
  273. -------
  274. values : ndarray with shape (k,) or float
  275. Values of the function. If `s` was 2-dimensional then ndarray is
  276. returned, otherwise float is returned.
  277. """
  278. if s.ndim == 1:
  279. Js = J.dot(s)
  280. q = np.dot(Js, Js)
  281. if diag is not None:
  282. q += np.dot(s * diag, s)
  283. else:
  284. Js = J.dot(s.T)
  285. q = np.sum(Js**2, axis=0)
  286. if diag is not None:
  287. q += np.sum(diag * s**2, axis=1)
  288. l = np.dot(s, g)
  289. return 0.5 * q + l
  290. # Utility functions to work with bound constraints.
  291. def in_bounds(x, lb, ub):
  292. """Check if a point lies within bounds."""
  293. return np.all((x >= lb) & (x <= ub))
  294. def step_size_to_bound(x, s, lb, ub):
  295. """Compute a min_step size required to reach a bound.
  296. The function computes a positive scalar t, such that x + s * t is on
  297. the bound.
  298. Returns
  299. -------
  300. step : float
  301. Computed step. Non-negative value.
  302. hits : ndarray of int with shape of x
  303. Each element indicates whether a corresponding variable reaches the
  304. bound:
  305. * 0 - the bound was not hit.
  306. * -1 - the lower bound was hit.
  307. * 1 - the upper bound was hit.
  308. """
  309. non_zero = np.nonzero(s)
  310. s_non_zero = s[non_zero]
  311. steps = np.empty_like(x)
  312. steps.fill(np.inf)
  313. with np.errstate(over='ignore'):
  314. steps[non_zero] = np.maximum((lb - x)[non_zero] / s_non_zero,
  315. (ub - x)[non_zero] / s_non_zero)
  316. min_step = np.min(steps)
  317. return min_step, np.equal(steps, min_step) * np.sign(s).astype(int)
  318. def find_active_constraints(x, lb, ub, rtol=1e-10):
  319. """Determine which constraints are active in a given point.
  320. The threshold is computed using `rtol` and the absolute value of the
  321. closest bound.
  322. Returns
  323. -------
  324. active : ndarray of int with shape of x
  325. Each component shows whether the corresponding constraint is active:
  326. * 0 - a constraint is not active.
  327. * -1 - a lower bound is active.
  328. * 1 - a upper bound is active.
  329. """
  330. active = np.zeros_like(x, dtype=int)
  331. if rtol == 0:
  332. active[x <= lb] = -1
  333. active[x >= ub] = 1
  334. return active
  335. lower_dist = x - lb
  336. upper_dist = ub - x
  337. lower_threshold = rtol * np.maximum(1, np.abs(lb))
  338. upper_threshold = rtol * np.maximum(1, np.abs(ub))
  339. lower_active = (np.isfinite(lb) &
  340. (lower_dist <= np.minimum(upper_dist, lower_threshold)))
  341. active[lower_active] = -1
  342. upper_active = (np.isfinite(ub) &
  343. (upper_dist <= np.minimum(lower_dist, upper_threshold)))
  344. active[upper_active] = 1
  345. return active
  346. def make_strictly_feasible(x, lb, ub, rstep=1e-10):
  347. """Shift a point to the interior of a feasible region.
  348. Each element of the returned vector is at least at a relative distance
  349. `rstep` from the closest bound. If ``rstep=0`` then `np.nextafter` is used.
  350. """
  351. x_new = x.copy()
  352. active = find_active_constraints(x, lb, ub, rstep)
  353. lower_mask = np.equal(active, -1)
  354. upper_mask = np.equal(active, 1)
  355. if rstep == 0:
  356. x_new[lower_mask] = np.nextafter(lb[lower_mask], ub[lower_mask])
  357. x_new[upper_mask] = np.nextafter(ub[upper_mask], lb[upper_mask])
  358. else:
  359. x_new[lower_mask] = (lb[lower_mask] +
  360. rstep * np.maximum(1, np.abs(lb[lower_mask])))
  361. x_new[upper_mask] = (ub[upper_mask] -
  362. rstep * np.maximum(1, np.abs(ub[upper_mask])))
  363. tight_bounds = (x_new < lb) | (x_new > ub)
  364. x_new[tight_bounds] = 0.5 * (lb[tight_bounds] + ub[tight_bounds])
  365. return x_new
  366. def CL_scaling_vector(x, g, lb, ub):
  367. """Compute Coleman-Li scaling vector and its derivatives.
  368. Components of a vector v are defined as follows:
  369. ::
  370. | ub[i] - x[i], if g[i] < 0 and ub[i] < np.inf
  371. v[i] = | x[i] - lb[i], if g[i] > 0 and lb[i] > -np.inf
  372. | 1, otherwise
  373. According to this definition v[i] >= 0 for all i. It differs from the
  374. definition in paper [1]_ (eq. (2.2)), where the absolute value of v is
  375. used. Both definitions are equivalent down the line.
  376. Derivatives of v with respect to x take value 1, -1 or 0 depending on a
  377. case.
  378. Returns
  379. -------
  380. v : ndarray with shape of x
  381. Scaling vector.
  382. dv : ndarray with shape of x
  383. Derivatives of v[i] with respect to x[i], diagonal elements of v's
  384. Jacobian.
  385. References
  386. ----------
  387. .. [1] M.A. Branch, T.F. Coleman, and Y. Li, "A Subspace, Interior,
  388. and Conjugate Gradient Method for Large-Scale Bound-Constrained
  389. Minimization Problems," SIAM Journal on Scientific Computing,
  390. Vol. 21, Number 1, pp 1-23, 1999.
  391. """
  392. v = np.ones_like(x)
  393. dv = np.zeros_like(x)
  394. mask = (g < 0) & np.isfinite(ub)
  395. v[mask] = ub[mask] - x[mask]
  396. dv[mask] = -1
  397. mask = (g > 0) & np.isfinite(lb)
  398. v[mask] = x[mask] - lb[mask]
  399. dv[mask] = 1
  400. return v, dv
  401. def reflective_transformation(y, lb, ub):
  402. """Compute reflective transformation and its gradient."""
  403. if in_bounds(y, lb, ub):
  404. return y, np.ones_like(y)
  405. lb_finite = np.isfinite(lb)
  406. ub_finite = np.isfinite(ub)
  407. x = y.copy()
  408. g_negative = np.zeros_like(y, dtype=bool)
  409. mask = lb_finite & ~ub_finite
  410. x[mask] = np.maximum(y[mask], 2 * lb[mask] - y[mask])
  411. g_negative[mask] = y[mask] < lb[mask]
  412. mask = ~lb_finite & ub_finite
  413. x[mask] = np.minimum(y[mask], 2 * ub[mask] - y[mask])
  414. g_negative[mask] = y[mask] > ub[mask]
  415. mask = lb_finite & ub_finite
  416. d = ub - lb
  417. t = np.remainder(y[mask] - lb[mask], 2 * d[mask])
  418. x[mask] = lb[mask] + np.minimum(t, 2 * d[mask] - t)
  419. g_negative[mask] = t > d[mask]
  420. g = np.ones_like(y)
  421. g[g_negative] = -1
  422. return x, g
  423. # Functions to display algorithm's progress.
  424. def print_header_nonlinear():
  425. print("{0:^15}{1:^15}{2:^15}{3:^15}{4:^15}{5:^15}"
  426. .format("Iteration", "Total nfev", "Cost", "Cost reduction",
  427. "Step norm", "Optimality"))
  428. def print_iteration_nonlinear(iteration, nfev, cost, cost_reduction,
  429. step_norm, optimality):
  430. if cost_reduction is None:
  431. cost_reduction = " " * 15
  432. else:
  433. cost_reduction = "{0:^15.2e}".format(cost_reduction)
  434. if step_norm is None:
  435. step_norm = " " * 15
  436. else:
  437. step_norm = "{0:^15.2e}".format(step_norm)
  438. print("{0:^15}{1:^15}{2:^15.4e}{3}{4}{5:^15.2e}"
  439. .format(iteration, nfev, cost, cost_reduction,
  440. step_norm, optimality))
  441. def print_header_linear():
  442. print("{0:^15}{1:^15}{2:^15}{3:^15}{4:^15}"
  443. .format("Iteration", "Cost", "Cost reduction", "Step norm",
  444. "Optimality"))
  445. def print_iteration_linear(iteration, cost, cost_reduction, step_norm,
  446. optimality):
  447. if cost_reduction is None:
  448. cost_reduction = " " * 15
  449. else:
  450. cost_reduction = "{0:^15.2e}".format(cost_reduction)
  451. if step_norm is None:
  452. step_norm = " " * 15
  453. else:
  454. step_norm = "{0:^15.2e}".format(step_norm)
  455. print("{0:^15}{1:^15.4e}{2}{3}{4:^15.2e}".format(
  456. iteration, cost, cost_reduction, step_norm, optimality))
  457. # Simple helper functions.
  458. def compute_grad(J, f):
  459. """Compute gradient of the least-squares cost function."""
  460. if isinstance(J, LinearOperator):
  461. return J.rmatvec(f)
  462. else:
  463. return J.T.dot(f)
  464. def compute_jac_scale(J, scale_inv_old=None):
  465. """Compute variables scale based on the Jacobian matrix."""
  466. if issparse(J):
  467. scale_inv = np.asarray(J.power(2).sum(axis=0)).ravel()**0.5
  468. else:
  469. scale_inv = np.sum(J**2, axis=0)**0.5
  470. if scale_inv_old is None:
  471. scale_inv[scale_inv == 0] = 1
  472. else:
  473. scale_inv = np.maximum(scale_inv, scale_inv_old)
  474. return 1 / scale_inv, scale_inv
  475. def left_multiplied_operator(J, d):
  476. """Return diag(d) J as LinearOperator."""
  477. J = aslinearoperator(J)
  478. def matvec(x):
  479. return d * J.matvec(x)
  480. def matmat(X):
  481. return d[:, np.newaxis] * J.matmat(X)
  482. def rmatvec(x):
  483. return J.rmatvec(x.ravel() * d)
  484. return LinearOperator(J.shape, matvec=matvec, matmat=matmat,
  485. rmatvec=rmatvec)
  486. def right_multiplied_operator(J, d):
  487. """Return J diag(d) as LinearOperator."""
  488. J = aslinearoperator(J)
  489. def matvec(x):
  490. return J.matvec(np.ravel(x) * d)
  491. def matmat(X):
  492. return J.matmat(X * d[:, np.newaxis])
  493. def rmatvec(x):
  494. return d * J.rmatvec(x)
  495. return LinearOperator(J.shape, matvec=matvec, matmat=matmat,
  496. rmatvec=rmatvec)
  497. def regularized_lsq_operator(J, diag):
  498. """Return a matrix arising in regularized least squares as LinearOperator.
  499. The matrix is
  500. [ J ]
  501. [ D ]
  502. where D is diagonal matrix with elements from `diag`.
  503. """
  504. J = aslinearoperator(J)
  505. m, n = J.shape
  506. def matvec(x):
  507. return np.hstack((J.matvec(x), diag * x))
  508. def rmatvec(x):
  509. x1 = x[:m]
  510. x2 = x[m:]
  511. return J.rmatvec(x1) + diag * x2
  512. return LinearOperator((m + n, n), matvec=matvec, rmatvec=rmatvec)
  513. def right_multiply(J, d, copy=True):
  514. """Compute J diag(d).
  515. If `copy` is False, `J` is modified in place (unless being LinearOperator).
  516. """
  517. if copy and not isinstance(J, LinearOperator):
  518. J = J.copy()
  519. if issparse(J):
  520. J.data *= d.take(J.indices, mode='clip') # scikit-learn recipe.
  521. elif isinstance(J, LinearOperator):
  522. J = right_multiplied_operator(J, d)
  523. else:
  524. J *= d
  525. return J
  526. def left_multiply(J, d, copy=True):
  527. """Compute diag(d) J.
  528. If `copy` is False, `J` is modified in place (unless being LinearOperator).
  529. """
  530. if copy and not isinstance(J, LinearOperator):
  531. J = J.copy()
  532. if issparse(J):
  533. J.data *= np.repeat(d, np.diff(J.indptr)) # scikit-learn recipe.
  534. elif isinstance(J, LinearOperator):
  535. J = left_multiplied_operator(J, d)
  536. else:
  537. J *= d[:, np.newaxis]
  538. return J
  539. def check_termination(dF, F, dx_norm, x_norm, ratio, ftol, xtol):
  540. """Check termination condition for nonlinear least squares."""
  541. ftol_satisfied = dF < ftol * F and ratio > 0.25
  542. xtol_satisfied = dx_norm < xtol * (xtol + x_norm)
  543. if ftol_satisfied and xtol_satisfied:
  544. return 4
  545. elif ftol_satisfied:
  546. return 2
  547. elif xtol_satisfied:
  548. return 3
  549. else:
  550. return None
  551. def scale_for_robust_loss_function(J, f, rho):
  552. """Scale Jacobian and residuals for a robust loss function.
  553. Arrays are modified in place.
  554. """
  555. J_scale = rho[1] + 2 * rho[2] * f**2
  556. J_scale[J_scale < EPS] = EPS
  557. J_scale **= 0.5
  558. f *= rho[1] / J_scale
  559. return left_multiply(J, J_scale, copy=False), f