_lbfgsb_py.py 23 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634
  1. """
  2. Functions
  3. ---------
  4. .. autosummary::
  5. :toctree: generated/
  6. fmin_l_bfgs_b
  7. """
  8. ## License for the Python wrapper
  9. ## ==============================
  10. ## Copyright (c) 2004 David M. Cooke <cookedm@physics.mcmaster.ca>
  11. ## Permission is hereby granted, free of charge, to any person obtaining a
  12. ## copy of this software and associated documentation files (the "Software"),
  13. ## to deal in the Software without restriction, including without limitation
  14. ## the rights to use, copy, modify, merge, publish, distribute, sublicense,
  15. ## and/or sell copies of the Software, and to permit persons to whom the
  16. ## Software is furnished to do so, subject to the following conditions:
  17. ## The above copyright notice and this permission notice shall be included in
  18. ## all copies or substantial portions of the Software.
  19. ## THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
  20. ## IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
  21. ## FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
  22. ## AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
  23. ## LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
  24. ## FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
  25. ## DEALINGS IN THE SOFTWARE.
  26. ## Modifications by Travis Oliphant and Enthought, Inc. for inclusion in SciPy
  27. import numpy as np
  28. from numpy import array, asarray, float64, zeros
  29. from . import _lbfgsb
  30. from ._optimize import (MemoizeJac, OptimizeResult, _call_callback_maybe_halt,
  31. _wrap_callback, _check_unknown_options,
  32. _prepare_scalar_function)
  33. from ._constraints import old_bound_to_new
  34. from scipy.sparse.linalg import LinearOperator
  35. from scipy._lib.deprecation import _NoValue
  36. import warnings
  37. __all__ = ['fmin_l_bfgs_b', 'LbfgsInvHessProduct']
  38. status_messages = {
  39. 0 : "START",
  40. 1 : "NEW_X",
  41. 2 : "RESTART",
  42. 3 : "FG",
  43. 4 : "CONVERGENCE",
  44. 5 : "STOP",
  45. 6 : "WARNING",
  46. 7 : "ERROR",
  47. 8 : "ABNORMAL"
  48. }
  49. task_messages = {
  50. 0 : "",
  51. 301 : "",
  52. 302 : "",
  53. 401 : "NORM OF PROJECTED GRADIENT <= PGTOL",
  54. 402 : "RELATIVE REDUCTION OF F <= FACTR*EPSMCH",
  55. 501 : "CPU EXCEEDING THE TIME LIMIT",
  56. 502 : "TOTAL NO. OF F,G EVALUATIONS EXCEEDS LIMIT",
  57. 503 : "PROJECTED GRADIENT IS SUFFICIENTLY SMALL",
  58. 504 : "TOTAL NO. OF ITERATIONS REACHED LIMIT",
  59. 505 : "CALLBACK REQUESTED HALT",
  60. 601 : "ROUNDING ERRORS PREVENT PROGRESS",
  61. 602 : "STP = STPMAX",
  62. 603 : "STP = STPMIN",
  63. 604 : "XTOL TEST SATISFIED",
  64. 701 : "NO FEASIBLE SOLUTION",
  65. 702 : "FACTR < 0",
  66. 703 : "FTOL < 0",
  67. 704 : "GTOL < 0",
  68. 705 : "XTOL < 0",
  69. 706 : "STP < STPMIN",
  70. 707 : "STP > STPMAX",
  71. 708 : "STPMIN < 0",
  72. 709 : "STPMAX < STPMIN",
  73. 710 : "INITIAL G >= 0",
  74. 711 : "M <= 0",
  75. 712 : "N <= 0",
  76. 713 : "INVALID NBD",
  77. }
  78. def fmin_l_bfgs_b(func, x0, fprime=None, args=(),
  79. approx_grad=0,
  80. bounds=None, m=10, factr=1e7, pgtol=1e-5,
  81. epsilon=1e-8,
  82. iprint=_NoValue, maxfun=15000, maxiter=15000, disp=_NoValue,
  83. callback=None, maxls=20):
  84. """
  85. Minimize a function func using the L-BFGS-B algorithm.
  86. Parameters
  87. ----------
  88. func : callable f(x,*args)
  89. Function to minimize.
  90. x0 : ndarray
  91. Initial guess.
  92. fprime : callable fprime(x,*args), optional
  93. The gradient of `func`. If None, then `func` returns the function
  94. value and the gradient (``f, g = func(x, *args)``), unless
  95. `approx_grad` is True in which case `func` returns only ``f``.
  96. args : sequence, optional
  97. Arguments to pass to `func` and `fprime`.
  98. approx_grad : bool, optional
  99. Whether to approximate the gradient numerically (in which case
  100. `func` returns only the function value).
  101. bounds : list, optional
  102. ``(min, max)`` pairs for each element in ``x``, defining
  103. the bounds on that parameter. Use None or +-inf for one of ``min`` or
  104. ``max`` when there is no bound in that direction.
  105. m : int, optional
  106. The maximum number of variable metric corrections
  107. used to define the limited memory matrix. (The limited memory BFGS
  108. method does not store the full hessian but uses this many terms in an
  109. approximation to it.)
  110. factr : float, optional
  111. The iteration stops when
  112. ``(f^k - f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= factr * eps``,
  113. where ``eps`` is the machine precision, which is automatically
  114. generated by the code. Typical values for `factr` are: 1e12 for
  115. low accuracy; 1e7 for moderate accuracy; 10.0 for extremely
  116. high accuracy. See Notes for relationship to `ftol`, which is exposed
  117. (instead of `factr`) by the `scipy.optimize.minimize` interface to
  118. L-BFGS-B.
  119. pgtol : float, optional
  120. The iteration will stop when
  121. ``max{|proj g_i | i = 1, ..., n} <= pgtol``
  122. where ``proj g_i`` is the i-th component of the projected gradient.
  123. epsilon : float, optional
  124. Step size used when `approx_grad` is True, for numerically
  125. calculating the gradient
  126. iprint : int, optional
  127. Deprecated option that previously controlled the text printed on the
  128. screen during the problem solution. Now the code does not emit any
  129. output and this keyword has no function.
  130. .. deprecated:: 1.15.0
  131. This keyword is deprecated and will be removed from SciPy 1.18.0.
  132. disp : int, optional
  133. Deprecated option that previously controlled the text printed on the
  134. screen during the problem solution. Now the code does not emit any
  135. output and this keyword has no function.
  136. .. deprecated:: 1.15.0
  137. This keyword is deprecated and will be removed from SciPy 1.18.0.
  138. maxfun : int, optional
  139. Maximum number of function evaluations. Note that this function
  140. may violate the limit because of evaluating gradients by numerical
  141. differentiation.
  142. maxiter : int, optional
  143. Maximum number of iterations.
  144. callback : callable, optional
  145. Called after each iteration, as ``callback(xk)``, where ``xk`` is the
  146. current parameter vector.
  147. maxls : int, optional
  148. Maximum number of line search steps (per iteration). Default is 20.
  149. Returns
  150. -------
  151. x : array_like
  152. Estimated position of the minimum.
  153. f : float
  154. Value of `func` at the minimum.
  155. d : dict
  156. Information dictionary.
  157. * d['warnflag'] is
  158. - 0 if converged,
  159. - 1 if too many function evaluations or too many iterations,
  160. - 2 if stopped for another reason, given in d['task']
  161. * d['grad'] is the gradient at the minimum (should be 0 ish)
  162. * d['funcalls'] is the number of function calls made.
  163. * d['nit'] is the number of iterations.
  164. See also
  165. --------
  166. minimize: Interface to minimization algorithms for multivariate
  167. functions. See the 'L-BFGS-B' `method` in particular. Note that the
  168. `ftol` option is made available via that interface, while `factr` is
  169. provided via this interface, where `factr` is the factor multiplying
  170. the default machine floating-point precision to arrive at `ftol`:
  171. ``ftol = factr * numpy.finfo(float).eps``.
  172. Notes
  173. -----
  174. SciPy uses a C-translated and modified version of the Fortran code,
  175. L-BFGS-B v3.0 (released April 25, 2011, BSD-3 licensed). Original Fortran
  176. version was written by Ciyou Zhu, Richard Byrd, Jorge Nocedal and,
  177. Jose Luis Morales.
  178. References
  179. ----------
  180. * R. H. Byrd, P. Lu and J. Nocedal. A Limited Memory Algorithm for Bound
  181. Constrained Optimization, (1995), SIAM Journal on Scientific and
  182. Statistical Computing, 16, 5, pp. 1190-1208.
  183. * C. Zhu, R. H. Byrd and J. Nocedal. L-BFGS-B: Algorithm 778: L-BFGS-B,
  184. FORTRAN routines for large scale bound constrained optimization (1997),
  185. ACM Transactions on Mathematical Software, 23, 4, pp. 550 - 560.
  186. * J.L. Morales and J. Nocedal. L-BFGS-B: Remark on Algorithm 778: L-BFGS-B,
  187. FORTRAN routines for large scale bound constrained optimization (2011),
  188. ACM Transactions on Mathematical Software, 38, 1.
  189. Examples
  190. --------
  191. Solve a linear regression problem via `fmin_l_bfgs_b`. To do this, first we
  192. define an objective function ``f(m, b) = (y - y_model)**2``, where `y`
  193. describes the observations and `y_model` the prediction of the linear model
  194. as ``y_model = m*x + b``. The bounds for the parameters, ``m`` and ``b``,
  195. are arbitrarily chosen as ``(0,5)`` and ``(5,10)`` for this example.
  196. >>> import numpy as np
  197. >>> from scipy.optimize import fmin_l_bfgs_b
  198. >>> X = np.arange(0, 10, 1)
  199. >>> M = 2
  200. >>> B = 3
  201. >>> Y = M * X + B
  202. >>> def func(parameters, *args):
  203. ... x = args[0]
  204. ... y = args[1]
  205. ... m, b = parameters
  206. ... y_model = m*x + b
  207. ... error = sum(np.power((y - y_model), 2))
  208. ... return error
  209. >>> initial_values = np.array([0.0, 1.0])
  210. >>> x_opt, f_opt, info = fmin_l_bfgs_b(func, x0=initial_values, args=(X, Y),
  211. ... approx_grad=True)
  212. >>> x_opt, f_opt
  213. array([1.99999999, 3.00000006]), 1.7746231151323805e-14 # may vary
  214. The optimized parameters in ``x_opt`` agree with the ground truth parameters
  215. ``m`` and ``b``. Next, let us perform a bound constrained optimization using
  216. the `bounds` parameter.
  217. >>> bounds = [(0, 5), (5, 10)]
  218. >>> x_opt, f_op, info = fmin_l_bfgs_b(func, x0=initial_values, args=(X, Y),
  219. ... approx_grad=True, bounds=bounds)
  220. >>> x_opt, f_opt
  221. array([1.65990508, 5.31649385]), 15.721334516453945 # may vary
  222. """
  223. # handle fprime/approx_grad
  224. if approx_grad:
  225. fun = func
  226. jac = None
  227. elif fprime is None:
  228. fun = MemoizeJac(func)
  229. jac = fun.derivative
  230. else:
  231. fun = func
  232. jac = fprime
  233. # build options
  234. callback = _wrap_callback(callback)
  235. opts = {'disp': disp,
  236. 'iprint': iprint,
  237. 'maxcor': m,
  238. 'ftol': factr * np.finfo(float).eps,
  239. 'gtol': pgtol,
  240. 'eps': epsilon,
  241. 'maxfun': maxfun,
  242. 'maxiter': maxiter,
  243. 'callback': callback,
  244. 'maxls': maxls}
  245. res = _minimize_lbfgsb(fun, x0, args=args, jac=jac, bounds=bounds,
  246. **opts)
  247. d = {'grad': res['jac'],
  248. 'task': res['message'],
  249. 'funcalls': res['nfev'],
  250. 'nit': res['nit'],
  251. 'warnflag': res['status']}
  252. f = res['fun']
  253. x = res['x']
  254. return x, f, d
  255. def _minimize_lbfgsb(fun, x0, args=(), jac=None, bounds=None,
  256. disp=_NoValue, maxcor=10, ftol=2.2204460492503131e-09,
  257. gtol=1e-5, eps=1e-8, maxfun=15000, maxiter=15000,
  258. iprint=_NoValue, callback=None, maxls=20,
  259. finite_diff_rel_step=None, workers=None,
  260. **unknown_options):
  261. """
  262. Minimize a scalar function of one or more variables using the L-BFGS-B
  263. algorithm.
  264. Options
  265. -------
  266. disp : None or int
  267. Deprecated option that previously controlled the text printed on the
  268. screen during the problem solution. Now the code does not emit any
  269. output and this keyword has no function.
  270. .. deprecated:: 1.15.0
  271. This keyword is deprecated and will be removed from SciPy 1.18.0.
  272. maxcor : int
  273. The maximum number of variable metric corrections used to
  274. define the limited memory matrix. (The limited memory BFGS
  275. method does not store the full hessian but uses this many terms
  276. in an approximation to it.)
  277. ftol : float
  278. The iteration stops when ``(f^k -
  279. f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= ftol``.
  280. gtol : float
  281. The iteration will stop when ``max{|proj g_i | i = 1, ..., n}
  282. <= gtol`` where ``proj g_i`` is the i-th component of the
  283. projected gradient.
  284. eps : float or ndarray
  285. If `jac is None` the absolute step size used for numerical
  286. approximation of the jacobian via forward differences.
  287. maxfun : int
  288. Maximum number of function evaluations before minimization terminates.
  289. Note that this function may violate the limit if the gradients
  290. are evaluated by numerical differentiation.
  291. maxiter : int
  292. Maximum number of algorithm iterations.
  293. iprint : int, optional
  294. Deprecated option that previously controlled the text printed on the
  295. screen during the problem solution. Now the code does not emit any
  296. output and this keyword has no function.
  297. .. deprecated:: 1.15.0
  298. This keyword is deprecated and will be removed from SciPy 1.18.0.
  299. maxls : int, optional
  300. Maximum number of line search steps (per iteration). Default is 20.
  301. finite_diff_rel_step : None or array_like, optional
  302. If ``jac in ['2-point', '3-point', 'cs']`` the relative step size to
  303. use for numerical approximation of the jacobian. The absolute step
  304. size is computed as ``h = rel_step * sign(x) * max(1, abs(x))``,
  305. possibly adjusted to fit into the bounds. For ``method='3-point'``
  306. the sign of `h` is ignored. If None (default) then step is selected
  307. automatically.
  308. workers : int, map-like callable, optional
  309. A map-like callable, such as `multiprocessing.Pool.map` for evaluating
  310. any numerical differentiation in parallel.
  311. This evaluation is carried out as ``workers(fun, iterable)``.
  312. .. versionadded:: 1.16.0
  313. Notes
  314. -----
  315. The option `ftol` is exposed via the `scipy.optimize.minimize` interface,
  316. but calling `scipy.optimize.fmin_l_bfgs_b` directly exposes `factr`. The
  317. relationship between the two is ``ftol = factr * numpy.finfo(float).eps``.
  318. I.e., `factr` multiplies the default machine floating-point precision to
  319. arrive at `ftol`.
  320. If the minimization is slow to converge the optimizer may halt if the
  321. total number of function evaluations exceeds `maxfun`, or the number of
  322. algorithm iterations has reached `maxiter` (whichever comes first). If
  323. this is the case then ``result.success=False``, and an appropriate
  324. error message is contained in ``result.message``.
  325. """
  326. _check_unknown_options(unknown_options)
  327. m = maxcor
  328. pgtol = gtol
  329. factr = ftol / np.finfo(float).eps
  330. x0 = asarray(x0).ravel()
  331. n, = x0.shape
  332. if disp is not _NoValue:
  333. warnings.warn("scipy.optimize: The `disp` and `iprint` options of the "
  334. "L-BFGS-B solver are deprecated and will be removed in "
  335. "SciPy 1.18.0.",
  336. DeprecationWarning, stacklevel=3)
  337. if iprint is not _NoValue:
  338. warnings.warn("scipy.optimize: The `disp` and `iprint` options of the "
  339. "L-BFGS-B solver are deprecated and will be removed in "
  340. "SciPy 1.18.0.",
  341. DeprecationWarning, stacklevel=3)
  342. # historically old-style bounds were/are expected by lbfgsb.
  343. # That's still the case but we'll deal with new-style from here on,
  344. # it's easier
  345. if bounds is None:
  346. pass
  347. elif len(bounds) != n:
  348. raise ValueError('length of x0 != length of bounds')
  349. else:
  350. bounds = np.array(old_bound_to_new(bounds))
  351. # check bounds
  352. if (bounds[0] > bounds[1]).any():
  353. raise ValueError(
  354. "LBFGSB - one of the lower bounds is greater than an upper bound."
  355. )
  356. # initial vector must lie within the bounds. Otherwise ScalarFunction and
  357. # approx_derivative will cause problems
  358. x0 = np.clip(x0, bounds[0], bounds[1])
  359. # _prepare_scalar_function can use bounds=None to represent no bounds
  360. sf = _prepare_scalar_function(fun, x0, jac=jac, args=args, epsilon=eps,
  361. bounds=bounds,
  362. finite_diff_rel_step=finite_diff_rel_step,
  363. workers=workers)
  364. func_and_grad = sf.fun_and_grad
  365. nbd = zeros(n, np.int32)
  366. low_bnd = zeros(n, float64)
  367. upper_bnd = zeros(n, float64)
  368. bounds_map = {(-np.inf, np.inf): 0,
  369. (1, np.inf): 1,
  370. (1, 1): 2,
  371. (-np.inf, 1): 3}
  372. if bounds is not None:
  373. for i in range(0, n):
  374. L, U = bounds[0, i], bounds[1, i]
  375. if not np.isinf(L):
  376. low_bnd[i] = L
  377. L = 1
  378. if not np.isinf(U):
  379. upper_bnd[i] = U
  380. U = 1
  381. nbd[i] = bounds_map[L, U]
  382. if not maxls > 0:
  383. raise ValueError('maxls must be positive.')
  384. x = array(x0, dtype=np.float64)
  385. f = array(0.0, dtype=np.float64)
  386. g = zeros((n,), dtype=np.float64)
  387. wa = zeros(2*m*n + 5*n + 11*m*m + 8*m, float64)
  388. iwa = zeros(3*n, dtype=np.int32)
  389. task = zeros(2, dtype=np.int32)
  390. ln_task = zeros(2, dtype=np.int32)
  391. lsave = zeros(4, dtype=np.int32)
  392. isave = zeros(44, dtype=np.int32)
  393. dsave = zeros(29, dtype=float64)
  394. n_iterations = 0
  395. while True:
  396. # g may become float32 if a user provides a function that calculates
  397. # the Jacobian in float32 (see gh-18730). The underlying code expects
  398. # float64, so upcast it
  399. g = g.astype(np.float64)
  400. # x, f, g, wa, iwa, task, csave, lsave, isave, dsave = \
  401. _lbfgsb.setulb(m, x, low_bnd, upper_bnd, nbd, f, g, factr, pgtol, wa,
  402. iwa, task, lsave, isave, dsave, maxls, ln_task)
  403. if task[0] == 3:
  404. # The minimization routine wants f and g at the current x.
  405. # Note that interruptions due to maxfun are postponed
  406. # until the completion of the current minimization iteration.
  407. # Overwrite f and g:
  408. f, g = func_and_grad(x)
  409. elif task[0] == 1:
  410. # new iteration
  411. n_iterations += 1
  412. intermediate_result = OptimizeResult(x=x, fun=f)
  413. if _call_callback_maybe_halt(callback, intermediate_result):
  414. task[0] = 5
  415. task[1] = 505
  416. if n_iterations >= maxiter:
  417. task[0] = 5
  418. task[1] = 504
  419. elif sf.nfev > maxfun:
  420. task[0] = 5
  421. task[1] = 502
  422. else:
  423. break
  424. if task[0] == 4:
  425. warnflag = 0
  426. elif sf.nfev > maxfun or n_iterations >= maxiter:
  427. warnflag = 1
  428. else:
  429. warnflag = 2
  430. # These two portions of the workspace are described in the mainlb
  431. # function docstring in "__lbfgsb.c", ws and wy arguments.
  432. s = wa[0: m*n].reshape(m, n)
  433. y = wa[m*n: 2*m*n].reshape(m, n)
  434. # isave(31) = the total number of BFGS updates prior the current iteration.
  435. n_bfgs_updates = isave[30]
  436. n_corrs = min(n_bfgs_updates, maxcor)
  437. hess_inv = LbfgsInvHessProduct(s[:n_corrs], y[:n_corrs])
  438. msg = status_messages[task[0]] + ": " + task_messages[task[1]]
  439. return OptimizeResult(fun=f, jac=g, nfev=sf.nfev,
  440. njev=sf.ngev,
  441. nit=n_iterations, status=warnflag, message=msg,
  442. x=x, success=(warnflag == 0), hess_inv=hess_inv)
  443. class LbfgsInvHessProduct(LinearOperator):
  444. """Linear operator for the L-BFGS approximate inverse Hessian.
  445. This operator computes the product of a vector with the approximate inverse
  446. of the Hessian of the objective function, using the L-BFGS limited
  447. memory approximation to the inverse Hessian, accumulated during the
  448. optimization.
  449. Objects of this class implement the ``scipy.sparse.linalg.LinearOperator``
  450. interface.
  451. Parameters
  452. ----------
  453. sk : array_like, shape=(n_corr, n)
  454. Array of `n_corr` most recent updates to the solution vector.
  455. (See [1]).
  456. yk : array_like, shape=(n_corr, n)
  457. Array of `n_corr` most recent updates to the gradient. (See [1]).
  458. References
  459. ----------
  460. .. [1] Nocedal, Jorge. "Updating quasi-Newton matrices with limited
  461. storage." Mathematics of computation 35.151 (1980): 773-782.
  462. """
  463. def __init__(self, sk, yk):
  464. """Construct the operator."""
  465. if sk.shape != yk.shape or sk.ndim != 2:
  466. raise ValueError('sk and yk must have matching shape, (n_corrs, n)')
  467. n_corrs, n = sk.shape
  468. super().__init__(dtype=np.float64, shape=(n, n))
  469. self.sk = sk
  470. self.yk = yk
  471. self.n_corrs = n_corrs
  472. self.rho = 1 / np.einsum('ij,ij->i', sk, yk)
  473. def _matvec(self, x):
  474. """Efficient matrix-vector multiply with the BFGS matrices.
  475. This calculation is described in Section (4) of [1].
  476. Parameters
  477. ----------
  478. x : ndarray
  479. An array with shape (n,) or (n,1).
  480. Returns
  481. -------
  482. y : ndarray
  483. The matrix-vector product
  484. """
  485. s, y, n_corrs, rho = self.sk, self.yk, self.n_corrs, self.rho
  486. q = np.array(x, dtype=self.dtype, copy=True)
  487. if q.ndim == 2 and q.shape[1] == 1:
  488. q = q.reshape(-1)
  489. alpha = np.empty(n_corrs)
  490. for i in range(n_corrs-1, -1, -1):
  491. alpha[i] = rho[i] * np.dot(s[i], q)
  492. q = q - alpha[i]*y[i]
  493. r = q
  494. for i in range(n_corrs):
  495. beta = rho[i] * np.dot(y[i], r)
  496. r = r + s[i] * (alpha[i] - beta)
  497. return r
  498. def _matmat(self, X):
  499. """Efficient matrix-matrix multiply with the BFGS matrices.
  500. This calculation is described in Section (4) of [1].
  501. Parameters
  502. ----------
  503. X : ndarray
  504. An array with shape (n,m)
  505. Returns
  506. -------
  507. Y : ndarray
  508. The matrix-matrix product
  509. Notes
  510. -----
  511. This implementation is written starting from _matvec and broadcasting
  512. all expressions along the second axis of X.
  513. """
  514. s, y, n_corrs, rho = self.sk, self.yk, self.n_corrs, self.rho
  515. Q = np.array(X, dtype=self.dtype, copy=True)
  516. alpha = np.empty((n_corrs, Q.shape[1]))
  517. for i in range(n_corrs-1, -1, -1):
  518. alpha[i] = rho[i] * np.dot(s[i], Q)
  519. Q -= alpha[i]*y[i][:, np.newaxis]
  520. R = Q
  521. for i in range(n_corrs):
  522. beta = rho[i] * np.dot(y[i], R)
  523. R += s[i][:, np.newaxis] * (alpha[i] - beta)
  524. return R
  525. def todense(self):
  526. """Return a dense array representation of this operator.
  527. Returns
  528. -------
  529. arr : ndarray, shape=(n, n)
  530. An array with the same shape and containing
  531. the same data represented by this `LinearOperator`.
  532. """
  533. I_arr = np.eye(*self.shape, dtype=self.dtype)
  534. return self._matmat(I_arr)