Surrogate Optimizers

When initializing a new MOOP object (see MOOP Classes), you must provide a surrogate optimization problem solver, which will be used to generate candidate solutions for each iteration.

from parmoo import optimizers

Note that when using a gradient-based technique, you must provide gradient evaluation options for all objective and constraint functions.

To implement a custom surrogate optimizer, import and extend the SurrogateOptimizer ABC.

from parmoo.optimizers.surrogate_optimizer import SurrogateOptimizer

The SurrogateOptimizer ABC and library of all available implementations in ParMOO are documented below.

SurrogateOptimizer

The abstract base class (ABC) for SurrogateOptimizer objects.

class optimizers.surrogate_optimizer.SurrogateOptimizer(o, lb, ub, hyperparams)

ABC describing surrogate optimization techniques.

This class contains the following methods.
  • setObjective(obj_func) (default implementation provided)

  • setSimulation(sim_func, sd_func) (default implementation provided)

  • setConstraints(constraint_func) (default implementation provided)

  • setPenalty(penaltyFunc, gradFunc) (default implementation provided)

  • setTrFunc(trFunc) (default implementation provided)

  • addAcquisition(*args) (default implementation provided)

  • returnResults(x, fx, sx, sdx)

  • solve(x)

  • save(filename)

  • load(filename)

abstract __init__(o, lb, ub, hyperparams)

Constructor for the SurrogateOptimizer class.

Parameters:
  • o (int) – The number of objectives.

  • lb (ndarray) – A 1D array of lower bounds for the design space.

  • ub (ndarray) – A 1D array of upper bounds for the design space.

  • hyperparams (dict) – A dictionary of hyperparameters for the optimization procedure.

Returns:

A new SurrogateOptimizer object.

Return type:

SurrogateOptimizer

setObjective(obj_func)

Add a vector-valued objective function that will be solved.

Parameters:

obj_func (function) – A vector-valued function that can be evaluated to solve the surrogate optimization problem.

setSimulation(sim_func, sd_func)

Add a vector-valued simulation function, used to calculate objs.

Parameters:
  • sim_func (function) – A vector-valued function that can be evaluated to determine the surrogate-predicted simulation outputs.

  • sd_func (function) – A vector-valued function that can be evaluated to determine the standard deviations of the surrogate predictions.

setPenalty(penalty_func)

Add a matrix-valued gradient function for obj_func.

Parameters:
  • penalty_func (function) – A vector-valued penalized objective that incorporates a penalty for violating constraints.

  • grad_func (function) – A matrix-valued function that can be evaluated to obtain the Jacobian matrix for obj_func.

setConstraints(constraint_func)

Add a constraint function that will be satisfied.

Parameters:

constraint_func (function) – A vector-valued function from the design space whose components correspond to constraint violations. If the problem has only bound constraints, this function returns zeros.

setTrFunc(trFunc)

Add a TR setter function for alerting surrogates.

Parameters:

trFunc (function) – A function with 2 inputs, which the optimizer must call prior to solving each surrogate optimization problem in order to set the trust-region center and radius.

returnResults(x, fx, sx, sdx)

This is a callback function to collect evaluation results.

Implement this function to receive the results of each true simulation evaluation from the MOOP class at runtime.

Parameters:
  • x (ndarray) – A 1D array with the design point evaluated.

  • fx (ndarray) – A 1D array with the objective function values at x.

  • sx (ndarray) – The simulation function values at x.

  • sdx (ndarray) – The standard deviation in the simulation prediction.

addAcquisition(*args)

Add an acquisition function for the surrogate optimizer.

Parameters:

*args (AcquisitionFunction) – Acquisition functions that are used to scalarize the list of objectives in order to solve the surrogate optimization problem.

abstract solve(x_k)

Solve the surrogate problem.

You may assume that the following internal attributes are defined and contain callable definitions of the objective, constraint, penalty, and simulation (surrogate) functions, respectively:

  • self.objectives,

  • self.constraints,

  • self.penalty_func, and

  • self.simulations.

Additionally, you may assume that:
  • self.acquisitions contains a list of one or more AcquisitionFunction object instances, each of whose acq.scalarize(f_vals, x_vals, s_vals_mean, s_vals_sd) is set and ready to call; and

  • self.setTR(x, r) can be called to set a trust-region centered at x with radius r (and re-fit the surrogates accordingly).

