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, by adding code to handle the optional ``der`` input.

def f(x, sx, der=0):
    # When using gradient-based solvers, define extra if-cases for
    # handling der=1 (calculate df/dx) and der=2 (caldculate df/dsx).

L-BFGS-B Variations (gradient-based)

Implementations of the SurrogateOptimizer class.

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:
  • LBFGSB – Limited-memory bound-constrained BFGS (L-BFGS-B) method

  • TR_LBFGSB – L-BFGS-B is applied within a trust region

class optimizers.lbfgsb.LBFGSB(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 LBFGSB 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 evaluation budget per solve (default: 1000).

    • opt_restarts (int): Number of multisolve 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.TR_LBFGSB(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 TR_LBFGSB 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 evaluation budget per solve (default: 1000).

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

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