pdfo.bobyqa#
- pdfo.bobyqa(fun, x0, args=(), bounds=None, options=None)[source]#
Bounded Optimization BY Quadratic Approximations.
Deprecated since version 1.3: Calling the BOBYQA solver via the
bobyqa
function is deprecated. The BOBYQA solver remains available in PDFO. Call thepdfo
function with the argumentmethod='bobyqa'
to use it.- Parameters:
- fun: callable
Objective function to be minimized.
fun(x, *args) -> float
where
x
is an array with shape (n,) and args is a tuple.- x0: ndarray, shape (n,)
Initial guess.
- args: tuple, optional
Parameters of the objective function. For example,
pdfo(fun, x0, args, ...)
is equivalent to
pdfo(lambda x: fun(x, *args), x0, ...)
- bounds: {Bounds, ndarray, shape (n, 2)}, optional
Bound constraints of the problem. It can be one of the cases below.
An instance of
Bounds
.An ndarray with shape (n, 2). The bound constraint for x[i] is
bounds[i, 0] <= x[i] <= bounds[i, 1]
. Setbounds[i, 0]
to \(-\infty\) orNone
if there is no lower bound, and setbounds[i, 1]
to \(\infty\) orNone
if there is no upper bound.
- options: dict, optional
The options passed to the solver. It contains optionally:
- rhobeg: float, optional
Initial value of the trust region radius, which should be a positive scalar. Typically,
options['rhobeg']
should be in the order of one tenth of the greatest expected change to a variable. By default, it ismin(1, min(ub - lb) / 4)
, and0.5
if the problem is scaled.- rhoend: float, optional
Final value of the trust region radius, which should be a positive scalar.
options['rhoend']
should indicate the accuracy required in the final values of the variables. Moreover,options['rhoend']
should be no more thanoptions['rhobeg']
and is by default1e-6
.- maxfev: int, optional
Upper bound of the number of calls of the objective function fun. Its value must be not less than
options['npt'] + 1
. By default, it is500 * n
.- npt: int, optional
Number of interpolation points of each model used in Powell’s Fortran code.
- ftarget: float, optional
Target value of the objective function. If a feasible iterate achieves an objective function value lower or equal to
`options['ftarget']
, the algorithm stops immediately. By default, it is \(-\infty\).- scale: bool, optional
Whether to scale the problem according to the bound constraints. By default, it is
False
. If the problem is to be scaled, thenrhobeg
andrhoend
will be used as the initial and final trust region radii for the scaled problem.- honour_x0: bool, optional
Whether to respect the user-defined
x0
. By default, it isFalse
.- quiet: bool, optional
Whether the interface is quiet. If it is set to
True
, the output message will not be printed. This flag does not interfere with the warning and error printing.- classical: bool, optional
Whether to call the classical Powell code or not. It is not encouraged in production. By default, it is
False
.- debug: bool, optional
Debugging flag. It is not encouraged in production. By default, it is
False
.- chkfunval: bool, optional
Flag used when debugging. If both
options['debug']
andoptions['chkfunval']
areTrue
, an extra function/constraint evaluation would be performed to check whether the returned values of objective function and constraint match the returnedx
. By default, it isFalse
.
- Returns:
- res: OptimizeResult
The results of the solver. Check
OptimizeResult
for a description of the attributes.
See also
References
[1]M. J. D. Powell. The BOBYQA algorithm for bound constrained optimization without derivatives. Technical Report DAMTP 2009/NA06, Department of Applied Mathematics and Theoretical Physics, University of Cambridge, 2009.
Examples
The following example shows how to solve a simple bound-constrained optimization problem. The problem considered below should be solved with a derivative-based method. It is used here only as an illustration.
We consider the 2-dimensional problem
\[\begin{split}\min_{x, y \in \R} \quad x^2 + y^2 \quad \text{s.t.} \quad \left\{ \begin{array}{l} 0 \le x \le 2,\\ 1 / 2 \le y \le 3. \end{array} \right.\end{split}\]We solve this problem using
bobyqa
starting from the initial guess \((x_0, y_0) = (0, 1)\) with at most 200 function evaluations.>>> from pdfo import Bounds, bobyqa >>> bounds = Bounds([0, 0.5], [2, 3]) >>> options = {'maxfev': 200} >>> res = bobyqa(lambda x: x[0]**2 + x[1]**2, [0, 1], bounds=bounds, options=options) >>> res.x array([0. , 0.5])
Note that
bobyqa
can also be used to solve unconstrained problems.