common.py 20 KB

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