qp_subproblem.py 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637
  1. """Equality-constrained quadratic programming solvers."""
  2. from scipy.sparse import (linalg, bmat, csc_matrix)
  3. from math import copysign
  4. import numpy as np
  5. from numpy.linalg import norm
  6. __all__ = [
  7. 'eqp_kktfact',
  8. 'sphere_intersections',
  9. 'box_intersections',
  10. 'box_sphere_intersections',
  11. 'inside_box_boundaries',
  12. 'modified_dogleg',
  13. 'projected_cg'
  14. ]
  15. # For comparison with the projected CG
  16. def eqp_kktfact(H, c, A, b):
  17. """Solve equality-constrained quadratic programming (EQP) problem.
  18. Solve ``min 1/2 x.T H x + x.t c`` subject to ``A x + b = 0``
  19. using direct factorization of the KKT system.
  20. Parameters
  21. ----------
  22. H : sparse matrix, shape (n, n)
  23. Hessian matrix of the EQP problem.
  24. c : array_like, shape (n,)
  25. Gradient of the quadratic objective function.
  26. A : sparse matrix
  27. Jacobian matrix of the EQP problem.
  28. b : array_like, shape (m,)
  29. Right-hand side of the constraint equation.
  30. Returns
  31. -------
  32. x : array_like, shape (n,)
  33. Solution of the KKT problem.
  34. lagrange_multipliers : ndarray, shape (m,)
  35. Lagrange multipliers of the KKT problem.
  36. """
  37. n, = np.shape(c) # Number of parameters
  38. m, = np.shape(b) # Number of constraints
  39. # Karush-Kuhn-Tucker matrix of coefficients.
  40. # Defined as in Nocedal/Wright "Numerical
  41. # Optimization" p.452 in Eq. (16.4).
  42. kkt_matrix = csc_matrix(bmat([[H, A.T], [A, None]]))
  43. # Vector of coefficients.
  44. kkt_vec = np.hstack([-c, -b])
  45. # TODO: Use a symmetric indefinite factorization
  46. # to solve the system twice as fast (because
  47. # of the symmetry).
  48. lu = linalg.splu(kkt_matrix)
  49. kkt_sol = lu.solve(kkt_vec)
  50. x = kkt_sol[:n]
  51. lagrange_multipliers = -kkt_sol[n:n+m]
  52. return x, lagrange_multipliers
  53. def sphere_intersections(z, d, trust_radius,
  54. entire_line=False):
  55. """Find the intersection between segment (or line) and spherical constraints.
  56. Find the intersection between the segment (or line) defined by the
  57. parametric equation ``x(t) = z + t*d`` and the ball
  58. ``||x|| <= trust_radius``.
  59. Parameters
  60. ----------
  61. z : array_like, shape (n,)
  62. Initial point.
  63. d : array_like, shape (n,)
  64. Direction.
  65. trust_radius : float
  66. Ball radius.
  67. entire_line : bool, optional
  68. When ``True``, the function returns the intersection between the line
  69. ``x(t) = z + t*d`` (``t`` can assume any value) and the ball
  70. ``||x|| <= trust_radius``. When ``False``, the function returns the intersection
  71. between the segment ``x(t) = z + t*d``, ``0 <= t <= 1``, and the ball.
  72. Returns
  73. -------
  74. ta, tb : float
  75. The line/segment ``x(t) = z + t*d`` is inside the ball for
  76. for ``ta <= t <= tb``.
  77. intersect : bool
  78. When ``True``, there is a intersection between the line/segment
  79. and the sphere. On the other hand, when ``False``, there is no
  80. intersection.
  81. """
  82. # Special case when d=0
  83. if norm(d) == 0:
  84. return 0, 0, False
  85. # Check for inf trust_radius
  86. if np.isinf(trust_radius):
  87. if entire_line:
  88. ta = -np.inf
  89. tb = np.inf
  90. else:
  91. ta = 0
  92. tb = 1
  93. intersect = True
  94. return ta, tb, intersect
  95. a = np.dot(d, d)
  96. b = 2 * np.dot(z, d)
  97. c = np.dot(z, z) - trust_radius**2
  98. discriminant = b*b - 4*a*c
  99. if discriminant < 0:
  100. intersect = False
  101. return 0, 0, intersect
  102. sqrt_discriminant = np.sqrt(discriminant)
  103. # The following calculation is mathematically
  104. # equivalent to:
  105. # ta = (-b - sqrt_discriminant) / (2*a)
  106. # tb = (-b + sqrt_discriminant) / (2*a)
  107. # but produce smaller round off errors.
  108. # Look at Matrix Computation p.97
  109. # for a better justification.
  110. aux = b + copysign(sqrt_discriminant, b)
  111. ta = -aux / (2*a)
  112. tb = -2*c / aux
  113. ta, tb = sorted([ta, tb])
  114. if entire_line:
  115. intersect = True
  116. else:
  117. # Checks to see if intersection happens
  118. # within vectors length.
  119. if tb < 0 or ta > 1:
  120. intersect = False
  121. ta = 0
  122. tb = 0
  123. else:
  124. intersect = True
  125. # Restrict intersection interval
  126. # between 0 and 1.
  127. ta = max(0, ta)
  128. tb = min(1, tb)
  129. return ta, tb, intersect
  130. def box_intersections(z, d, lb, ub,
  131. entire_line=False):
  132. """Find the intersection between segment (or line) and box constraints.
  133. Find the intersection between the segment (or line) defined by the
  134. parametric equation ``x(t) = z + t*d`` and the rectangular box
  135. ``lb <= x <= ub``.
  136. Parameters
  137. ----------
  138. z : array_like, shape (n,)
  139. Initial point.
  140. d : array_like, shape (n,)
  141. Direction.
  142. lb : array_like, shape (n,)
  143. Lower bounds to each one of the components of ``x``. Used
  144. to delimit the rectangular box.
  145. ub : array_like, shape (n, )
  146. Upper bounds to each one of the components of ``x``. Used
  147. to delimit the rectangular box.
  148. entire_line : bool, optional
  149. When ``True``, the function returns the intersection between the line
  150. ``x(t) = z + t*d`` (``t`` can assume any value) and the rectangular
  151. box. When ``False``, the function returns the intersection between the segment
  152. ``x(t) = z + t*d``, ``0 <= t <= 1``, and the rectangular box.
  153. Returns
  154. -------
  155. ta, tb : float
  156. The line/segment ``x(t) = z + t*d`` is inside the box for
  157. for ``ta <= t <= tb``.
  158. intersect : bool
  159. When ``True``, there is a intersection between the line (or segment)
  160. and the rectangular box. On the other hand, when ``False``, there is no
  161. intersection.
  162. """
  163. # Make sure it is a numpy array
  164. z = np.asarray(z)
  165. d = np.asarray(d)
  166. lb = np.asarray(lb)
  167. ub = np.asarray(ub)
  168. # Special case when d=0
  169. if norm(d) == 0:
  170. return 0, 0, False
  171. # Get values for which d==0
  172. zero_d = (d == 0)
  173. # If the boundaries are not satisfied for some coordinate
  174. # for which "d" is zero, there is no box-line intersection.
  175. if (z[zero_d] < lb[zero_d]).any() or (z[zero_d] > ub[zero_d]).any():
  176. intersect = False
  177. return 0, 0, intersect
  178. # Remove values for which d is zero
  179. not_zero_d = np.logical_not(zero_d)
  180. z = z[not_zero_d]
  181. d = d[not_zero_d]
  182. lb = lb[not_zero_d]
  183. ub = ub[not_zero_d]
  184. # Find a series of intervals (t_lb[i], t_ub[i]).
  185. t_lb = (lb-z) / d
  186. t_ub = (ub-z) / d
  187. # Get the intersection of all those intervals.
  188. ta = max(np.minimum(t_lb, t_ub))
  189. tb = min(np.maximum(t_lb, t_ub))
  190. # Check if intersection is feasible
  191. if ta <= tb:
  192. intersect = True
  193. else:
  194. intersect = False
  195. # Checks to see if intersection happens within vectors length.
  196. if not entire_line:
  197. if tb < 0 or ta > 1:
  198. intersect = False
  199. ta = 0
  200. tb = 0
  201. else:
  202. # Restrict intersection interval between 0 and 1.
  203. ta = max(0, ta)
  204. tb = min(1, tb)
  205. return ta, tb, intersect
  206. def box_sphere_intersections(z, d, lb, ub, trust_radius,
  207. entire_line=False,
  208. extra_info=False):
  209. """Find the intersection between segment (or line) and box/sphere constraints.
  210. Find the intersection between the segment (or line) defined by the
  211. parametric equation ``x(t) = z + t*d``, the rectangular box
  212. ``lb <= x <= ub`` and the ball ``||x|| <= trust_radius``.
  213. Parameters
  214. ----------
  215. z : array_like, shape (n,)
  216. Initial point.
  217. d : array_like, shape (n,)
  218. Direction.
  219. lb : array_like, shape (n,)
  220. Lower bounds to each one of the components of ``x``. Used
  221. to delimit the rectangular box.
  222. ub : array_like, shape (n, )
  223. Upper bounds to each one of the components of ``x``. Used
  224. to delimit the rectangular box.
  225. trust_radius : float
  226. Ball radius.
  227. entire_line : bool, optional
  228. When ``True``, the function returns the intersection between the line
  229. ``x(t) = z + t*d`` (``t`` can assume any value) and the constraints.
  230. When ``False``, the function returns the intersection between the segment
  231. ``x(t) = z + t*d``, ``0 <= t <= 1`` and the constraints.
  232. extra_info : bool, optional
  233. When ``True``, the function returns ``intersect_sphere`` and ``intersect_box``.
  234. Returns
  235. -------
  236. ta, tb : float
  237. The line/segment ``x(t) = z + t*d`` is inside the rectangular box and
  238. inside the ball for ``ta <= t <= tb``.
  239. intersect : bool
  240. When ``True``, there is a intersection between the line (or segment)
  241. and both constraints. On the other hand, when ``False``, there is no
  242. intersection.
  243. sphere_info : dict, optional
  244. Dictionary ``{ta, tb, intersect}`` containing the interval ``[ta, tb]``
  245. for which the line intercepts the ball. And a boolean value indicating
  246. whether the sphere is intersected by the line.
  247. box_info : dict, optional
  248. Dictionary ``{ta, tb, intersect}`` containing the interval ``[ta, tb]``
  249. for which the line intercepts the box. And a boolean value indicating
  250. whether the box is intersected by the line.
  251. """
  252. ta_b, tb_b, intersect_b = box_intersections(z, d, lb, ub,
  253. entire_line)
  254. ta_s, tb_s, intersect_s = sphere_intersections(z, d,
  255. trust_radius,
  256. entire_line)
  257. ta = np.maximum(ta_b, ta_s)
  258. tb = np.minimum(tb_b, tb_s)
  259. if intersect_b and intersect_s and ta <= tb:
  260. intersect = True
  261. else:
  262. intersect = False
  263. if extra_info:
  264. sphere_info = {'ta': ta_s, 'tb': tb_s, 'intersect': intersect_s}
  265. box_info = {'ta': ta_b, 'tb': tb_b, 'intersect': intersect_b}
  266. return ta, tb, intersect, sphere_info, box_info
  267. else:
  268. return ta, tb, intersect
  269. def inside_box_boundaries(x, lb, ub):
  270. """Check if lb <= x <= ub."""
  271. return (lb <= x).all() and (x <= ub).all()
  272. def reinforce_box_boundaries(x, lb, ub):
  273. """Return clipped value of x"""
  274. return np.minimum(np.maximum(x, lb), ub)
  275. def modified_dogleg(A, Y, b, trust_radius, lb, ub):
  276. """Approximately minimize ``1/2*|| A x + b ||^2`` inside trust-region.
  277. Approximately solve the problem of minimizing ``1/2*|| A x + b ||^2``
  278. subject to ``||x|| < Delta`` and ``lb <= x <= ub`` using a modification
  279. of the classical dogleg approach.
  280. Parameters
  281. ----------
  282. A : LinearOperator (or sparse matrix or ndarray), shape (m, n)
  283. Matrix ``A`` in the minimization problem. It should have
  284. dimension ``(m, n)`` such that ``m < n``.
  285. Y : LinearOperator (or sparse matrix or ndarray), shape (n, m)
  286. LinearOperator that apply the projection matrix
  287. ``Q = A.T inv(A A.T)`` to the vector. The obtained vector
  288. ``y = Q x`` being the minimum norm solution of ``A y = x``.
  289. b : array_like, shape (m,)
  290. Vector ``b``in the minimization problem.
  291. trust_radius: float
  292. Trust radius to be considered. Delimits a sphere boundary
  293. to the problem.
  294. lb : array_like, shape (n,)
  295. Lower bounds to each one of the components of ``x``.
  296. It is expected that ``lb <= 0``, otherwise the algorithm
  297. may fail. If ``lb[i] = -Inf``, the lower
  298. bound for the ith component is just ignored.
  299. ub : array_like, shape (n, )
  300. Upper bounds to each one of the components of ``x``.
  301. It is expected that ``ub >= 0``, otherwise the algorithm
  302. may fail. If ``ub[i] = Inf``, the upper bound for the ith
  303. component is just ignored.
  304. Returns
  305. -------
  306. x : array_like, shape (n,)
  307. Solution to the problem.
  308. Notes
  309. -----
  310. Based on implementations described in pp. 885-886 from [1]_.
  311. References
  312. ----------
  313. .. [1] Byrd, Richard H., Mary E. Hribar, and Jorge Nocedal.
  314. "An interior point algorithm for large-scale nonlinear
  315. programming." SIAM Journal on Optimization 9.4 (1999): 877-900.
  316. """
  317. # Compute minimum norm minimizer of 1/2*|| A x + b ||^2.
  318. newton_point = -Y.dot(b)
  319. # Check for interior point
  320. if inside_box_boundaries(newton_point, lb, ub) \
  321. and norm(newton_point) <= trust_radius:
  322. x = newton_point
  323. return x
  324. # Compute gradient vector ``g = A.T b``
  325. g = A.T.dot(b)
  326. # Compute Cauchy point
  327. # `cauchy_point = g.T g / (g.T A.T A g)``.
  328. A_g = A.dot(g)
  329. cauchy_point = -np.dot(g, g) / np.dot(A_g, A_g) * g
  330. # Origin
  331. origin_point = np.zeros_like(cauchy_point)
  332. # Check the segment between cauchy_point and newton_point
  333. # for a possible solution.
  334. z = cauchy_point
  335. p = newton_point - cauchy_point
  336. _, alpha, intersect = box_sphere_intersections(z, p, lb, ub,
  337. trust_radius)
  338. if intersect:
  339. x1 = z + alpha*p
  340. else:
  341. # Check the segment between the origin and cauchy_point
  342. # for a possible solution.
  343. z = origin_point
  344. p = cauchy_point
  345. _, alpha, _ = box_sphere_intersections(z, p, lb, ub,
  346. trust_radius)
  347. x1 = z + alpha*p
  348. # Check the segment between origin and newton_point
  349. # for a possible solution.
  350. z = origin_point
  351. p = newton_point
  352. _, alpha, _ = box_sphere_intersections(z, p, lb, ub,
  353. trust_radius)
  354. x2 = z + alpha*p
  355. # Return the best solution among x1 and x2.
  356. if norm(A.dot(x1) + b) < norm(A.dot(x2) + b):
  357. return x1
  358. else:
  359. return x2
  360. def projected_cg(H, c, Z, Y, b, trust_radius=np.inf,
  361. lb=None, ub=None, tol=None,
  362. max_iter=None, max_infeasible_iter=None,
  363. return_all=False):
  364. """Solve EQP problem with projected CG method.
  365. Solve equality-constrained quadratic programming problem
  366. ``min 1/2 x.T H x + x.t c`` subject to ``A x + b = 0`` and,
  367. possibly, to trust region constraints ``||x|| < trust_radius``
  368. and box constraints ``lb <= x <= ub``.
  369. Parameters
  370. ----------
  371. H : LinearOperator (or sparse matrix or ndarray), shape (n, n)
  372. Operator for computing ``H v``.
  373. c : array_like, shape (n,)
  374. Gradient of the quadratic objective function.
  375. Z : LinearOperator (or sparse matrix or ndarray), shape (n, n)
  376. Operator for projecting ``x`` into the null space of A.
  377. Y : LinearOperator, sparse matrix, ndarray, shape (n, m)
  378. Operator that, for a given a vector ``b``, compute smallest
  379. norm solution of ``A x + b = 0``.
  380. b : array_like, shape (m,)
  381. Right-hand side of the constraint equation.
  382. trust_radius : float, optional
  383. Trust radius to be considered. By default, uses ``trust_radius=inf``,
  384. which means no trust radius at all.
  385. lb : array_like, shape (n,), optional
  386. Lower bounds to each one of the components of ``x``.
  387. If ``lb[i] = -Inf`` the lower bound for the i-th
  388. component is just ignored (default).
  389. ub : array_like, shape (n, ), optional
  390. Upper bounds to each one of the components of ``x``.
  391. If ``ub[i] = Inf`` the upper bound for the i-th
  392. component is just ignored (default).
  393. tol : float, optional
  394. Tolerance used to interrupt the algorithm.
  395. max_iter : int, optional
  396. Maximum algorithm iterations. Where ``max_inter <= n-m``.
  397. By default, uses ``max_iter = n-m``.
  398. max_infeasible_iter : int, optional
  399. Maximum infeasible (regarding box constraints) iterations the
  400. algorithm is allowed to take.
  401. By default, uses ``max_infeasible_iter = n-m``.
  402. return_all : bool, optional
  403. When ``true``, return the list of all vectors through the iterations.
  404. Returns
  405. -------
  406. x : array_like, shape (n,)
  407. Solution of the EQP problem.
  408. info : Dict
  409. Dictionary containing the following:
  410. - niter : Number of iterations.
  411. - stop_cond : Reason for algorithm termination:
  412. 1. Iteration limit was reached;
  413. 2. Reached the trust-region boundary;
  414. 3. Negative curvature detected;
  415. 4. Tolerance was satisfied.
  416. - allvecs : List containing all intermediary vectors (optional).
  417. - hits_boundary : True if the proposed step is on the boundary
  418. of the trust region.
  419. Notes
  420. -----
  421. Implementation of Algorithm 6.2 on [1]_.
  422. In the absence of spherical and box constraints, for sufficient
  423. iterations, the method returns a truly optimal result.
  424. In the presence of those constraints, the value returned is only
  425. a inexpensive approximation of the optimal value.
  426. References
  427. ----------
  428. .. [1] Gould, Nicholas IM, Mary E. Hribar, and Jorge Nocedal.
  429. "On the solution of equality constrained quadratic
  430. programming problems arising in optimization."
  431. SIAM Journal on Scientific Computing 23.4 (2001): 1376-1395.
  432. """
  433. CLOSE_TO_ZERO = 1e-25
  434. n, = np.shape(c) # Number of parameters
  435. m, = np.shape(b) # Number of constraints
  436. # Initial Values
  437. x = Y.dot(-b)
  438. r = Z.dot(H.dot(x) + c)
  439. g = Z.dot(r)
  440. p = -g
  441. # Store ``x`` value
  442. if return_all:
  443. allvecs = [x]
  444. # Values for the first iteration
  445. H_p = H.dot(p)
  446. rt_g = norm(g)**2 # g.T g = r.T Z g = r.T g (ref [1]_ p.1389)
  447. # If x > trust-region the problem does not have a solution.
  448. tr_distance = trust_radius - norm(x)
  449. if tr_distance < 0:
  450. raise ValueError("Trust region problem does not have a solution.")
  451. # If x == trust_radius, then x is the solution
  452. # to the optimization problem, since x is the
  453. # minimum norm solution to Ax=b.
  454. elif tr_distance < CLOSE_TO_ZERO:
  455. info = {'niter': 0, 'stop_cond': 2, 'hits_boundary': True}
  456. if return_all:
  457. allvecs.append(x)
  458. info['allvecs'] = allvecs
  459. return x, info
  460. # Set default tolerance
  461. if tol is None:
  462. tol = max(min(0.01 * np.sqrt(rt_g), 0.1 * rt_g), CLOSE_TO_ZERO)
  463. # Set default lower and upper bounds
  464. if lb is None:
  465. lb = np.full(n, -np.inf)
  466. if ub is None:
  467. ub = np.full(n, np.inf)
  468. # Set maximum iterations
  469. if max_iter is None:
  470. max_iter = n-m
  471. max_iter = min(max_iter, n-m)
  472. # Set maximum infeasible iterations
  473. if max_infeasible_iter is None:
  474. max_infeasible_iter = n-m
  475. hits_boundary = False
  476. stop_cond = 1
  477. counter = 0
  478. last_feasible_x = np.zeros_like(x)
  479. k = 0
  480. for i in range(max_iter):
  481. # Stop criteria - Tolerance : r.T g < tol
  482. if rt_g < tol:
  483. stop_cond = 4
  484. break
  485. k += 1
  486. # Compute curvature
  487. pt_H_p = H_p.dot(p)
  488. # Stop criteria - Negative curvature
  489. if pt_H_p <= 0:
  490. if np.isinf(trust_radius):
  491. raise ValueError("Negative curvature not allowed "
  492. "for unrestricted problems.")
  493. else:
  494. # Find intersection with constraints
  495. _, alpha, intersect = box_sphere_intersections(
  496. x, p, lb, ub, trust_radius, entire_line=True)
  497. # Update solution
  498. if intersect:
  499. x = x + alpha*p
  500. # Reinforce variables are inside box constraints.
  501. # This is only necessary because of roundoff errors.
  502. x = reinforce_box_boundaries(x, lb, ub)
  503. # Attribute information
  504. stop_cond = 3
  505. hits_boundary = True
  506. break
  507. # Get next step
  508. alpha = rt_g / pt_H_p
  509. x_next = x + alpha*p
  510. # Stop criteria - Hits boundary
  511. if np.linalg.norm(x_next) >= trust_radius:
  512. # Find intersection with box constraints
  513. _, theta, intersect = box_sphere_intersections(x, alpha*p, lb, ub,
  514. trust_radius)
  515. # Update solution
  516. if intersect:
  517. x = x + theta*alpha*p
  518. # Reinforce variables are inside box constraints.
  519. # This is only necessary because of roundoff errors.
  520. x = reinforce_box_boundaries(x, lb, ub)
  521. # Attribute information
  522. stop_cond = 2
  523. hits_boundary = True
  524. break
  525. # Check if ``x`` is inside the box and start counter if it is not.
  526. if inside_box_boundaries(x_next, lb, ub):
  527. counter = 0
  528. else:
  529. counter += 1
  530. # Whenever outside box constraints keep looking for intersections.
  531. if counter > 0:
  532. _, theta, intersect = box_sphere_intersections(x, alpha*p, lb, ub,
  533. trust_radius)
  534. if intersect:
  535. last_feasible_x = x + theta*alpha*p
  536. # Reinforce variables are inside box constraints.
  537. # This is only necessary because of roundoff errors.
  538. last_feasible_x = reinforce_box_boundaries(last_feasible_x,
  539. lb, ub)
  540. counter = 0
  541. # Stop after too many infeasible (regarding box constraints) iteration.
  542. if counter > max_infeasible_iter:
  543. break
  544. # Store ``x_next`` value
  545. if return_all:
  546. allvecs.append(x_next)
  547. # Update residual
  548. r_next = r + alpha*H_p
  549. # Project residual g+ = Z r+
  550. g_next = Z.dot(r_next)
  551. # Compute conjugate direction step d
  552. rt_g_next = norm(g_next)**2 # g.T g = r.T g (ref [1]_ p.1389)
  553. beta = rt_g_next / rt_g
  554. p = - g_next + beta*p
  555. # Prepare for next iteration
  556. x = x_next
  557. g = g_next
  558. r = g_next
  559. rt_g = norm(g)**2 # g.T g = r.T Z g = r.T g (ref [1]_ p.1389)
  560. H_p = H.dot(p)
  561. if not inside_box_boundaries(x, lb, ub):
  562. x = last_feasible_x
  563. hits_boundary = True
  564. info = {'niter': k, 'stop_cond': stop_cond,
  565. 'hits_boundary': hits_boundary}
  566. if return_all:
  567. info['allvecs'] = allvecs
  568. return x, info