Note: If implementing your own solver, try to jit (or re-jit) any of the objective, constraint, penalty, simulation surrogate, and/or acquisition functions after each call to self.setTR. Additionally, if provided by the user, the objectives, constraints, penalty, and acq.scalarize, functions should all be differentiable by importing and calling jax.jacrev().

Parameters:

x_k (ndarray) – A 2D array containing a list of current iterates.

Returns:

A 2D array matching the shape of x_k specifying x_{k+1}.

Return type:

ndarray

save(filename)

Save important data from this class so that it can be reloaded.

Note: If this function is left unimplemented, ParMOO will reinitialize a fresh instance after a save/load. If this is the desired behavior, then this method and the load method need not be implemented.

Parameters:

filename (string) – The relative or absolute path to the file where all reload data should be saved.

load(filename)

Reload important data into this class after a previous save.

Note: If this function is left unimplemented, ParMOO will reinitialize a fresh instance after a save/load. If this is the desired behavior, then this method and the save method need not be implemented.

Parameters:

filename (string) – The relative or absolute path to the file where all reload data has been saved.

Pattern Search Techniques (gradient-free)

This module contains implementations of the SurrogateOptimizer ABC, which are based on pattern search.

Note that these strategies are all gradient-free, and therefore do not require any gradients to be defined. This makes them friendly for getting started.

The classes include:
  • LocalSurrogate_PS – A multistart pattern search (PS) algorithm

  • GlobalSurrogate_PS – global random search, followed by LocalPS

class optimizers.pattern_search.GlobalSurrogate_PS(o, lb, ub, hyperparams)

Use randomized search globally followed by a local PS.

Use RandomSearch to globally search the design space followed LocalSurrogate_PS to refine the potentially efficient solutions.

__init__(o, lb, ub, hyperparams)

Constructor for the GlobalPS class.

Parameters:
  • o (int) – The number of objectives.

  • lb (numpy.ndarray) – A 1d array of lower bounds for the design region. The number of design variables is inferred from the dimension of lb.

  • ub (numpy.ndarray) – A 1d array of upper bounds for the design region. The dimension must match ub.

  • hyperparams (dict) –

    A dictionary of hyperparameters for the optimization procedure. It may contain the following:

    • opt_budget (int): The function evaluation budget (default: 1500)

    • ps_budget (int): The number of the total opt_budget evaluations that will be used by PS (default: 2/3 of opt_budget).

Returns:

A new SurrogateOptimizer object.

Return type:

SurrogateOptimizer

solve(x)

Solve the surrogate problem by using random search followed by PS.

Parameters:

x (np.ndarray) – A 2d array containing a list of feasible design points used to warm start the search.

Returns:

A 2d numpy.ndarray containing a list of potentially efficient design points that were found by the optimizers.

Return type:

np.ndarray

class optimizers.pattern_search.LocalSurrogate_PS(o, lb, ub, hyperparams)

Use multistart Pattern Search to solve surrogate problem locally.

Applies PS to the surrogate problem, in order to identify design points that are locally Pareto optimal, with respect to the surrogate problem. Sorts poll directions by most recently used and attempts to step in promising directions in late iterations.

__init__(o, lb, ub, hyperparams)

Constructor for the LocalSurrogate_PS class.

Parameters:
  • o (int) – The number of objectives.

  • lb (numpy.ndarray) – A 1d array of lower bounds for the design region. The number of design variables is inferred from the dimension of lb.

  • ub (numpy.ndarray) – A 1d array of upper bounds for the design region. The dimension must match ub.

  • hyperparams (dict) –

    A dictionary of hyperparameters for the optimization procedure. It may contain the following:

    • opt_budget (int): The PS iteration limit (default: 1000).

    • opt_restarts (int): Number of multistart restarts per scalarization (default: n+1).

Returns:

A new SurrogateOptimizer object.

Return type:

SurrogateOptimizer

returnResults(x, fx, sx, sdx)

Collect the results of a function evaluation.

Parameters:
  • x (np.ndarray) – The design point evaluated.

  • fx (np.ndarray) – The objective function values at x.

  • sx (np.ndarray) – The simulation function values at x.

  • sdx (np.ndarray) – The standard deviation in the simulation outputs at x.

solve(x)

Solve the surrogate problem in a trust region via pattern search.

Parameters:

