Constrained Optimization

conOpt: Constraint optimization

This subpackage contains algorithms for solving constraint optimization problems. Equality and inequality constraints problems

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\begin{eqnarray} e(x) & = 0 \\ g(x) & \leq 0 \end{eqnarray}

can be solved with

  1. Penalty method
  2. Augmented Lagrangian method
  3. SQP with a BFGS update strategy (at the moment only equality constraints).

Note, the penalty and the augmented lagrangian method works with an inner optimization method which can be chosen from the methods from the unConopt subpackage.

For problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]
\[\text{s. t. }\quad x_a \leq x \leq x_b\]

one can use

  1. Projected gradient method
  2. The L-BFGS-B method
  3. Projected Newton-Krylov method (if you can provide the action of the second derivative).

The projected Newton-Krylov method uses either the oppy implementation of cg or gmres to solve the corresponding linear systems.

The following line search methods for the box constraint solvers are available:

  1. Armijo Goldstein
  2. Projected Newton

In a future release we will add solvers for a combination of equality, inequality and box constraints.

Documentation is available in the docstrings and online here.

Methods

oppy.conOpt.armijo_goldstein module

Armijo-Goldstein linesearch

This linesearch method can be used in optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]
\[\text{s. t. } l\leq x \leq u\]

for finding an optimal stepsize length \(t\) in every optimization iteration.

oppy.conOpt.armijo_goldstein.armijo_goldstein(fhandle, dfhandle, x, d, low, up, options=None)

Armijo-Goldstein stepsize algorithm is a linesearch method for box constrained optimization problems. The linesearch terminates if the sufficient decrease condition

\[f(x(t)) \leq f(x)-c_1 \langle \nabla f(x),(x-x(t)) \rangle\]

is violated. Here, \(\langle \cdot,\cdot\rangle\) denotes the inner product and \(x(t) = P(x + t\cdot d)\), where the projector \(P\) is given by:

\[\begin{split}P(y)=\begin{cases} y \quad &\mathrm{,if}\quad l \leq y \leq u\\ l \quad &\mathrm{,if}\quad l >y\\ u \quad &\mathrm{,if}\quad u <y\\ \end{cases}\end{split}\]

where \(l\) is the lower and \(u\) the upper bound.

Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Current iteration point. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • d (numpy.ndarray, shape (n,)) – New direction. Array of real elements of size (n,), where ‘n’ is the number of independent variables. low : ndarray, shape (n,) Lower bound for the projection. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • up (numpy.ndarray, shape (n,)) – Upper bound for the projection. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • low (numpy.ndarray, shape (n,)) – Lower bound for the projection. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • options (Options, optional) –

    Options for armijo_goldstein(), represented as a Options object. The default is None. Attributes containing

    • initstepfloat
      Initial stepsize, initialized to \(1\).
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • c1float
      Weighting parameter, initialized to \(10^{-2}\).
    • betafloat
      Parameter to reduce stepsize, initialized to \(1/2\).
    • nmaxint or float
      Maximum number of iterations, initialized with \(30\).
    • inner_prodcallable()
      Callable function that computes the inner product of the difference \(x - x(t)\). The default is None and therefore initialized to the standard \(\mathbb{R}^n\) scalarproduct, i.e.
      \[\mathrm{inner\_prod} : (x,x(t)) \mapsto\]
      \[(x-x(t))^T(x-x(t))\]

      Note

      The norm function in the optimizer has to be changed appropriately to the inner_prod function in the linesearch algorithm, i.e.

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}\]
Returns:

  • t (float) – Final steplength

  • flag (int) – The flag indicated, whether the linesearch was successful, possible outputs are:

    1: Final steplength complies with sufficient decrease.
    0: Algorithm did not terminate and ended with nmax.

Examples

Minimizing problem of the rosenbrock() function

>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock

The l_bfgs_b() method was chosen as the minimizing method

>>> from oppy.conOpt import l_bfgs_b
>>> from oppy.conOpt import armijo_goldstein
>>> import numpy as np

Defining the armijo_goldstein() linesearch with default options:

>>> def linesearch(x, d, low, up):
>>>     t, flag = armijo_goldstein(fhandle, dfhandle, x, d, low, up)
>>>     return t, flag

The Options for the l_bfgs_b() are set as default as well. Input data for the l_bfgs_b() method:

