upoqa package
UPOQA Solver
UPOQA
UPOQA is a derivative-free model-based optimization method for partially separable problems, whose name is the abbreviation for Unconstrained Partially-separable Optimization by Quadratic Approximation.
- upoqa.minimize(fun, x0, coords={}, maxiter=None, maxfev={}, weights={}, xforms={}, xform_bounds={}, extra_fun=None, npt=None, radius_init=1.0, radius_final=1e-06, noise_level=0, seek_global_minimum=False, f_target=None, tr_shape='structured', callback=None, disp=True, verbose=False, debug=False, return_internals=False, options={}, **kwargs)[source]
Minimize an unconstrained objective function that has a partially-separable structure using a derivative-free model-based trust-region method.
Partial separability means the objective can be expressed as:
\[\min_{x\in\mathbb{R}^n} \quad \sum_{i=1}^q f_i(U_i x),\]where \(f_i:\mathbb{R}^{|\mathcal{I}_i|} \to \mathbb{R}\) are black-box functions (termed as element function or element) whose gradients and Hessians are unavailable, \(U_i:\mathbb{R}^n \to \mathbb{R}^{|\mathcal{I}_i|}\) are projection operators, and \(\mathcal{I}_i \subset [i]\) is an index set satisfying \(|\mathcal{I}_i| < n\) (element functions depend on small subsets of variables).
The solver also supports a more general objective form:
\[\min_{x\in\mathbb{R}^n} \quad f_0(x) + \sum_{i=1}^q w_i h_i\left(f_i(U_i x)\right),\]where:
\(f_0: \mathbb{R}^n \to \mathbb{R}\) is a white-box component with computable gradients and Hessians
\(w_i \in \mathbb{R}\) are element weights
\(h_i: \mathbb{R}\to\mathbb{R}\) are smooth transformations on the outputs of \(f_i\), with computable gradients/Hessians
- Parameters:
fun (list of callables, dict of callables, or callable) –
The objective element function(s). Can be:
A dictionary mapping element names to element functions
A list containing all element functions
A single function (not recommended; treated as a single element without exploiting structure)
When using dictionaries, element names can be any hashable type. Dictionary input is recommended for better readability, especially when elements have different contexts or lack inherent ordering.
For list input, the solver assigns default names as
"fun[0]","fun[1]", etc.x0 (ndarray, shape (n,)) – Initial guess for the optimization variables.
coords (dict of array_like or list of array-like, optional) –
Mapping of element names to variable indices. Each value specifies the variables an element depends on.
- Example for dictionary input:
{'f1': [0, 1], 'f2': [2, 3]}indicates:Element f1 depends on variables 0 and 1
Element f2 depends on variables 2 and 3
- Example for list input:
[[0, 1], [2, 3]]assigns indices to elements in list order.
If not provided, all elements default to depending on all variables.
maxiter (int, optional) – Maximum number of iterations. Default:
1000 * n.maxfev (int, dict of int or list of int, optional) –
Maximum function evaluations:
dict: Per-element limits (keys match
fun)list: Per-element limits (order matches
fun)int: Uniform limit for all elements
weights (dict or list, optional) –
Nonnegative weights \(w_i\) for elements. Formats:
dict: Per-element weights (keys match
fun)list: Per-element weights (order matches
fun)scalar: Uniform weight for all elements
xforms (dict or list, optional) –
Transformations \(h_i\) for element outputs. Each entry must be a list of length 3 containing [function, gradient, Hessian] callables.
Example for two simple transformations:
xforms = { 'ele1': [lambda x: x**2, lambda x: 2*x, lambda x: 2], 'ele2': [ lambda x: np.sqrt(x), lambda x: 1/(2*np.sqrt(x)), lambda x: -1/(4*x**(3/2)) ] }
Note: Appropriate use of transformations can accelerate convergence speed but may introduce instability. For simple cases, consider incorporating them directly into element functions.
xform_bounds (dict of tuple, list of tuple, optional) –
Bounds for the input domain of transformation functions.
Each entry must be a tuple of length 2 (lower_bound, upper_bound) defining the valid input range for a transformation. When provided, the solver automatically performs bound checking on transformation inputs, ensuring they remain within specified limits.
Example:
{'log_transform': (1e-8, np.inf)}prevents log(0) errors.Note: Even if element function outputs naturally stay within bounds, surrogate model outputs may exceed valid ranges during optimization. If your transformation functions cannot gracefully handle out-of-bound inputs, it is strongly recommended to provide this parameter to prevent algorithm crashes due to invalid transformation evaluations.
extra_fun (list of callables, optional) – White-box component \(f_0\) of the objective function, which is a list of length 3 specified as
[function_value, gradient, Hessian]callables.npt (int or dict or list, optional) – Number of interpolation points per element. Formats similar to
maxfev.radius_init (float, default=1.0) – Initial trust region radius
radius_final (float, default=1e-6) – Final (minimum) trust region radius
noise_level ({0, 1, 2}, default=0) –
Estimated noise level in objective function:
noise_level=0: No noisenoise_level=1: Moderate noise. Restart mechanisms enabled and the default interpolation points are set to\[\max\{\min\{0.8 n_i ^{1.5} + n_i + 1, (n_i + 1)(n_i + 2) / 2\}, 2 n_i + 1\}\]where \(n_i\) are the dimension of the \(i\)-th element.
noise_level=2: Significant noise. Restart mechanisms enabled and the default interpolation points are set to\[(n_i + 1)(n_i + 2) / 2\]for better noise resilience.
seek_global_minimum (bool, default=False) – Whether to enable global optimization restarts. When True, activates restart mechanisms with parameters optimized for global search (similar to Py-BOBYQA [1]).
f_target (float, optional) – Target on the objective function value. The optimization procedure is terminated when the objective function value at the iterate is less than or equal to this target.
tr_shape ({'structured', 'spherical'}, default='structured') –
Geometry of the trust region:
'structured': Complex region formed by intersecting cylinders with element-specific radii. In simple cases (e.g., two orthogonal cylinders), resembles a Steinmetz solid. Allows larger steps in well-modeled directions while restricting poorly-modeled ones.'spherical': Traditional spherical trust region. Simpler but less adaptive to element-wise model accuracy. Recommended if the'structured'option incurs excessive runtime overhead (typically when the number of elements is very large).
callback (callable, optional) –
A callback executed at each objective function evaluation. The method terminates if a
StopIterationexception is raised by the callback function. Its signature can be one of the following:callback(intermediate_result)where
intermediate_resultis a keyword parameter that contains an instance ofOptimizeResult, with attributesx,funand other useful information, being the point at which the objective function is evaluated and the value of the objective function, respectively. The name of the parameter must beintermediate_resultfor the callback to be passed an instance ofOptimizeResult.Alternatively, the callback function can have the signature:
callback(xk)where
xkis the point at which the objective function is evaluated. Introspection is used to determine which of the signatures to invoke.disp (bool, default=True) – Whether to display progress.
verbose (bool, default=False) – Whether to display detailed progress information.
debug (bool, default=False) –
Enable debugging features:
debug=1: NaN checking in function evaluationsdebug=2: Include solver internals (manager, interp_set, model) in resultdebug=3: Enhanced verbosity when both disp and verbose are True
return_internals (bool, default=False) – If True, the returned result will contain the fields
manager,interp_setandmodelfor debugging purposes even ifdebug=False.options (dict) – Advanced algorithm parameters. See the source code of
UPOQAParameterListfor available options.
- Returns:
Result of the optimization procedure, with the following fields:
- messagestr
Description of the cause of the termination.
- successbool
Whether the optimization procedure terminated successfully.
- funfloat
Value of objective function at the solution.
- funsdict or list
Values of the element functions evaluated at the solution. The type matches the input
fun:If
funwas a dictionary, return a dictionary mapping element names to their function values.If
funwas a list, return a list of function values in the same order as the input list.
- extra_funfloat
Value of extra function (the white-box component \(f_0\)) at the solution
- xndarray
The solution of the optimization.
- jac, hessndarray
Values of objective function’s Jacobian and its Hessian at
x. The Hessian is an approximation.- nitint
Number of iterations performed by the optimizer.
- nfevdict or list
Number of evaluations of element functions.
The type matches the input
fun:If
funwas a dictionary, return a dictionary mapping element names to their numbers of function evaluations.If
funwas a list, return a list of evaluation numbers in the same order as the input list.
- max_nfevint
Maximum number of evaluations of element functions.
- avg_nfevfloat
Average number of evaluations of element functions.
- nrunint
Number of runs (when enabling restarts).
If
debug=Trueorreturn_internals=True, the result also has the following fields:- manager
UPOQAManager Algorithm manager.
- interp_set
OverallInterpSet Interpolation point set.
- model
OverallSurrogate Surrogate model.
If an exception is raised during the optimization that cannot be handled by the algorithm, the result will also contain the following fields:
- exceptionException
Exception raised during the optimization that caused the termination, if any.
- tracebackstr
Traceback of the exception, if any.
- Return type:
References
[1]C. Cartis, J. Fiala, B. Marteau, and L. Roberts. Improving the flexibility and robustness of model-based derivative-free optimization solvers. ACM Trans. Math. Softw. 45, 3 (2019), 1-41. doi:10.1145/3338517.
Examples
Here is a simple example of using the solver to minimize a function with two elements. Consider the following problem:
\[\min_{x,y,z\in \mathbb{R}} \quad x^2 + 2y^2 + z^2 + 2xy - (y + 1)z\]let’s replace \([x,y,z]\) with \(\mathbf{x} = [x_1,x_2,x_3]\), and rewrite this problem into
\[\min_{\mathbf{x}\in\mathbb{R}^3} \quad f_1(x_1, x_2) + f_2(x_2, x_3),\]where
\[\begin{split}\begin{align*} &f_1(x_1, x_2) = x_1^2 + x_2^2 + 2x_1 x_2, \\ &f_2(x_2, x_3) = x_2^2 + x_3^2 - (x_2 + 1)x_3, \end{align*}\end{split}\]then the objective function can be implemented as:
>>> def f1(x): # f1(x,y) ... return x[0] ** 2 + x[1] ** 2 + 2 * x[0] * x[1] # x^2 + y^2 + 2xy >>> def f2(x): # f2(y,z) ... return x[0] ** 2 + x[1] ** 2 - (x[0] + 1) * x[1] # y^2 + z^2 - (y+1)z >>> fun = {'xy': f1, 'yz': f2 } >>> coords = {'xy': [0, 1], 'yz': [1, 2]}
The problem can be solved by:
>>> x0 = [0, 0, 0] >>> res = minimize(fun, x0, coords = coords, disp = False) >>> res.x array([-0.33333332, 0.3333333 , 0.66666665]) >>> res.fun np.float64(-0.33333333333333215)
UPOQA Utils
- upoqa.utils package
- upoqa.utils.interp_set module
- Interpolation Set
InterpSetInterpSet.nInterpSet.nptInterpSet.interp_set_YInterpSet.interp_set_fvalInterpSet.x_opt_idxInterpSet.x_anchor_idxInterpSet.x_optInterpSet.x_anchorInterpSet.f_optInterpSet.f_anchorInterpSet.f_anchorInterpSet.f_optInterpSet.get_anchor()InterpSet.get_interp_set()InterpSet.get_opt()InterpSet.init_interp_set()InterpSet.init_interp_set_Y()InterpSet.set_interp_set_fval()InterpSet.update_point_on_idx()InterpSet.update_x_opt()InterpSet.x_anchorInterpSet.x_opt
OverallInterpSetOverallInterpSet.nOverallInterpSet.ele_interp_setsOverallInterpSet.proj_onto_eleOverallInterpSet.ele_namesOverallInterpSet.ele_idxsOverallInterpSet.nptsOverallInterpSet.interp_set_YOverallInterpSet.interp_set_fvalOverallInterpSet.interp_set_fval_elesOverallInterpSet.interp_set_extra_fvalOverallInterpSet.x_opt_idxOverallInterpSet.enableOverallInterpSet.ele_numOverallInterpSet.set_sizeOverallInterpSet.disableOverallInterpSet.capacityOverallInterpSet.unallocated_idxOverallInterpSet.allocated_idxOverallInterpSet.x_optOverallInterpSet.f_optOverallInterpSet.allocated_idxOverallInterpSet.append()OverallInterpSet.capacityOverallInterpSet.clear()OverallInterpSet.clear_invalid_point()OverallInterpSet.clear_with_only_one_left()OverallInterpSet.delete_points()OverallInterpSet.disableOverallInterpSet.ele_numOverallInterpSet.f_optOverallInterpSet.get_ele_fvals()OverallInterpSet.get_interp_set()OverallInterpSet.get_opt()OverallInterpSet.set_sizeOverallInterpSet.unallocated_idxOverallInterpSet.update_fval_with_new_weights_and_xforms()OverallInterpSet.update_point_on_idx()OverallInterpSet.update_x_opt()OverallInterpSet.x_opt
- upoqa.utils.manager module
- Algorithm Manager
ExitInfoExitStatusMaxEvalNumReachedUPOQAManagerUPOQAManager.nUPOQAManager.nptsUPOQAManager.resolutionUPOQAManager.radiiUPOQAManager.ele_namesUPOQAManager.ele_numUPOQAManager.proj_onto_eleUPOQAManager.interp_setUPOQAManager.modelUPOQAManager.ele_modelsUPOQAManager.xformsUPOQAManager.weightsUPOQAManager.extra_funUPOQAManager.paramsUPOQAManager.maxiterUPOQAManager.maxfevUPOQAManager.nfevUPOQAManager.nrunUPOQAManager.dispUPOQAManager.x_optUPOQAManager.f_optUPOQAManager.average_radiusUPOQAManager.build_fval()UPOQAManager.calc_detailed_decrease()UPOQAManager.check_slow_iteration()UPOQAManager.couple_with_utils()UPOQAManager.decide_to_restart_due_to_exploding_coeff()UPOQAManager.ele_interp_setsUPOQAManager.ele_modelsUPOQAManager.extra_funUPOQAManager.f_optUPOQAManager.find_best_start_point()UPOQAManager.get_geometry_step()UPOQAManager.get_index_to_remove()UPOQAManager.get_reduction_ratio()UPOQAManager.get_update_requests()UPOQAManager.has_reached_maxfev()UPOQAManager.init_coord_overlap_relation()UPOQAManager.max_radiusUPOQAManager.nUPOQAManager.nptsUPOQAManager.prepare_for_a_restart()UPOQAManager.reduce_resolution()UPOQAManager.reset_restart_track_info()UPOQAManager.set_radius()UPOQAManager.set_resolution()UPOQAManager.set_weights()UPOQAManager.set_xforms()UPOQAManager.shift_center_to()UPOQAManager.soft_restart()UPOQAManager.update_radius()UPOQAManager.update_restart_track_info()UPOQAManager.update_weights()UPOQAManager.update_xforms()UPOQAManager.weightsUPOQAManager.wrap_xforms_with_bounds()UPOQAManager.x_optUPOQAManager.xforms
- upoqa.utils.model module
- Interpolation Surrogate Model
OverallSurrogateOverallSurrogate.nOverallSurrogate.interp_setOverallSurrogate.proj_onto_eleOverallSurrogate.coordsOverallSurrogate.ele_modelsOverallSurrogate.xformsOverallSurrogate.weightsOverallSurrogate.extra_funOverallSurrogate.paramsOverallSurrogate.ele_namesOverallSurrogate.model_centerOverallSurrogate.fun_eval()OverallSurrogate.grad_eval()OverallSurrogate.hess_eval()OverallSurrogate.hess_operator()OverallSurrogate.set_weights()OverallSurrogate.set_xforms()OverallSurrogate.shift_center_to()OverallSurrogate.update()OverallSurrogate.update_weights()OverallSurrogate.update_xforms()
QuadSurrogateQuadSurrogate.nQuadSurrogate.interp_setQuadSurrogate.nptQuadSurrogate.model_centerQuadSurrogate.model_consQuadSurrogate.model_gradQuadSurrogate.model_hess_explicitQuadSurrogate.model_hess_implicitQuadSurrogate.alt_grad_eval()QuadSurrogate.fun_eval()QuadSurrogate.get_determinant_ratio()QuadSurrogate.grad_eval()QuadSurrogate.hess_eval()QuadSurrogate.hess_operator()QuadSurrogate.inv_kkt_matrix_dot()QuadSurrogate.inv_kkt_matrix_partial_dot()QuadSurrogate.model_hessQuadSurrogate.reinit()QuadSurrogate.reset()QuadSurrogate.reset_surrogate_coeff()QuadSurrogate.shift_center_to()QuadSurrogate.solve_interpolation_KKT()QuadSurrogate.update()
SurrogateLinAlgError
- upoqa.utils.params module
- upoqa.utils.projection module
- upoqa.utils.result module
- Optimization Result
OptimizeResultOptimizeResult.messageOptimizeResult.successOptimizeResult.exceptionOptimizeResult.tracebackOptimizeResult.funOptimizeResult.funsOptimizeResult.extra_funOptimizeResult.xOptimizeResult.nitOptimizeResult.nfevOptimizeResult.max_nfevOptimizeResult.avg_nfevOptimizeResult.nrunOptimizeResult.managerOptimizeResult.interp_setOptimizeResult.model
- upoqa.utils.sub_solver module
- upoqa.utils.interp_set module