framework.py 38 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240
  1. import warnings
  2. import numpy as np
  3. from scipy.optimize import lsq_linear
  4. from .models import Models, Quadratic
  5. from .settings import Options, Constants
  6. from .subsolvers import (
  7. cauchy_geometry,
  8. spider_geometry,
  9. normal_byrd_omojokun,
  10. tangential_byrd_omojokun,
  11. constrained_tangential_byrd_omojokun,
  12. )
  13. from .subsolvers.optim import qr_tangential_byrd_omojokun
  14. from .utils import get_arrays_tol
  15. TINY = np.finfo(float).tiny
  16. EPS = np.finfo(float).eps
  17. class TrustRegion:
  18. """
  19. Trust-region framework.
  20. """
  21. def __init__(self, pb, options, constants):
  22. """
  23. Initialize the trust-region framework.
  24. Parameters
  25. ----------
  26. pb : `cobyqa.problem.Problem`
  27. Problem to solve.
  28. options : dict
  29. Options of the solver.
  30. constants : dict
  31. Constants of the solver.
  32. Raises
  33. ------
  34. `cobyqa.utils.MaxEvalError`
  35. If the maximum number of evaluations is reached.
  36. `cobyqa.utils.TargetSuccess`
  37. If a nearly feasible point has been found with an objective
  38. function value below the target.
  39. `cobyqa.utils.FeasibleSuccess`
  40. If a feasible point has been found for a feasibility problem.
  41. `numpy.linalg.LinAlgError`
  42. If the initial interpolation system is ill-defined.
  43. """
  44. # Set the initial penalty parameter.
  45. self._penalty = 0.0
  46. # Initialize the models.
  47. self._pb = pb
  48. self._models = Models(self._pb, options, self.penalty)
  49. self._constants = constants
  50. # Set the index of the best interpolation point.
  51. self._best_index = 0
  52. self.set_best_index()
  53. # Set the initial Lagrange multipliers.
  54. self._lm_linear_ub = np.zeros(self.m_linear_ub)
  55. self._lm_linear_eq = np.zeros(self.m_linear_eq)
  56. self._lm_nonlinear_ub = np.zeros(self.m_nonlinear_ub)
  57. self._lm_nonlinear_eq = np.zeros(self.m_nonlinear_eq)
  58. self.set_multipliers(self.x_best)
  59. # Set the initial trust-region radius and the resolution.
  60. self._resolution = options[Options.RHOBEG]
  61. self._radius = self.resolution
  62. @property
  63. def n(self):
  64. """
  65. Number of variables.
  66. Returns
  67. -------
  68. int
  69. Number of variables.
  70. """
  71. return self._pb.n
  72. @property
  73. def m_linear_ub(self):
  74. """
  75. Number of linear inequality constraints.
  76. Returns
  77. -------
  78. int
  79. Number of linear inequality constraints.
  80. """
  81. return self._pb.m_linear_ub
  82. @property
  83. def m_linear_eq(self):
  84. """
  85. Number of linear equality constraints.
  86. Returns
  87. -------
  88. int
  89. Number of linear equality constraints.
  90. """
  91. return self._pb.m_linear_eq
  92. @property
  93. def m_nonlinear_ub(self):
  94. """
  95. Number of nonlinear inequality constraints.
  96. Returns
  97. -------
  98. int
  99. Number of nonlinear inequality constraints.
  100. """
  101. return self._pb.m_nonlinear_ub
  102. @property
  103. def m_nonlinear_eq(self):
  104. """
  105. Number of nonlinear equality constraints.
  106. Returns
  107. -------
  108. int
  109. Number of nonlinear equality constraints.
  110. """
  111. return self._pb.m_nonlinear_eq
  112. @property
  113. def radius(self):
  114. """
  115. Trust-region radius.
  116. Returns
  117. -------
  118. float
  119. Trust-region radius.
  120. """
  121. return self._radius
  122. @radius.setter
  123. def radius(self, radius):
  124. """
  125. Set the trust-region radius.
  126. Parameters
  127. ----------
  128. radius : float
  129. New trust-region radius.
  130. """
  131. self._radius = radius
  132. if (
  133. self.radius
  134. <= self._constants[Constants.DECREASE_RADIUS_THRESHOLD]
  135. * self.resolution
  136. ):
  137. self._radius = self.resolution
  138. @property
  139. def resolution(self):
  140. """
  141. Resolution of the trust-region framework.
  142. The resolution is a lower bound on the trust-region radius.
  143. Returns
  144. -------
  145. float
  146. Resolution of the trust-region framework.
  147. """
  148. return self._resolution
  149. @resolution.setter
  150. def resolution(self, resolution):
  151. """
  152. Set the resolution of the trust-region framework.
  153. Parameters
  154. ----------
  155. resolution : float
  156. New resolution of the trust-region framework.
  157. """
  158. self._resolution = resolution
  159. @property
  160. def penalty(self):
  161. """
  162. Penalty parameter.
  163. Returns
  164. -------
  165. float
  166. Penalty parameter.
  167. """
  168. return self._penalty
  169. @property
  170. def models(self):
  171. """
  172. Models of the objective function and constraints.
  173. Returns
  174. -------
  175. `cobyqa.models.Models`
  176. Models of the objective function and constraints.
  177. """
  178. return self._models
  179. @property
  180. def best_index(self):
  181. """
  182. Index of the best interpolation point.
  183. Returns
  184. -------
  185. int
  186. Index of the best interpolation point.
  187. """
  188. return self._best_index
  189. @property
  190. def x_best(self):
  191. """
  192. Best interpolation point.
  193. Its value is interpreted as relative to the origin, not the base point.
  194. Returns
  195. -------
  196. `numpy.ndarray`
  197. Best interpolation point.
  198. """
  199. return self.models.interpolation.point(self.best_index)
  200. @property
  201. def fun_best(self):
  202. """
  203. Value of the objective function at `x_best`.
  204. Returns
  205. -------
  206. float
  207. Value of the objective function at `x_best`.
  208. """
  209. return self.models.fun_val[self.best_index]
  210. @property
  211. def cub_best(self):
  212. """
  213. Values of the nonlinear inequality constraints at `x_best`.
  214. Returns
  215. -------
  216. `numpy.ndarray`, shape (m_nonlinear_ub,)
  217. Values of the nonlinear inequality constraints at `x_best`.
  218. """
  219. return self.models.cub_val[self.best_index, :]
  220. @property
  221. def ceq_best(self):
  222. """
  223. Values of the nonlinear equality constraints at `x_best`.
  224. Returns
  225. -------
  226. `numpy.ndarray`, shape (m_nonlinear_eq,)
  227. Values of the nonlinear equality constraints at `x_best`.
  228. """
  229. return self.models.ceq_val[self.best_index, :]
  230. def lag_model(self, x):
  231. """
  232. Evaluate the Lagrangian model at a given point.
  233. Parameters
  234. ----------
  235. x : `numpy.ndarray`, shape (n,)
  236. Point at which the Lagrangian model is evaluated.
  237. Returns
  238. -------
  239. float
  240. Value of the Lagrangian model at `x`.
  241. """
  242. return (
  243. self.models.fun(x)
  244. + self._lm_linear_ub
  245. @ (self._pb.linear.a_ub @ x - self._pb.linear.b_ub)
  246. + self._lm_linear_eq
  247. @ (self._pb.linear.a_eq @ x - self._pb.linear.b_eq)
  248. + self._lm_nonlinear_ub @ self.models.cub(x)
  249. + self._lm_nonlinear_eq @ self.models.ceq(x)
  250. )
  251. def lag_model_grad(self, x):
  252. """
  253. Evaluate the gradient of the Lagrangian model at a given point.
  254. Parameters
  255. ----------
  256. x : `numpy.ndarray`, shape (n,)
  257. Point at which the gradient of the Lagrangian model is evaluated.
  258. Returns
  259. -------
  260. `numpy.ndarray`, shape (n,)
  261. Gradient of the Lagrangian model at `x`.
  262. """
  263. return (
  264. self.models.fun_grad(x)
  265. + self._lm_linear_ub @ self._pb.linear.a_ub
  266. + self._lm_linear_eq @ self._pb.linear.a_eq
  267. + self._lm_nonlinear_ub @ self.models.cub_grad(x)
  268. + self._lm_nonlinear_eq @ self.models.ceq_grad(x)
  269. )
  270. def lag_model_hess(self):
  271. """
  272. Evaluate the Hessian matrix of the Lagrangian model at a given point.
  273. Returns
  274. -------
  275. `numpy.ndarray`, shape (n, n)
  276. Hessian matrix of the Lagrangian model at `x`.
  277. """
  278. hess = self.models.fun_hess()
  279. if self.m_nonlinear_ub > 0:
  280. hess += self._lm_nonlinear_ub @ self.models.cub_hess()
  281. if self.m_nonlinear_eq > 0:
  282. hess += self._lm_nonlinear_eq @ self.models.ceq_hess()
  283. return hess
  284. def lag_model_hess_prod(self, v):
  285. """
  286. Evaluate the right product of the Hessian matrix of the Lagrangian
  287. model with a given vector.
  288. Parameters
  289. ----------
  290. v : `numpy.ndarray`, shape (n,)
  291. Vector with which the Hessian matrix of the Lagrangian model is
  292. multiplied from the right.
  293. Returns
  294. -------
  295. `numpy.ndarray`, shape (n,)
  296. Right product of the Hessian matrix of the Lagrangian model with
  297. `v`.
  298. """
  299. return (
  300. self.models.fun_hess_prod(v)
  301. + self._lm_nonlinear_ub @ self.models.cub_hess_prod(v)
  302. + self._lm_nonlinear_eq @ self.models.ceq_hess_prod(v)
  303. )
  304. def lag_model_curv(self, v):
  305. """
  306. Evaluate the curvature of the Lagrangian model along a given direction.
  307. Parameters
  308. ----------
  309. v : `numpy.ndarray`, shape (n,)
  310. Direction along which the curvature of the Lagrangian model is
  311. evaluated.
  312. Returns
  313. -------
  314. float
  315. Curvature of the Lagrangian model along `v`.
  316. """
  317. return (
  318. self.models.fun_curv(v)
  319. + self._lm_nonlinear_ub @ self.models.cub_curv(v)
  320. + self._lm_nonlinear_eq @ self.models.ceq_curv(v)
  321. )
  322. def sqp_fun(self, step):
  323. """
  324. Evaluate the objective function of the SQP subproblem.
  325. Parameters
  326. ----------
  327. step : `numpy.ndarray`, shape (n,)
  328. Step along which the objective function of the SQP subproblem is
  329. evaluated.
  330. Returns
  331. -------
  332. float
  333. Value of the objective function of the SQP subproblem along `step`.
  334. """
  335. return step @ (
  336. self.models.fun_grad(self.x_best)
  337. + 0.5 * self.lag_model_hess_prod(step)
  338. )
  339. def sqp_cub(self, step):
  340. """
  341. Evaluate the linearization of the nonlinear inequality constraints.
  342. Parameters
  343. ----------
  344. step : `numpy.ndarray`, shape (n,)
  345. Step along which the linearization of the nonlinear inequality
  346. constraints is evaluated.
  347. Returns
  348. -------
  349. `numpy.ndarray`, shape (m_nonlinear_ub,)
  350. Value of the linearization of the nonlinear inequality constraints
  351. along `step`.
  352. """
  353. return (
  354. self.models.cub(self.x_best)
  355. + self.models.cub_grad(self.x_best) @ step
  356. )
  357. def sqp_ceq(self, step):
  358. """
  359. Evaluate the linearization of the nonlinear equality constraints.
  360. Parameters
  361. ----------
  362. step : `numpy.ndarray`, shape (n,)
  363. Step along which the linearization of the nonlinear equality
  364. constraints is evaluated.
  365. Returns
  366. -------
  367. `numpy.ndarray`, shape (m_nonlinear_ub,)
  368. Value of the linearization of the nonlinear equality constraints
  369. along `step`.
  370. """
  371. return (
  372. self.models.ceq(self.x_best)
  373. + self.models.ceq_grad(self.x_best) @ step
  374. )
  375. def merit(self, x, fun_val=None, cub_val=None, ceq_val=None):
  376. """
  377. Evaluate the merit function at a given point.
  378. Parameters
  379. ----------
  380. x : `numpy.ndarray`, shape (n,)
  381. Point at which the merit function is evaluated.
  382. fun_val : float, optional
  383. Value of the objective function at `x`. If not provided, the
  384. objective function is evaluated at `x`.
  385. cub_val : `numpy.ndarray`, shape (m_nonlinear_ub,), optional
  386. Values of the nonlinear inequality constraints. If not provided,
  387. the nonlinear inequality constraints are evaluated at `x`.
  388. ceq_val : `numpy.ndarray`, shape (m_nonlinear_eq,), optional
  389. Values of the nonlinear equality constraints. If not provided,
  390. the nonlinear equality constraints are evaluated at `x`.
  391. Returns
  392. -------
  393. float
  394. Value of the merit function at `x`.
  395. """
  396. if fun_val is None or cub_val is None or ceq_val is None:
  397. fun_val, cub_val, ceq_val = self._pb(x, self.penalty)
  398. m_val = fun_val
  399. if self._penalty > 0.0:
  400. c_val = self._pb.violation(x, cub_val=cub_val, ceq_val=ceq_val)
  401. if np.count_nonzero(c_val):
  402. m_val += self._penalty * np.linalg.norm(c_val)
  403. return m_val
  404. def get_constraint_linearizations(self, x):
  405. """
  406. Get the linearizations of the constraints at a given point.
  407. Parameters
  408. ----------
  409. x : `numpy.ndarray`, shape (n,)
  410. Point at which the linearizations of the constraints are evaluated.
  411. Returns
  412. -------
  413. `numpy.ndarray`, shape (m_linear_ub + m_nonlinear_ub, n)
  414. Left-hand side matrix of the linearized inequality constraints.
  415. `numpy.ndarray`, shape (m_linear_ub + m_nonlinear_ub,)
  416. Right-hand side vector of the linearized inequality constraints.
  417. `numpy.ndarray`, shape (m_linear_eq + m_nonlinear_eq, n)
  418. Left-hand side matrix of the linearized equality constraints.
  419. `numpy.ndarray`, shape (m_linear_eq + m_nonlinear_eq,)
  420. Right-hand side vector of the linearized equality constraints.
  421. """
  422. aub = np.block(
  423. [
  424. [self._pb.linear.a_ub],
  425. [self.models.cub_grad(x)],
  426. ]
  427. )
  428. bub = np.block(
  429. [
  430. self._pb.linear.b_ub - self._pb.linear.a_ub @ x,
  431. -self.models.cub(x),
  432. ]
  433. )
  434. aeq = np.block(
  435. [
  436. [self._pb.linear.a_eq],
  437. [self.models.ceq_grad(x)],
  438. ]
  439. )
  440. beq = np.block(
  441. [
  442. self._pb.linear.b_eq - self._pb.linear.a_eq @ x,
  443. -self.models.ceq(x),
  444. ]
  445. )
  446. return aub, bub, aeq, beq
  447. def get_trust_region_step(self, options):
  448. """
  449. Get the trust-region step.
  450. The trust-region step is computed by solving the derivative-free
  451. trust-region SQP subproblem using a Byrd-Omojokun composite-step
  452. approach. For more details, see Section 5.2.3 of [1]_.
  453. Parameters
  454. ----------
  455. options : dict
  456. Options of the solver.
  457. Returns
  458. -------
  459. `numpy.ndarray`, shape (n,)
  460. Normal step.
  461. `numpy.ndarray`, shape (n,)
  462. Tangential step.
  463. References
  464. ----------
  465. .. [1] T. M. Ragonneau. *Model-Based Derivative-Free Optimization
  466. Methods and Software*. PhD thesis, Department of Applied
  467. Mathematics, The Hong Kong Polytechnic University, Hong Kong, China,
  468. 2022. URL: https://theses.lib.polyu.edu.hk/handle/200/12294.
  469. """
  470. # Evaluate the linearizations of the constraints.
  471. aub, bub, aeq, beq = self.get_constraint_linearizations(self.x_best)
  472. xl = self._pb.bounds.xl - self.x_best
  473. xu = self._pb.bounds.xu - self.x_best
  474. # Evaluate the normal step.
  475. radius = self._constants[Constants.BYRD_OMOJOKUN_FACTOR] * self.radius
  476. normal_step = normal_byrd_omojokun(
  477. aub,
  478. bub,
  479. aeq,
  480. beq,
  481. xl,
  482. xu,
  483. radius,
  484. options[Options.DEBUG],
  485. **self._constants,
  486. )
  487. if options[Options.DEBUG]:
  488. tol = get_arrays_tol(xl, xu)
  489. if (np.any(normal_step + tol < xl)
  490. or np.any(xu < normal_step - tol)):
  491. warnings.warn(
  492. "the normal step does not respect the bound constraint.",
  493. RuntimeWarning,
  494. 2,
  495. )
  496. if np.linalg.norm(normal_step) > 1.1 * radius:
  497. warnings.warn(
  498. "the normal step does not respect the trust-region "
  499. "constraint.",
  500. RuntimeWarning,
  501. 2,
  502. )
  503. # Evaluate the tangential step.
  504. radius = np.sqrt(self.radius**2.0 - normal_step @ normal_step)
  505. xl -= normal_step
  506. xu -= normal_step
  507. bub = np.maximum(bub - aub @ normal_step, 0.0)
  508. g_best = self.models.fun_grad(self.x_best) + self.lag_model_hess_prod(
  509. normal_step
  510. )
  511. if self._pb.type in ["unconstrained", "bound-constrained"]:
  512. tangential_step = tangential_byrd_omojokun(
  513. g_best,
  514. self.lag_model_hess_prod,
  515. xl,
  516. xu,
  517. radius,
  518. options[Options.DEBUG],
  519. **self._constants,
  520. )
  521. else:
  522. tangential_step = constrained_tangential_byrd_omojokun(
  523. g_best,
  524. self.lag_model_hess_prod,
  525. xl,
  526. xu,
  527. aub,
  528. bub,
  529. aeq,
  530. radius,
  531. options["debug"],
  532. **self._constants,
  533. )
  534. if options[Options.DEBUG]:
  535. tol = get_arrays_tol(xl, xu)
  536. if np.any(tangential_step + tol < xl) or np.any(
  537. xu < tangential_step - tol
  538. ):
  539. warnings.warn(
  540. "The tangential step does not respect the bound "
  541. "constraints.",
  542. RuntimeWarning,
  543. 2,
  544. )
  545. if (
  546. np.linalg.norm(normal_step + tangential_step)
  547. > 1.1 * np.sqrt(2.0) * self.radius
  548. ):
  549. warnings.warn(
  550. "The trial step does not respect the trust-region "
  551. "constraint.",
  552. RuntimeWarning,
  553. 2,
  554. )
  555. return normal_step, tangential_step
  556. def get_geometry_step(self, k_new, options):
  557. """
  558. Get the geometry-improving step.
  559. Three different geometry-improving steps are computed and the best one
  560. is returned. For more details, see Section 5.2.7 of [1]_.
  561. Parameters
  562. ----------
  563. k_new : int
  564. Index of the interpolation point to be modified.
  565. options : dict
  566. Options of the solver.
  567. Returns
  568. -------
  569. `numpy.ndarray`, shape (n,)
  570. Geometry-improving step.
  571. Raises
  572. ------
  573. `numpy.linalg.LinAlgError`
  574. If the computation of a determinant fails.
  575. References
  576. ----------
  577. .. [1] T. M. Ragonneau. *Model-Based Derivative-Free Optimization
  578. Methods and Software*. PhD thesis, Department of Applied
  579. Mathematics, The Hong Kong Polytechnic University, Hong Kong, China,
  580. 2022. URL: https://theses.lib.polyu.edu.hk/handle/200/12294.
  581. """
  582. if options[Options.DEBUG]:
  583. assert (
  584. k_new != self.best_index
  585. ), "The index `k_new` must be different from the best index."
  586. # Build the k_new-th Lagrange polynomial.
  587. coord_vec = np.squeeze(np.eye(1, self.models.npt, k_new))
  588. lag = Quadratic(
  589. self.models.interpolation,
  590. coord_vec,
  591. options[Options.DEBUG],
  592. )
  593. g_lag = lag.grad(self.x_best, self.models.interpolation)
  594. # Compute a simple constrained Cauchy step.
  595. xl = self._pb.bounds.xl - self.x_best
  596. xu = self._pb.bounds.xu - self.x_best
  597. step = cauchy_geometry(
  598. 0.0,
  599. g_lag,
  600. lambda v: lag.curv(v, self.models.interpolation),
  601. xl,
  602. xu,
  603. self.radius,
  604. options[Options.DEBUG],
  605. )
  606. sigma = self.models.determinants(self.x_best + step, k_new)
  607. # Compute the solution on the straight lines joining the interpolation
  608. # points to the k-th one, and choose it if it provides a larger value
  609. # of the determinant of the interpolation system in absolute value.
  610. xpt = (
  611. self.models.interpolation.xpt
  612. - self.models.interpolation.xpt[:, self.best_index, np.newaxis]
  613. )
  614. xpt[:, [0, self.best_index]] = xpt[:, [self.best_index, 0]]
  615. step_alt = spider_geometry(
  616. 0.0,
  617. g_lag,
  618. lambda v: lag.curv(v, self.models.interpolation),
  619. xpt[:, 1:],
  620. xl,
  621. xu,
  622. self.radius,
  623. options[Options.DEBUG],
  624. )
  625. sigma_alt = self.models.determinants(self.x_best + step_alt, k_new)
  626. if abs(sigma_alt) > abs(sigma):
  627. step = step_alt
  628. sigma = sigma_alt
  629. # Compute a Cauchy step on the tangent space of the active constraints.
  630. if self._pb.type in [
  631. "linearly constrained",
  632. "nonlinearly constrained",
  633. ]:
  634. aub, bub, aeq, beq = (
  635. self.get_constraint_linearizations(self.x_best))
  636. tol_bd = get_arrays_tol(xl, xu)
  637. tol_ub = get_arrays_tol(bub)
  638. free_xl = xl <= -tol_bd
  639. free_xu = xu >= tol_bd
  640. free_ub = bub >= tol_ub
  641. # Compute the Cauchy step.
  642. n_act, q = qr_tangential_byrd_omojokun(
  643. aub,
  644. aeq,
  645. free_xl,
  646. free_xu,
  647. free_ub,
  648. )
  649. g_lag_proj = q[:, n_act:] @ (q[:, n_act:].T @ g_lag)
  650. norm_g_lag_proj = np.linalg.norm(g_lag_proj)
  651. if 0 < n_act < self._pb.n and norm_g_lag_proj > TINY * self.radius:
  652. step_alt = (self.radius / norm_g_lag_proj) * g_lag_proj
  653. if lag.curv(step_alt, self.models.interpolation) < 0.0:
  654. step_alt = -step_alt
  655. # Evaluate the constraint violation at the Cauchy step.
  656. cbd = np.block([xl - step_alt, step_alt - xu])
  657. cub = aub @ step_alt - bub
  658. ceq = aeq @ step_alt - beq
  659. maxcv_val = max(
  660. np.max(array, initial=0.0)
  661. for array in [cbd, cub, np.abs(ceq)]
  662. )
  663. # Accept the new step if it is nearly feasible and do not
  664. # drastically worsen the determinant of the interpolation
  665. # system in absolute value.
  666. tol = np.max(np.abs(step_alt[~free_xl]), initial=0.0)
  667. tol = np.max(np.abs(step_alt[~free_xu]), initial=tol)
  668. tol = np.max(np.abs(aub[~free_ub, :] @ step_alt), initial=tol)
  669. tol = min(10.0 * tol, 1e-2 * np.linalg.norm(step_alt))
  670. if maxcv_val <= tol:
  671. sigma_alt = self.models.determinants(
  672. self.x_best + step_alt, k_new
  673. )
  674. if abs(sigma_alt) >= 0.1 * abs(sigma):
  675. step = np.clip(step_alt, xl, xu)
  676. if options[Options.DEBUG]:
  677. tol = get_arrays_tol(xl, xu)
  678. if np.any(step + tol < xl) or np.any(xu < step - tol):
  679. warnings.warn(
  680. "The geometry step does not respect the bound "
  681. "constraints.",
  682. RuntimeWarning,
  683. 2,
  684. )
  685. if np.linalg.norm(step) > 1.1 * self.radius:
  686. warnings.warn(
  687. "The geometry step does not respect the "
  688. "trust-region constraint.",
  689. RuntimeWarning,
  690. 2,
  691. )
  692. return step
  693. def get_second_order_correction_step(self, step, options):
  694. """
  695. Get the second-order correction step.
  696. Parameters
  697. ----------
  698. step : `numpy.ndarray`, shape (n,)
  699. Trust-region step.
  700. options : dict
  701. Options of the solver.
  702. Returns
  703. -------
  704. `numpy.ndarray`, shape (n,)
  705. Second-order correction step.
  706. """
  707. # Evaluate the linearizations of the constraints.
  708. aub, bub, aeq, beq = self.get_constraint_linearizations(self.x_best)
  709. xl = self._pb.bounds.xl - self.x_best
  710. xu = self._pb.bounds.xu - self.x_best
  711. radius = np.linalg.norm(step)
  712. soc_step = normal_byrd_omojokun(
  713. aub,
  714. bub,
  715. aeq,
  716. beq,
  717. xl,
  718. xu,
  719. radius,
  720. options[Options.DEBUG],
  721. **self._constants,
  722. )
  723. if options[Options.DEBUG]:
  724. tol = get_arrays_tol(xl, xu)
  725. if np.any(soc_step + tol < xl) or np.any(xu < soc_step - tol):
  726. warnings.warn(
  727. "The second-order correction step does not "
  728. "respect the bound constraints.",
  729. RuntimeWarning,
  730. 2,
  731. )
  732. if np.linalg.norm(soc_step) > 1.1 * radius:
  733. warnings.warn(
  734. "The second-order correction step does not "
  735. "respect the trust-region constraint.",
  736. RuntimeWarning,
  737. 2,
  738. )
  739. return soc_step
  740. def get_reduction_ratio(self, step, fun_val, cub_val, ceq_val):
  741. """
  742. Get the reduction ratio.
  743. Parameters
  744. ----------
  745. step : `numpy.ndarray`, shape (n,)
  746. Trust-region step.
  747. fun_val : float
  748. Objective function value at the trial point.
  749. cub_val : `numpy.ndarray`, shape (m_nonlinear_ub,)
  750. Nonlinear inequality constraint values at the trial point.
  751. ceq_val : `numpy.ndarray`, shape (m_nonlinear_eq,)
  752. Nonlinear equality constraint values at the trial point.
  753. Returns
  754. -------
  755. float
  756. Reduction ratio.
  757. """
  758. merit_old = self.merit(
  759. self.x_best,
  760. self.fun_best,
  761. self.cub_best,
  762. self.ceq_best,
  763. )
  764. merit_new = self.merit(self.x_best + step, fun_val, cub_val, ceq_val)
  765. merit_model_old = self.merit(
  766. self.x_best,
  767. 0.0,
  768. self.models.cub(self.x_best),
  769. self.models.ceq(self.x_best),
  770. )
  771. merit_model_new = self.merit(
  772. self.x_best + step,
  773. self.sqp_fun(step),
  774. self.sqp_cub(step),
  775. self.sqp_ceq(step),
  776. )
  777. if abs(merit_model_old - merit_model_new) > TINY * abs(
  778. merit_old - merit_new
  779. ):
  780. return (merit_old - merit_new) / abs(
  781. merit_model_old - merit_model_new
  782. )
  783. else:
  784. return -1.0
  785. def increase_penalty(self, step):
  786. """
  787. Increase the penalty parameter.
  788. Parameters
  789. ----------
  790. step : `numpy.ndarray`, shape (n,)
  791. Trust-region step.
  792. """
  793. aub, bub, aeq, beq = self.get_constraint_linearizations(self.x_best)
  794. viol_diff = max(
  795. np.linalg.norm(
  796. np.block(
  797. [
  798. np.maximum(0.0, -bub),
  799. beq,
  800. ]
  801. )
  802. )
  803. - np.linalg.norm(
  804. np.block(
  805. [
  806. np.maximum(0.0, aub @ step - bub),
  807. aeq @ step - beq,
  808. ]
  809. )
  810. ),
  811. 0.0,
  812. )
  813. sqp_val = self.sqp_fun(step)
  814. threshold = np.linalg.norm(
  815. np.block(
  816. [
  817. self._lm_linear_ub,
  818. self._lm_linear_eq,
  819. self._lm_nonlinear_ub,
  820. self._lm_nonlinear_eq,
  821. ]
  822. )
  823. )
  824. if abs(viol_diff) > TINY * abs(sqp_val):
  825. threshold = max(threshold, sqp_val / viol_diff)
  826. best_index_save = self.best_index
  827. if (
  828. self._penalty
  829. <= self._constants[Constants.PENALTY_INCREASE_THRESHOLD]
  830. * threshold
  831. ):
  832. self._penalty = max(
  833. self._constants[Constants.PENALTY_INCREASE_FACTOR] * threshold,
  834. 1.0,
  835. )
  836. self.set_best_index()
  837. return best_index_save == self.best_index
  838. def decrease_penalty(self):
  839. """
  840. Decrease the penalty parameter.
  841. """
  842. self._penalty = min(self._penalty, self._get_low_penalty())
  843. self.set_best_index()
  844. def set_best_index(self):
  845. """
  846. Set the index of the best point.
  847. """
  848. best_index = self.best_index
  849. m_best = self.merit(
  850. self.x_best,
  851. self.models.fun_val[best_index],
  852. self.models.cub_val[best_index, :],
  853. self.models.ceq_val[best_index, :],
  854. )
  855. r_best = self._pb.maxcv(
  856. self.x_best,
  857. self.models.cub_val[best_index, :],
  858. self.models.ceq_val[best_index, :],
  859. )
  860. tol = (
  861. 10.0
  862. * EPS
  863. * max(self.models.n, self.models.npt)
  864. * max(abs(m_best), 1.0)
  865. )
  866. for k in range(self.models.npt):
  867. if k != self.best_index:
  868. x_val = self.models.interpolation.point(k)
  869. m_val = self.merit(
  870. x_val,
  871. self.models.fun_val[k],
  872. self.models.cub_val[k, :],
  873. self.models.ceq_val[k, :],
  874. )
  875. r_val = self._pb.maxcv(
  876. x_val,
  877. self.models.cub_val[k, :],
  878. self.models.ceq_val[k, :],
  879. )
  880. if m_val < m_best or (m_val < m_best + tol and r_val < r_best):
  881. best_index = k
  882. m_best = m_val
  883. r_best = r_val
  884. self._best_index = best_index
  885. def get_index_to_remove(self, x_new=None):
  886. """
  887. Get the index of the interpolation point to remove.
  888. If `x_new` is not provided, the index returned should be used during
  889. the geometry-improvement phase. Otherwise, the index returned is the
  890. best index for included `x_new` in the interpolation set.
  891. Parameters
  892. ----------
  893. x_new : `numpy.ndarray`, shape (n,), optional
  894. New point to be included in the interpolation set.
  895. Returns
  896. -------
  897. int
  898. Index of the interpolation point to remove.
  899. float
  900. Distance between `x_best` and the removed point.
  901. Raises
  902. ------
  903. `numpy.linalg.LinAlgError`
  904. If the computation of a determinant fails.
  905. """
  906. dist_sq = np.sum(
  907. (
  908. self.models.interpolation.xpt
  909. - self.models.interpolation.xpt[:, self.best_index, np.newaxis]
  910. )
  911. ** 2.0,
  912. axis=0,
  913. )
  914. if x_new is None:
  915. sigma = 1.0
  916. weights = dist_sq
  917. else:
  918. sigma = self.models.determinants(x_new)
  919. weights = (
  920. np.maximum(
  921. 1.0,
  922. dist_sq
  923. / max(
  924. self._constants[Constants.LOW_RADIUS_FACTOR]
  925. * self.radius,
  926. self.resolution,
  927. )
  928. ** 2.0,
  929. )
  930. ** 3.0
  931. )
  932. weights[self.best_index] = -1.0 # do not remove the best point
  933. k_max = np.argmax(weights * np.abs(sigma))
  934. return k_max, np.sqrt(dist_sq[k_max])
  935. def update_radius(self, step, ratio):
  936. """
  937. Update the trust-region radius.
  938. Parameters
  939. ----------
  940. step : `numpy.ndarray`, shape (n,)
  941. Trust-region step.
  942. ratio : float
  943. Reduction ratio.
  944. """
  945. s_norm = np.linalg.norm(step)
  946. if ratio <= self._constants[Constants.LOW_RATIO]:
  947. self.radius *= self._constants[Constants.DECREASE_RADIUS_FACTOR]
  948. elif ratio <= self._constants[Constants.HIGH_RATIO]:
  949. self.radius = max(
  950. self._constants[Constants.DECREASE_RADIUS_FACTOR]
  951. * self.radius,
  952. s_norm,
  953. )
  954. else:
  955. self.radius = min(
  956. self._constants[Constants.INCREASE_RADIUS_FACTOR]
  957. * self.radius,
  958. max(
  959. self._constants[Constants.DECREASE_RADIUS_FACTOR]
  960. * self.radius,
  961. self._constants[Constants.INCREASE_RADIUS_THRESHOLD]
  962. * s_norm,
  963. ),
  964. )
  965. def enhance_resolution(self, options):
  966. """
  967. Enhance the resolution of the trust-region framework.
  968. Parameters
  969. ----------
  970. options : dict
  971. Options of the solver.
  972. """
  973. if (
  974. self._constants[Constants.LARGE_RESOLUTION_THRESHOLD]
  975. * options[Options.RHOEND]
  976. < self.resolution
  977. ):
  978. self.resolution *= self._constants[
  979. Constants.DECREASE_RESOLUTION_FACTOR
  980. ]
  981. elif (
  982. self._constants[Constants.MODERATE_RESOLUTION_THRESHOLD]
  983. * options[Options.RHOEND]
  984. < self.resolution
  985. ):
  986. self.resolution = np.sqrt(self.resolution
  987. * options[Options.RHOEND])
  988. else:
  989. self.resolution = options[Options.RHOEND]
  990. # Reduce the trust-region radius.
  991. self._radius = max(
  992. self._constants[Constants.DECREASE_RADIUS_FACTOR] * self._radius,
  993. self.resolution,
  994. )
  995. def shift_x_base(self, options):
  996. """
  997. Shift the base point to `x_best`.
  998. Parameters
  999. ----------
  1000. options : dict
  1001. Options of the solver.
  1002. """
  1003. self.models.shift_x_base(np.copy(self.x_best), options)
  1004. def set_multipliers(self, x):
  1005. """
  1006. Set the Lagrange multipliers.
  1007. This method computes and set the Lagrange multipliers of the linear and
  1008. nonlinear constraints to be the QP multipliers.
  1009. Parameters
  1010. ----------
  1011. x : `numpy.ndarray`, shape (n,)
  1012. Point at which the Lagrange multipliers are computed.
  1013. """
  1014. # Build the constraints of the least-squares problem.
  1015. incl_linear_ub = self._pb.linear.a_ub @ x >= self._pb.linear.b_ub
  1016. incl_nonlinear_ub = self.cub_best >= 0.0
  1017. incl_xl = self._pb.bounds.xl >= x
  1018. incl_xu = self._pb.bounds.xu <= x
  1019. m_linear_ub = np.count_nonzero(incl_linear_ub)
  1020. m_nonlinear_ub = np.count_nonzero(incl_nonlinear_ub)
  1021. m_xl = np.count_nonzero(incl_xl)
  1022. m_xu = np.count_nonzero(incl_xu)
  1023. if (
  1024. m_linear_ub + m_nonlinear_ub + self.m_linear_eq
  1025. + self.m_nonlinear_eq > 0
  1026. ):
  1027. identity = np.eye(self._pb.n)
  1028. c_jac = np.r_[
  1029. -identity[incl_xl, :],
  1030. identity[incl_xu, :],
  1031. self._pb.linear.a_ub[incl_linear_ub, :],
  1032. self.models.cub_grad(x, incl_nonlinear_ub),
  1033. self._pb.linear.a_eq,
  1034. self.models.ceq_grad(x),
  1035. ]
  1036. # Solve the least-squares problem.
  1037. g_best = self.models.fun_grad(x)
  1038. xl_lm = np.full(c_jac.shape[0], -np.inf)
  1039. xl_lm[: m_xl + m_xu + m_linear_ub + m_nonlinear_ub] = 0.0
  1040. res = lsq_linear(
  1041. c_jac.T,
  1042. -g_best,
  1043. bounds=(xl_lm, np.inf),
  1044. method="bvls",
  1045. )
  1046. # Extract the Lagrange multipliers.
  1047. self._lm_linear_ub[incl_linear_ub] = res.x[
  1048. m_xl + m_xu:m_xl + m_xu + m_linear_ub
  1049. ]
  1050. self._lm_linear_ub[~incl_linear_ub] = 0.0
  1051. self._lm_nonlinear_ub[incl_nonlinear_ub] = res.x[
  1052. m_xl
  1053. + m_xu
  1054. + m_linear_ub:m_xl
  1055. + m_xu
  1056. + m_linear_ub
  1057. + m_nonlinear_ub
  1058. ]
  1059. self._lm_nonlinear_ub[~incl_nonlinear_ub] = 0.0
  1060. self._lm_linear_eq[:] = res.x[
  1061. m_xl
  1062. + m_xu
  1063. + m_linear_ub
  1064. + m_nonlinear_ub:m_xl
  1065. + m_xu
  1066. + m_linear_ub
  1067. + m_nonlinear_ub
  1068. + self.m_linear_eq
  1069. ]
  1070. self._lm_nonlinear_eq[:] = res.x[
  1071. m_xl + m_xu + m_linear_ub + m_nonlinear_ub + self.m_linear_eq:
  1072. ]
  1073. def _get_low_penalty(self):
  1074. r_val_ub = np.c_[
  1075. (
  1076. self.models.interpolation.x_base[np.newaxis, :]
  1077. + self.models.interpolation.xpt.T
  1078. )
  1079. @ self._pb.linear.a_ub.T
  1080. - self._pb.linear.b_ub[np.newaxis, :],
  1081. self.models.cub_val,
  1082. ]
  1083. r_val_eq = (
  1084. self.models.interpolation.x_base[np.newaxis, :]
  1085. + self.models.interpolation.xpt.T
  1086. ) @ self._pb.linear.a_eq.T - self._pb.linear.b_eq[np.newaxis, :]
  1087. r_val_eq = np.block(
  1088. [
  1089. r_val_eq,
  1090. -r_val_eq,
  1091. self.models.ceq_val,
  1092. -self.models.ceq_val,
  1093. ]
  1094. )
  1095. r_val = np.block([r_val_ub, r_val_eq])
  1096. c_min = np.nanmin(r_val, axis=0)
  1097. c_max = np.nanmax(r_val, axis=0)
  1098. indices = (
  1099. c_min
  1100. < self._constants[Constants.THRESHOLD_RATIO_CONSTRAINTS] * c_max
  1101. )
  1102. if np.any(indices):
  1103. f_min = np.nanmin(self.models.fun_val)
  1104. f_max = np.nanmax(self.models.fun_val)
  1105. c_min_neg = np.minimum(0.0, c_min[indices])
  1106. c_diff = np.min(c_max[indices] - c_min_neg)
  1107. if c_diff > TINY * (f_max - f_min):
  1108. penalty = (f_max - f_min) / c_diff
  1109. else:
  1110. penalty = np.inf
  1111. else:
  1112. penalty = 0.0
  1113. return penalty