>>> x0 = np.array([0.5, 0.5])
>>> low = np.array([-1.5, -0.5])
>>> up = np.array([1.5, 2.5])

Execute l_bfgs_b() method with armijo_goldstein() linesearch:

>>> result = l_bfgs_b(rosenbrock, gradient_rosenbrock, x0, low, up)
>>> print(result.x)

The result is [1., 1.], which is the exact solution [1, 1].

oppy.conOpt.augmented_lagrangian module

Augmented Lagrangian method

This method solves optimization problems with equality and inequality constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\begin{eqnarray} e(x) & = 0 \\ g(x) & \leq 0 \end{eqnarray}

using augmented_lagrangian() method.

oppy.conOpt.augmented_lagrangian.augmented_lagrangian(fhandle, dfhandle, x, ehandle=None, dehandle=None, ghandle=None, dghandle=None, options=None)

Augmented Lagrangian.

This method solves optimization problems with equality and inequality constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\begin{eqnarray} e(x) & = 0 \\ g(x) & \leq 0. \end{eqnarray}
Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and fargs is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and fargs is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Initial guess. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • ehandle (callable(), optional) –

    The equality restriction function with \(\texttt{ehandle(x)=0}\).

    ehandle(x, *eargs) -> array_like, shape (me,)

    where x is an 1-D array with shape (n,) and eargs is a tuple of the fixed parameters needed to completely specify the function.

    Note

    If you pass ehandle, you must also pass the derivative dehandle. Also, you must pass either ehandle or ghandle.

  • dehandle (callable(), optional) –

    The gradient of the equality restriction function.

    dehandle(x, *eargs) -> array_like, shape (n,me).

    where x is an 1-D array with shape (n,) and eargs is a tuple of the fixed parameters needed to completely specify the function.

  • ghandle (callable(), optional) –

    The inequality restriction function with \(\texttt{ghandle(x)\leq 0}\).

    ghandle(x, *gargs) -> array_like, shape (mg,)

    where x is an 1-D array with shape (n,) and gargs is a tuple of the fixed parameters needed to completely specify the function.

    Note

    If you pass ghandle, you must also pass the derivative dghandle. Also, you must pass either ehandle or ghandle.

  • dghandle (callable(), optional) –

    The gradient of the inequality restriction function.

    dghandle(x, *gargs) -> array_like, shape (n,mg).

    where x is an 1-D array with shape (n,) and gargs is a tuple of the fixed parameters needed to completely specify the function.

  • options (Options, optional) –

    Options for augmented_lagrangian() method, represented as a Options object. The default is None. Possible options are:

    • tol_absfloat
      Absolute tolerance for termination condition, initialized with \(10^{-7}\).
    • tol_relfloat
      Relative tolerance for termination condition, initialized with \(10^{-7}\).

    Note

    The algorithm stops, if the norm of the current residual res[k] satisfies the inequality

    \[\mathrm{res}[k] < \mathrm{stopping\_criteria},\]

    where stopping_criteria is given by

    \[\mathrm{stopping\_criteria} = \mathrm{tol\_rel} \cdot \mathrm{res}[0] + \mathrm{tol\_abs},\]

    and res[0] is the norm of the first residual. If the value of stopping_criteria is bigger than \(0.1\), e.g. if the norm of the first residual is big and the relative tolerance is small, then this value is set to \(0.1\).

    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the functions \(e\) and \(g\) in a vector \(x\). The default is None and therefore initialized to the standard norm in \(\mathbb{R}^n\), i.e.
      \[\mathrm{norm}: x \mapsto ||x||_2.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • eargstuple
      Tuple containing additional arguments for ehandle.
    • gargstuple
      Tuple containing additional arguments for ghandle.
    • lamnumpy.ndarray, shape (me,)
      Initial lagrange mulitplicator for the equality constraints, initialized with zeros.
    • munumpy.ndarray, shape (mg,)
      Initial lagrange mulitplicator for the inequality constraints, initialized with zeros.
    • taufloat
      Intitial tolerance for the inner optimization method, intialized with \(10^{-3}\).
    • c0float
      Intitial penalty factor, initialized with \(0.1\).
    • betafloat
      Update parameter for the penalty factor, intialized with \(4\).
    • optimeth: callable()
      Callable function that computes a optimum, initialized with gradmethod().
    • resultsstr or list
      String or list to specify the attributes for the Results object, initialized with 'normal'.
Returns:

results – Instance of the class Results. Represents the optimization results. Attributes of the object depend on the choice of the options.results parameter. With default choice the object has the following attributes:

  • xnumpy.ndarray, shape (n,)

    Solution vector of the optimization problem.

  • iterationsint

    Number of iterations the algorithm needed.

  • reslist

    List with the norm of the residuals.

Return type:

Results

References

  1. Bertsekas. Nonlinear Programming. Second Edition, see Bertsekas [Ber99].

Examples

Minimizing problem of the rosenbrock() function subject to

\begin{eqnarray} (x_1-1)^3-x_2+1& \leq 0 \\ x_1+x_2-2 & \leq 0 \end{eqnarray}

Importing the augmented_lagrangian() method and defining the problem:

>>> import numpy as np
>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock
>>> from oppy.conOpt import augmented_lagrangian
>>> g = (lambda x: np.array([(x[0]-1)**3-x[1]+1, x[0]+x[1]-2]))
>>> dg = (lambda x: np.array([[3*(x[0]-1)**2, -1],
>>>                           [1, 1]]))
>>> f = rosenbrock
>>> df = gradient_rosenbrock
>>> x0 = np.array([0.5, -0.5])

Execute augmented_lagrangian() method:

>>> result = augmented_lagrangian(fhandle=f, dfhandle=df, x=x0, ghandle=g,
>>>                               dghandle=dg)
>>> print(result.x)

The result is [0.99967779 1.0001933 ], which is close to the exact solution [1, 1].

oppy.conOpt.l_bfgs_b module

L_BFGS_B method

This method solves optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\[l \leq x \leq u\]

using the l_bfgs_b() method.

oppy.conOpt.l_bfgs_b.bfgsrecb(nt, sstore, ystore, d, activeset)

Help method bfgsrecb for L_BFGS_B.

Recursively compute action of an approximated Hessian on a vector using stored information of the history of the iteration

Parameters:
  • nt (int) – Number of recursive steps or vectors stored
  • sstore/ystore (numpy.ndarray, shape (n,nt)) – Matrix containing stored information
  • d (numpy.ndarray, shape (n,)) – Vector with current direction
  • activeset (numpy.ndarray, shape (n,)) – Array containing the active indices.
Returns:

dnew – Vector with new direction

Return type:

numpy.ndarray, shape (n,)

oppy.conOpt.l_bfgs_b.l_bfgs_b(fhandle, dfhandle, x, low, up, options=None)

Limited-memory Broyden-Fletcher–Goldfarb-Shanno Box constraints.

This method solves optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\[l \leq x \leq u.\]
Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Initial guess. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • low (numpy.ndarray, shape (n,)) – Lower bound. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • up (numpy.ndarray, shape (n,)) – Upper bound. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • options (Options, optional) –

    Options for l_bfgs_b() method, represented as a Options object. The default is None. Possible options are:

    • tol_absfloat
      Absolute tolerance for termination condition, initialized with \(10^{-7}\).
    • tol_relfloat
      Relative tolerance for termination condition, initialized with \(10^{-7}\).

    Note

    The algorithm stops, if the norm of the current residual res[k] satisfies the inequality

    \[\mathrm{res}[k] < \mathrm{stopping\_criteria},\]

    where stopping_criteria is given by

    \[\mathrm{stopping\_criteria} = \mathrm{tol\_rel} \cdot \mathrm{res}[0] + \mathrm{tol\_abs},\]

    and res[0] is the norm of the first residual. If the value of stopping_criteria is bigger than \(0.1\), e.g. if the norm of the first residual is big and the relative tolerance is small, then this value is set to \(0.1\).

    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • linesearchcallable,
      Callable function that computes a steplength, initialized with armijo_goldstein().
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the functions gradient \(\nabla g\) in a vector \(x\). The default is None and therefore initialized to the standard norm in \(\mathbb{R}^n\), i.e.
      \[\mathrm{norm}: x \mapsto ||x||_2.\]

      Note

      The inner_prod function in the linesearch algorithm has to be changed appropriately to the norm function, i.e.

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • fargs:class`:tuple
      Tuple containing additional arguments for fhandle.
    • resultsstr or list
      String or list to specify the attributes for the Results object, initialized with 'normal'.
    • nsint
      Minimal number when storage is refreshed, initialized with \(0\) .
    • nsmaxint
      Maximal number when storage is refreshed, initialized with \(50\).
Returns:

results – Instance of the class Results. Represents the optimization results. Attributes of the object depend on the choice of the options.results parameter. With default choice the object has the following attributes:

  • xnumpy.ndarray, shape (n,)

    Solution vector of the optimization problem.

  • iterationsint

    Number of iterations the algorithm needed.

  • reslist

    List with the norm of the residuals.

Return type:

Results

References

    1. Kelley. Iterative Methods for Optimization see Kelley [Kel99].

Examples

Minimizing problem of the rosenbrock() function subject to

\begin{eqnarray} -1.5\leq x_1 \leq 1.5 \\ -0.5 \leq x_2 \leq 2.5 \end{eqnarray}

Importing the l_bfgs_b() method and defining the problem:

>>> import numpy as np
>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock
>>> from oppy.conOpt import l_bfgs_b
>>> low = np.array([-1.5, -0.5])
>>> up = np.array([1.5, 2.5])
>>> x0 = np.array([0.5, 0.5])

Execute l_bfgs_b() method:

>>> result = l_bfgs_b(rosenbrock, gradient_rosenbrock, x0, low, up)
>>> print(result.x)

The result is [1. 1.], which is the exact solution [1, 1].

oppy.conOpt.penalty module

Penalty method

This method solves optimization problems with equality and inequality constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\begin{eqnarray} e(x) & = 0 \\ g(x) & \leq 0 \end{eqnarray}

using penalty() method.

oppy.conOpt.penalty.penalty(fhandle, dfhandle, x, ehandle=None, dehandle=None, ghandle=None, dghandle=None, phandle=None, dphandle=None, options=None)

Penalty method.

This method solves optimization problems with equality and inequality constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\begin{eqnarray} e(x) & = 0 \\ g(x) & \leq 0. \end{eqnarray}
Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and fargs is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and fargs is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Initial guess. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • ehandle (callable(), optional) –

    The equality restriction function with \(\texttt{ehandle(x)=0}\).

    ehandle(x, *eargs) -> array_like, shape (me,)

    where x is an 1-D array with shape (n,) and eargs is a tuple of the fixed parameters needed to completely specify the function.

    Note

    If you pass ehandle, you must also pass the derivative dehandle. Also, you must pass either ehandle or ghandle.

  • dehandle (callable(), optional) –

    The gradient of the equality restriction function.

    dehandle(x, *eargs) -> array_like, shape (n,me).

    where x is an 1-D array with shape (n,) and eargs is a tuple of the fixed parameters needed to completely specify the function.

  • ghandle (callable(), optional) –

    The inequality restriction function with \(\texttt{ghandle(x)\leq 0}\).

    ghandle(x, *gargs) -> array_like, shape (mg,)

    where x is an 1-D array with shape (n,) and gargs is a tuple of the fixed parameters needed to completely specify the function.

    Note

    If you pass ghandle, you must also pass the derivative dghandle. Also, you must pass either ehandle or ghandle.

  • dghandle (callable(), optional) –

    The gradient of the inequality restriction function.

    dghandle(x, *gargs) -> array_like, shape (n,mg).

    where x is an 1-D array with shape (n,) and gargs is a tuple of the fixed parameters needed to completely specify the function.

  • options (Options, optional) –

    Options for penalty() method, represented as a Options object. The default is None. Possible options are:

    • tol_absfloat
      Absolute tolerance for termination condition, initialized with \(10^{-7}\).
    • tol_relfloat
      Relative tolerance for termination condition, initialized with \(10^{-7}\).

    Note

    The algorithm stops, if the norm of the current residual res[k] satisfies the inequality

    \[\mathrm{res}[k] < \mathrm{stopping\_criteria},\]

    where stopping_criteria is given by

    \[\mathrm{stopping\_criteria} = \mathrm{tol\_rel} \cdot \mathrm{res}[0] + \mathrm{tol\_abs},\]

    and res[0] is the norm of the first residual. If the value of stopping_criteria is bigger than \(0.1\), e.g. if the norm of the first residual is big and the relative tolerance is small, then this value is set to \(0.1\).

    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the functions \(e\) and \(g\) in a vector \(x\). The default is None and therefore initialized to the standard norm in \(\mathbb{R}^n\), i.e.
      \[\mathrm{norm}: x \mapsto ||x||_2.}.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • eargstuple
      Tuple containing additional arguments for ehandle.
    • gargstuple
      Tuple containing additional arguments for ghandle.
    • taufloat
      Intitial tolerance for the inner optimization method, intialized with \(10^{-3}\).
    • betafloat
      Update parameter for the penalty factor, intialized with \(2\).
    • c0float
      Intitial penalty factor, initialized with \(0.1\).
    • optimeth: callable()
      Callable function that computes a optimum, initialized with gradmethod.
    • resultsstr or list
      String or list to specify the attributes for the Results object, initialized with 'normal'.
Returns:

results – Instance of the class Results. Represents the optimization results. Attributes of the object depend on the choice of the options.results parameter. With default choice the object has the following attributes:

  • xnumpy.ndarray, shape (n,)

    Solution vector of the optimization problem.

  • iterationsint

    Number of iterations the algorithm needed.

  • reslist

    List with the norm of the residuals.

Return type:

Results

Examples

Minimizing problem of the rosenbrock() function subject to

\begin{eqnarray} (x_1-1)^3-x_2+1& \leq 0 \\ x_1+x_2-2 & \leq 0 \end{eqnarray}

Importing the penalty() method and defining the problem:

>>> import numpy as np
>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock
>>> from oppy.conOpt import penalty
>>> g = (lambda x: np.array([(x[0]-1)**3-x[1]+1, x[0]+x[1]-2]))
>>> dg = (lambda x: np.array([[3*(x[0]-1)**2, -1],
                            [1, 1]]))
>>> f = rosenbrock
>>> df = gradient_rosenbrock
>>> x0 = np.array([0.5, -0.5])

Execute penalty() method:

>>> result = penalty(fhandle=f, dfhandle=df, x=x0, ghandle=g,
                                  dghandle=dg)
>>> print(result.x)

The result is [0.99999989 0.99999986], which is close to the exact solution [1, 1].

oppy.conOpt.projected_armijo module

Projected Armijo linesearch

This linesearch method can be used in optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]
\[\text{s. t. } l\leq x \leq u\]

for finding an optimal stepsize length \(t\) in every optimization iteration.

oppy.conOpt.projected_armijo.projected_armijo(fhandle, x, d, low, up, options=None)

Projected Armijo stepsize algorithm is a linesearch method for box constrained optimization problems. The linesearch terminates if the sufficient decrease condition

\[f(x(t)) \leq f(x)-c_1/t \langle (x-x(t)),(x-x(t)) \rangle\]

is violated. Here, \(\langle \cdot,\cdot\rangle\) denotes the inner product and \(x(t) = P(x + t\cdot d)\), where the projector \(P\) is given by:

\[\begin{split}P(y)=\begin{cases} y \quad &\mathrm{,if}\quad l \leq y \leq u\\ l \quad &\mathrm{,if}\quad l >y\\ u \quad &\mathrm{,if}\quad u <y\\ \end{cases}\end{split}\]

where \(l\) is the lower and \(u\) the upper bound.

Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Current iteration point. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • d (numpy.ndarray, shape (n,)) – New direction. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • low (numpy.ndarray, shape (n,)) –
    Lower bound for the projection. Array of real elements of size (n,),
    where ‘n’ is the number of independent variables.
    upnumpy.ndarray, shape (n,)
    Upper bound for the projection. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • options (Options, optional) –

    Options for newton method, represented as a Options object. The default is None. Possible options are:

    • initstepfloat
      Initial stepsize, initialized to \(1\).
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • c1float
      Weighting parameter, initialized to \(10^{-2}\).
    • betafloat
      Parameter to reduce stepsize, initialized to \(1/2\).
    • nmaxint or float
      Maximum number of iterations, initialized with \(30\).
    • inner_prodcallable()
      Callable function that computes the inner product of the difference \(x - x(t)\). The default is None and therefore initialized to the standard \(\mathbb{R}^n\) scalarproduct, i.e.
      \[\mathrm{inner\_prod} : (x,x(t)) \mapsto\]
      \[(x-x(t))^T(x-x(t))\]

      Note

      The norm function in the optimizer has to be changed appropriately to the inner_prod function in the linesearch algorithm, i.e.

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}\]
Returns:

  • t (float) – Final steplength

  • flag (int) – The flag indicated, whether the linesearch was successful, possible outputs are:

    1: Final steplength complies with sufficient decrease.
    0: Algorithm did not terminate and ended with nmax.

Examples

Minimizing problem of the rosenbrock() function

>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock

The l_bfgs_b() method was chosen as the minimizing method

>>> from oppy.conOpt import l_bfgs_b
>>> from oppy.conOpt import armijo_goldstein
>>> import numpy as np

Defining the armijo_goldstein() linesearch with default options:

>>> def linesearch(x, d, low, up):
>>>     t, flag = armijo_goldstein(fhandle, dfhandle, x, d, low, up)
>>>     return t, flag

The options for the l_bfgs_b() are set as default as well. Input data for the l_bfgs_b() method:

>>> x0 = np.array([0.5, 0.5])
>>> low = np.array([-1.5, -0.5])
>>> up = np.array([1.5, 2.5])

Execute l_bfgs_b() method with armijo_goldstein() linesearch:

>>> result = l_bfgs_b(rosenbrock, gradient_rosenbrock, x0, low, up)
>>> print(result.x)

The result is [1., 1.], which is the exact solution [1, 1].

oppy.conOpt.projected_gradient module

Projected Gradient method

This method solves optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\[l \leq x \leq u\]

using projected_gradient() method.

oppy.conOpt.projected_gradient.projected_gradient(fhandle, dfhandle, x, low, up, options=None)

Projected Gradient method.

This method solves optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\[l \leq x \leq u.\]
Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Initial guess. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • low (numpy.ndarray, shape (n,)) – Lower bound. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • up (numpy.ndarray, shape (n,)) – Upper bound. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • options (Options, optional) –

    Options for projected_gradient() method, represented as a Options object. The default is None. Possible options are:

    • tol_absfloat
      Absolute tolerance for termination condition, initialized with \(10^{-7}\).
    • tol_relfloat
      Relative tolerance for termination condition, initialized with \(10^{-7}\).

    Note

    The algorithm stops, if the norm of the current residual res[k] satisfies the inequality

    \[\mathrm{res}[k] < \mathrm{stopping\_criteria},\]

    where stopping_criteria is given by

    \[\mathrm{stopping\_criteria} = \mathrm{tol\_rel} \cdot \mathrm{res}[0] + \mathrm{tol\_abs},\]

    and res[0] is the norm of the first residual. If the value of stopping_criteria is bigger than \(0.1\), e.g. if the norm of the first residual is big and the relative tolerance is small, then this value is set to \(0.1\).

    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • linesearchcallable(),
      Callable function that computes a steplength, initialized with armijo_goldstein().
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the functions gradient \(\nabla g\) in a vector \(x\). The default is None and therefore initialized to the standard norm in \(\mathbb{R}^n\), i.e.
      \[\mathrm{norm}: x \mapsto ||x||_2.\]

      Note

      The inner_prod function in the linesearch algorithm has to be changed appropriately to the norm function, i.e.

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • resultsstr or list
      String or list to specify the attributes for the Results object, initialized with 'normal'.
Returns:

results – Instance of the class Results. Represents the optimization results. Attributes of the object depend on the choice of the options.results parameter. With default choice the object has the following attributes:

  • xnumpy.ndarray, shape (n,)

    Solution vector of the optimization problem.

  • iterationsint

    Number of iterations the algorithm needed.

  • reslist

    List with the norm of the residuals.

Return type:

Results

References

    1. Kelley. Iterative Methods for Optimization, see Kelley [Kel99].

Examples

Minimizing problem of the rosenbrock() function subject to

\begin{eqnarray} -1.5\leq x_1 \leq 1.5 \\ -0.5 \leq x_2 \leq 2.5 \end{eqnarray}

Importing the projected_gradient() method and defining the problem:

>>> import numpy as np
>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock
>>> from oppy.conOpt import projected_gradient
>>> x0 = np.array([0.5, 0.5])
>>> low = np.array([-1.5, -0.5])
>>> up = np.array([1.5, 2.5])

Execute projected_gradient() method:

>>> result = projected_gradient(rosenbrock, gradient_rosenbrock, x0,
>>>                             low, up)
>>> print(result.x)

The maximum number was reached. Therefore the algorithm stopped and the result is [0.76556834, 0.58547492], which is not close to the exact solution [1, 1].

oppy.conOpt.projected_newton module

Projected Newton-CG (GMRES)

This method solves optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\[l \leq x \leq u\]

using projected_newton() method.

oppy.conOpt.projected_newton.projected_newton(fhandle, dfhandle, ddfhandle, x, low, up, options=None)

Projected Newton method.

This method solves optimization problems with box constraints, e.g.

\[\min_{x \in \mathbb{R}^n} f(x)\]

subject to

\[l \leq x \leq u.\]
Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • ddfhandle (callable()) –

    The action of the Hessian of the objective function in a point x in direction u.

    ddfhandle(u, x, Active, Inactive) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,), u is an 1-D array with shape (n,) and Active and Inactive are are 1-D array with shape (n,) as well.

    Note

    If you have an analytical represantation of the Hessian you can calculate the the action of the Hessian as

    \[ddfhandle(u, x, Act, Inact) = Act*u + Inact*((Hessian(x) \cdot (Inact*u)))\]

    where \(*\) represents to dotwise multiplication and \(\cdot\) the Matrix-Vecotr multiplication.

  • x (numpy.ndarray, shape (n,)) – Initial guess. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • low (numpy.ndarray, shape (n,)) – Lower bound. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • up (numpy.ndarray, shape (n,)) – Upper bound. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • options (Options, optional) –

    Options for projected_newton() method, represented as a Options object. The default is None. Possible options are:

    • tol_absfloat
      Absolute tolerance for termination condition, initialized with \(10^{-7}\).
    • tol_relfloat
      Relative tolerance for termination condition, initialized with \(10^{-7}\).

    Note

    The algorithm stops, if the norm of the current residual res[k] satisfies the inequality

    \[\mathrm{res}[k] < \mathrm{stopping\_criteria},\]

    where stopping_criteria is given by

    \[\mathrm{stopping\_criteria} = \mathrm{tol\_rel} \cdot \mathrm{res}[0] + \mathrm{tol\_abs},\]

    and res[0] is the norm of the first residual. If the value of stopping_criteria is bigger than \(0.1\), e.g. if the norm of the first residual is big and the relative tolerance is small, then this value is set to \(0.1\).

    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • linesearchcallable,
      Callable function that computes a steplength, initialized with armijo_goldstein().
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the functions gradient \(\nabla g\) in a vector \(x\). The default is None and therefore initialized to the standard norm in \(\mathbb{R}^n\), i.e.
      \[\mathrm{norm}: x \mapsto ||x||_2.\]

      Note

      The inner_prod function in the linesearch algorithm has to be changed appropriately to the norm function, i.e.

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • use:class: str
      String to specify how the newton equation shall be solved. You can decide between cg() or gmres().
    • resultsstr or list
      String or list to specify the attributes for the Results object, initialized with 'normal'.
Returns:

results – Instance of the class Results. Represents the optimization results. Attributes of the object depend on the choice of the options.results parameter. With default choice the object has the following attributes:

  • xnumpy.ndarray, shape (n,)

    Solution vector of the optimization problem.

  • iterationsint

    Number of iterations the algorithm needed.

  • reslist

    List with the norm of the residuals.

Return type:

Results

References

    1. Kelley. Iterative Methods for Optimization, see Kelley [Kel99].

Examples

Minimizing problem of the rosenbrock() function subject to

\begin{eqnarray} -1.5\leq x_1 \leq 1.5 \\ -0.5 \leq x_2 \leq 2.5 \end{eqnarray}

Importing the projected_newton() method and defining the problem:

>>> import numpy as np
>>> from oppy.tests.costfunctions import (rosenbrock, gradient_rosenbrock,
>>>                                  hessian_rosenbrock)
>>> from oppy.conOpt import projected_newton
>>> x0 = np.array([0.5, 0.5])
>>> low = np.array([-1.5, -0.5])
>>> up = np.array([1.5, 2.5])

As written in the description of the method we need the action of the Hessian.

>>> action_rosenbrock = lambda u, x, Act, Inact: Act*u + Inact*((hessian_rosenbrock(x)@(Inact*u)))

Execute projected_newton() method:

>>> result = projected_newton(rosenbrock, gradient_rosenbrock,
>>>                           action_rosenbrock, x0, low, up)
>>> print(result.x)

The result is [1. 1.], which is the exact solution [1, 1].

oppy.conOpt.sqp_bfgs module

Sequential Quadratic Programming with damped BFGS.

This method solves optimization problems with box constraints, e.g.

\[\min f(x)\]

subject to

\[e(x) = 0\]

using a sequential quadratic qrogramming method with damped BFGS (sqp_bfgs()) update to aproximate the second derivative.

oppy.conOpt.sqp_bfgs.sqp_bfgs(fhandle, dfhandle, x, ehandle=None, dehandle=None, options=None)

Sequential Quadratic Programming with damped BFGS.

With this method you can solve problems of the form

\[\min f(x)\]

subject to

\[e(x) = 0\]

using a sequential quadratic qrogramming method with damped BFGS update to aproximate the second derivative.

Parameters:
  • fhandle (callable()) –

    The objective function to be minimized.

    fhandle(x, *fargs) -> float

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • dfhandle (callable()) –

    The gradient of the objective function.

    dfhandle(x, *fargs) -> array_like, shape (n,).

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • x (numpy.ndarray, shape (n,)) – Initial guess. Array of real elements of size (n,), where ‘n’ is the number of independent variables.
  • ehandle (callable(), optional) –

    The restriction function.

    ehandle(x, *eargs) -> array_like, shape (m,)

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

    Note

    In a future realease, inequality restrictions will be allowed. Currently, the restrictions are optional, but they can only be omitted when the inequality restrictions are implemented.

  • dehandle (callable(), optional) –

    The gradient of the restriction function.

    dehandle(x, *eargs) -> array_like, shape (n,m).

    where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

  • options (Options, optional) –

    Options for sqp_bfgs(), represented as a Options object. The default is None. Class Attributes containing

    • tol_absfloat
      Absolute tolerance for termination condition, initialized with \(10^{-7}\).
    • tol_relfloat
      Relative tolerance for termination condition, initialized with \(10^{-7}\).

    Note

    The algorithm stops, if the norm of the current residual res[k] satisfies the inequality

    \[\mathrm{res}[k] < \mathrm{stopping\_criteria},\]

    where stopping_criteria is given by

    \[\mathrm{stopping\_criteria} = \mathrm{tol\_rel} \cdot \mathrm{res}[0] + \mathrm{tol\_abs},\]

    and res[0] is the norm of the first residual. If the value of stopping_criteria is bigger than \(0.1\), e.g. if the norm of the first residual is big and the relative tolerance is small, then this value is set to \(0.1\).

    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the first order criterium, e.g.
      \[\Vert \nabla f(x) + \lambda^\top \nabla e(x) \Vert + \Vert e(x) \Vert\]

      in a vector \(x\). The default is None and therefore initialized to the standard norm in \(\mathbb{R}^n\), i.e.

      \[\mathrm{norm}: x \mapsto ||x||_2\]

      Note

      The inner_prod function in the linesearch algorithm has to be changed appropriately to the norm function, i.e.

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • eargstuple
      Tuple containing additional arguments for ehandle.
    • lamnumpy.ndarray, shape (me,)
      Initial lagrange mulitplicator for the equality constraints, initialized with zeros.
    • resultsstr or list
      String or list to specify the attributes for the Results object, initialized with 'normal'.
Returns:

results – Instance of the class Results. Represents the optimization results. Attributes of the object depend on the choice of the options.results parameter. With default choice the object has

  • xnumpy.ndarray, shape (n,)

    Solution vector of the optimization problem.

  • iterationsint

    Number of iterations the algorithm needed.

  • reslist

    List with the norm of the residuals.

Return type:

Results

Examples

Minimizing problem of the rosenbrock() function subject to

\begin{eqnarray} (x_1-1)^3-x_2+1& \leq 0 \\ x_1+x_2-2 & \leq 0 \end{eqnarray}

Importing the sqp_bfgs() method and defining the problem:

>>> import numpy as np
>>> from oppy.tests.costfunctions import rosenbrock, gradient_rosenbrock
>>> from oppy.conOpt import sqp_bfgs
>>> ehandle = (lambda x: np.array([(x[0]-1)**3-x[1]+1, x[0]+x[1]-2]))
>>> dehandle = (lambda x: np.array([[3*(x[0]-1)**2, -1],
>>>                                 [1, 1]]))
>>> fhandle = rosenbrock
>>> dfhandle = gradient_rosenbrock
>>> x0 = np.array([0.5, -0.5])

Execute sqp_bfgs():

>>> result = sqp_bfgs(fhandle, dfhandle, x0, ehandle, dehandle)
>>> print(result.x)

The result is [1., 1.], which is the exact solution [1, 1].