x (np.ndarray) – A 2d array containing a list of feasible design points used to warm start the search.

Returns:

A 2d numpy.ndarray of potentially efficient design points that were found by the PS optimizer.

Return type:

np.ndarray

save(filename)

Save important data from this class so that it can be reloaded.

Parameters:

filename (string) – The relative or absolute path to the file where all reload data should be saved.

load(filename)

Reload important data into this class after a previous save.

Parameters:

filename (string) – The relative or absolute path to the file where all reload data has been saved.

L-BFGS-B Variations (gradient-based)

Optimization methods based on limited-memory BFGS-B (L-BFGS-B).

This module contains implementations of the SurrogateOptimizer ABC, which are based on the L-BFGS-B quasi-Newton algorithm.

Note that all of these methods are gradient based, and therefore require objective, constraint, and surrogate gradient methods to be defined.

The classes include:
  • GlobalSurrogate_BFGS – Minimize the surrogate globally via L-BFGS-B

  • LocalSurrogate_BFGS – Minimize surrogate in trust region via L-BFGS-B

class optimizers.lbfgsb.GlobalSurrogate_BFGS(o, lb, ub, hyperparams)

Use L-BFGS-B and gradients to identify local solutions.

Applies L-BFGS-B to the surrogate problem, in order to identify design points that are locally Pareto optimal with respect to the surrogate problem.

__init__(o, lb, ub, hyperparams)

Constructor for the GlobalSurrogate_BFGS class.

Parameters:
  • o (int) – The number of objectives.

  • lb (numpy.ndarray) – A 1d array of lower bounds for the design region. The number of design variables is inferred from the dimension of lb.

  • ub (numpy.ndarray) – A 1d array of upper bounds for the design region. The dimension must match ub.

  • hyperparams (dict) –

    A dictionary of hyperparameters for the optimization procedure. It may contain the following:

    • opt_budget (int): The iteration limit per solve (default: 100).

    • opt_restarts (int): Number of multistart restarts per scalarization (default: n+1).

Returns:

A new SurrogateOptimizer object.

Return type:

SurrogateOptimizer

solve(x)

Solve the surrogate problem using L-BFGS-B.

Parameters:

x (np.ndarray) – A 2d array containing a list of design points used to warm start the search.

Returns:

A 2d numpy.ndarray of potentially efficient design points that were found by L-BFGS-B.

Return type:

np.ndarray

class optimizers.lbfgsb.LocalSurrogate_BFGS(o, lb, ub, hyperparams)

Use L-BFGS-B and gradients to identify solutions within a trust region.

Applies L-BFGS-B to the surrogate problem, in order to identify design points that are locally Pareto optimal with respect to the surrogate problem.

__init__(o, lb, ub, hyperparams)

Constructor for the LocalSurrogate_BFGS class.

Parameters:
  • o (int) – The number of objectives.

  • lb (numpy.ndarray) – A 1d array of lower bounds for the design region. The number of design variables is inferred from the dimension of lb.

  • ub (numpy.ndarray) – A 1d array of upper bounds for the design region. The dimension must match ub.

  • hyperparams (dict) –

    A dictionary of hyperparameters for the optimization procedure. It may contain the following:

    • opt_budget (int): The iteration limit per solve (default: 500).

    • opt_restarts (int): Number of multistart restarts per scalarization (default: 2).

Returns:

A new SurrogateOptimizer object.

Return type:

SurrogateOptimizer

returnResults(x, fx, sx, sdx)

Collect the results of a function evaluation.

Parameters:
  • x (np.ndarray) – The design point evaluated.

  • fx (np.ndarray) – The objective function values at x.

  • sx (np.ndarray) – The simulation function values at x.

  • sdx (np.ndarray) – The standard deviation in the simulation outputs at x.

solve(x)

Solve the surrogate problem using L-BFGS-B.

Parameters:

x (np.ndarray) – A 2d array containing a list of design points used to warm start the search.

Returns:

A 2d numpy.ndarray of potentially efficient design points that were found by L-BFGS-B.

Return type:

np.ndarray

save(filename)

Save important data from this class so that it can be reloaded.

Parameters:

filename (string) – The relative or absolute path to the file where all reload data should be saved.

load(filename)

Reload important data into this class after a previous save.

Parameters:

filename (string) – The relative or absolute path to the file where all reload data has been saved.