Unconstrained Optimization

unconOpt: unconstrained optimization

Subpackage which provide some methods for unconstrained optimization, e.g:

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

Right now we can solve this kind of problems with line search based first- and second-order methods.

  • Line Search Methods
    1. Armijo
    2. Wolfe-Powell
  • Optimization Methods
    1. Gradient Method
    2. Newton Method
    3. Nonlinear CG (with different strategies like Fletcher-Reves)
    4. Quasi-Newton Methods (with different strategies like BFGS, Broyden, DFP, …)

Documentation is available in the docstrings and online here.

Methods

oppy.unconOpt.armijo module

Armijo Linesearch.

Given a starting point \(x\) and a search direction \(d\) the stepsize algorithm armijo() computes a stepsize \(t\) which satisfies the sufficient decrease condition:

\[f(x+td) \leq f(x)+c_1 t \nabla f(x)^T d\]
oppy.unconOpt.armijo.armijo(fhandle, dfhandle, x, d, options=None)

Compute a stepsize satisfying the Armijo conditon.

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.
  • options (Options, optional) –

    Options for armijo, represented as a Options object. The default is None. Possible options are

    • initstepfloat,
      Initial stepsize, initialized to \(1\).
    • fargsfloat,
      Tuple containing additional arguments for fhandle.
    • c1float,
      Weightingparameter, initialized to \(10^{-2}\).
    • betafloat,
      Parameter to reduce stepsize, initialized to \(1/2\).
    • nmaxint or float,
      Maximum number of iterations, initialized to \(30\).
    • inner_prodcallable(),
      Callable function that computes the scalarproduct between a search direction \(d \in \mathbb{R}^n\) and the gradient \(\nabla f\) in a vector \(x \in \mathbb{R}^n\). The default is None and therefore initialized to the standard \(\mathbb{R}^n\) scalarproduct, i.e.
      \[\mathrm{inner\_prod} : (x,d) \mapsto x^T d\]

      Note

      Please note, that if you want to change the inner_prod function in the linesearch algorithm, that you have to change the norm function in your optimizer to

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

  • t (float) – Final steplength

  • flag (int) – The flag indicates, 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 negative cosine function

>>> from oppy.tests.costfunctions import negativeCosine,\
... gradient_negativeCosine, hessian_negativeCosine

The newton() method was choosen as the minimizing method

>>> from oppy.unconOpt.armijo import armijo
>>> from oppy.unconOpt.newton import newton
>>> import numpy as np
>>> from oppy.options.options import Options

Set all Options manually to their default value:

>>> options_armijo = Options()
>>> options_armijo.initstep = 1
>>> options_armijo.c1 = 1e-2
>>> options_armijo.beta = 1/2
>>> options_armijo.nmax = 30

Defining the armijo() linesearch:

>>> def linesearch(x,d):
...     t, flag = armijo(negativeCosine, gradient_negativeCosine, x, d,\
...                      options_armijo)
...     return(t, flag)

Override default linesearch of the newton() method with manually defined linesearch method:

>>> options_newton = Options()
>>> options_newton.linesearch = linesearch
>>> options_newton.results = "all"

Set an initial guess for the newton() method:

>>> x0 = np.array([1.655])
>>> results_newton = newton(negativeCosine, gradient_negativeCosine,\
...                     hessian_negativeCosine,x0, options_newton)

Show optimization results:

>>> print(results_newton)
          X: array([[1.65500000e+00,  6.58543035e-01, -1.15229962e-01,
                     5.12729027e-04, -4.49306297e-11]])
 additional: [{'linesearch': <function linesearch at 0x7fc8352c8940>,
               'results': 'all', 'tol_abs': 1e-06, 'tol_rel': 1e-06,
               'nmax': 100, 'disp': False, 'inner_prod': None, 'fargs': (),
               'parent_results': None, 'alpha1': 0.1, 'alpha2': 0.1}]
 iterations: 4
 linesearch: array([0, 1, 1, 1, 1])
        res: array([9.96456965e-01, 6.11965211e-01, 1.14975128e-01,
                    5.12729004e-04, 4.49306297e-11])
          x: array([-4.49306297e-11])

oppy.unconOpt.gradmethod module

Gradientscheme.

This module contains the function gradmethod() to solve the minimzation problem

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

using the gradient descent method.

oppy.unconOpt.gradmethod.gradmethod(fhandle, dfhandle, x, options=None)

Gradient method with linesearch algorithm.

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.
  • options (Options, optional) –

    Options for gradmethod, 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 optimization 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 wolfe().
    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the gradient \(\nabla f\) 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

      If you want to change the norm, please remember, that you have to change the inner_prod function in your linesearch algorithm, such that norm is computed by

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • usestr
      If use == 'scaled', the descent direction is computed by \(d = -\frac{\mathrm{dfhandle}(x)} {\mathrm{norm}(\mathrm{dfhandle}(x))}\) otherwise by \(d = -\mathrm{dfhandle}(x)\).
    • 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 negative cosine function

>>> from oppy.tests.costfunctions import negativeCosine, \
...      gradient_negativeCosine

The minimizing method was chosen as the gradmethod() method

>>> from oppy.unconOpt import gradmethod
>>> from oppy.unconOpt import armijo
>>> import numpy as np
>>> from oppy.options import Options

Define a custom linesearch method as armijo() linesearch

>>> def linesearch(x,d):
...     t, flag = armijo(negativeCosine, gradient_negativeCosine, x, d)
...     return(t, flag)

Set options for the gradmethod() method

>>> options_gradmethod = Options()
>>> options_gradmethod.linesearch = linesearch
>>> options_gradmethod.nmax = 100
>>> options_gradmethod.tol_abs = 1e-6
>>> options_gradmethod.tol_rel = 1e-6

Input data for the gradmethod() method:

>>> x0 = np.array([1.1655])

Execute gradmethod() algorithm with armijo() linesearch:

>>> results_gradmethod = gradmethod(negativeCosine,
...                                 gradient_negativeCosine,
...                                 x0, options_gradmethod)
>>> print(results_gradmethod)
  iterations: 12
        res: array([9.18985598e-01, 1.64745520e-01, 8.43994774e-02,
                    4.04889292e-02, 2.19982254e-02, 9.24986809e-03,
                    6.37495682e-03, 1.43749950e-03, 5.15624977e-04,
                    4.60937484e-04, 2.73437500e-05, 3.17382812e-06,
        6.40869141e-07])
          x: array([-6.40869141e-07])

oppy.unconOpt.newton module

Newton’s method.

This module contains the function newton() to solve the minimzation problem

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

using the (globalized) Newton’s method.

oppy.unconOpt.newton.newton(fhandle, dfhandle, ddfhandle, x, options=None)

Newton’s method with linesearch algorithm.

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 Hessian of the objective function.

    ddfhandle(x, *fargs) -> array_like, shape (n,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.
  • options (Options, optional) –

    Options for Newton’s 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 optimization 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
      Maximum number of iterations, initialized with 100.
    • linesearchcallable()
      Callable function that computes a steplength initialized with wolfe().
    • dispbool
      If disp == True, some details about the progress will shown during the optimization
    • normcallable()
      Callable function that computes the norm of the gradient \(\nabla f\) 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

      If you want to change the norm, please remember, that you have to change the inner_prod function in your linesearch algorithm, such that norm is computed by

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • fargstuple
      Tuple containing additional arguments for fhandle.
    • alpha1float
      Parameter for the globalization strategy, initialized with \(10^{-1}\).
    • alpha2float
      Parameter for the globalization strategy, initialized with \(10^{-1}\).
    • 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 negative cosine function

>>> from oppy.tests.costfunctions import negativeCosine,\
... gradient_negativeCosine, hessian_negativeCosine
>>> from oppy.options.options import Options

The minimizing method was chosen as the newton() method

>>> from oppy.unconOpt.newton import newton
>>> from oppy.unconOpt.armijo import armijo
>>> import numpy as np

Define a custom linesearch method as armijo() linesearch

>>> def linesearch(x, d):
...      t, flag = armijo(negativeCosine, gradient_negativeCosine, x, d)
...      return(t, flag)

The Options for the newton() method are set as

>>> options_newton = Options()
>>> options_newton.linesearch = linesearch
>>> options_newton.nmax = 100
>>> options_newton.tol_rel = 1e-6
>>> options_newton.tol_abs = 1e-6

Input data for the newton() method

>>> x0 = np.array([1.1655])
Execute newton() method with armijo() linesearch
>>> results_newton = newton(negativeCosine, gradient_negativeCosine,\
...                         hessian_negativeCosine, x0, options_newton)
>>> print(results_newton)
 iterations: 2
        res: array([9.18985598e-01, 1.35623556e-04, 8.31541894e-13])
          x: array([-8.31541894e-13])

oppy.unconOpt.nonlinear_cg module

Nonlinear Conjugated Gradient (with restart).

This module contains the function nonlinear_cg() and all helper functions to solve the minimzation problem

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

using the nonlinear conjugated gradient method (with restart).

oppy.unconOpt.nonlinear_cg.nonlinear_cg(fhandle, dfhandle, x, options=None)

Nonlinear Conjugated Gradient (with restart).

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.
  • options (Options, optional) –

    Options for nonlinear_cg 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 optimization 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 wolfe().

      Note

      A restart is currently only executed if the linesearch algorithm failes (e.g. if the returned step width does not satisfy the desired condition(s)).

    • dispbool
      If disp == True, some details about the progress will be shown during the optimization.
    • normcallable()
      Callable function that computes the norm of the gradient \(\nabla f\) 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

      If you want to change the norm, please remember, that you have to change the inner_prod function in your linesearch algorithm, such that norm is computed by

      \[\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'.
    • updatestr
      Update rule for the subsequent conjugate directions. Default is set to 'fr_pr'.
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 negative cosine function

>>> from oppy.tests.costfunctions import negativeCosine, \
...           gradient_negativeCosine

The minimizing method was chosen as the nonlinear CG method

>>> from oppy.unconOpt.nonlinear_cg import nonlinear_cg
>>> from oppy.unconOpt.armijo import armijo
>>> import numpy as np
>>> from oppy.options.options import Options

Define a custom linesearch method as armijo() linesearch

>>> def linesearch(x,d):
...      t, flag = armijo(negativeCosine, gradient_negativeCosine,x,d)
...      return(t,flag)

Set Options for the nonlinear_cg() algorithm

>>> options_frpr = Options()
>>> options_frpr.linesearch = linesearch
>>> options_frpr.tol_abs = 1e-6
>>> options_frpr.tol_rel = 1e-6
>>> options_frpr.nmax = 100

Input data for the nonlinear_cg() method:

>>> x0 = np.array([1.1655])

Execute nonlinear_cg() method with armijo() linesearch

>>> results_frpr = nonlinear_cg(negativeCosine,gradient_negativeCosine,x0,
...                      options_frpr)
>>> print(results_frpr)
 iterations: 8
        res: array([9.18985598e-01, 2.44025224e-01, 6.72362857e-02,
                    1.36567083e-02, 2.21296590e-03, 3.00499737e-04,
                    3.52641138e-05, 3.65266310e-06, 3.39154601e-07])
          x: array([3.39154601e-07])
oppy.unconOpt.nonlinear_cg.fr_pr(grad_old, grad_new)

Compute the beta parameter for nonlinear_cg().

Parameters:
  • grad_old (numpy.ndarray) – Gradient of the function evaluated at the previous iterate.
  • grad_new (numpy.ndarray) – Gradient of the function evaluated at the current iterate.
Returns:

beta\(\beta\) parameter for the subsequent conjugate direction.

Return type:

float

oppy.unconOpt.quasinewton module

Globalised quasi-Newton.

This module contains the function quasinewton() and all helper functions (computation of the approximated hessian matrix) to solve the minimization problem

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

using the (globalized) quasi-Newton method.

oppy.unconOpt.quasinewton.quasinewton(fhandle, dfhandle, x, options=None)

Quasi-Newton method with linesearch algorithm.

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.
  • options (Options, optional) –

    Options for quasi-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 optimization 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\).

    • fargstuple
      Tuple containing additional arguments for fhandle.
    • nmaxint or float
      Maximum number of iterations, initialized with 100.
    • linesearchcallable(),
      Callable function that computes a steplength initialized with wolfe().
    • dispbool
      If disp == True, some details about the progress will shown during the optimization
    • normcallable()
      Callable function that computes the norm of the gradient \(\nabla f\) 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

      If you want to change the norm, please remember, that you have to change the inner_prod function in your linesearch algorithm, such that norm is computed by

      \[\mathrm{norm} : x \mapsto \sqrt{\mathrm{inner\_prod}(x,x)}.\]
    • updatestr
      Which quasi-Newton to update the approx hessian, initialized with bfgs (also broyden, dfp and sr1 are available).
    • 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 negative cosine function

>>> from oppy.tests.costfunctions import negativeCosine,\
... gradient_negativeCosine
>>> from oppy.options.options import Options

The minimizing method was chosen as the quasinewton() method

>>> from oppy.unconOpt.quasinewton import quasinewton
>>> from oppy.unconOpt.armijo import armijo
>>> import numpy as np

Define a custom linesearch method as armijo() linesearch

>>> def linesearch(x, d):
...     t, flag = armijo(negativeCosine, gradient_negativeCosine, x, d)
...     return(t, flag)

The Options for the quasinewton() method are set as

>>> options_quasinewton = Options()
>>> options_quasinewton.linesearch = linesearch
>>> options_quasinewton.nmax = 100
>>> options_quasinewton.tol_rel = 1e-6
>>> options_quasinewton.tol_abs = 1e-6

Input data for the quasinewton() method

>>> x0 = np.array([1.1655])

Execute quasinewton() method with armijo() linesearch

>>> result_quasinewton = quasinewton(negativeCosine,\
...                                  gradient_negativeCosine, x0,\
...                                  options_quasinewton)
>>> print(result_quasinewton)
 iterations: 4
        res: array([9.18985598e-01, 2.44025224e-01, 8.56307228e-02,
                    5.68860736e-04, 6.92877097e-07])
          x: array([-6.92877097e-07])
oppy.unconOpt.quasinewton.bfgs(s, y, H)

BFGS update for quasinewton().

Calculating the bfgs update during the formula:

\[H_{k+1} = H_k + yy^T/y^Ts - H_ks(H_ks)^T/s^TH_ks\]
Parameters:
  • s (numpy.ndarray, shape(n,)) – Represents \(x_{k+1} - x_k\).
  • y (numpy.ndarray, shape(n,)) – Represents \(\nabla f(x_{k+1}) - \nabla f(x_k)\).
  • H (numpy.ndarray, shape(n,n)) – Represents the old approximate Hessian \(H_k\).
Returns:

Represents the new approxiamte Hessian \(H_{k+1}\).

Return type:

numpy.ndarray, shape(n,n)

oppy.unconOpt.quasinewton.broyden(s, y, H)

Broyden update formula for quasinewton().

Calculating the bfgs update during the formula:

Parameters:
  • s (numpy.ndarray, shape(n,)) – Represents \(x_{k+1} - x_k\).
  • y (numpy.ndarray, shape(n,)) – Represents \(\nabla f(x_{k+1}) - \nabla f(x_k)\).
  • H (numpy.ndarray, shape(n,n)) – Represents the old approximate Hessian \(H_{k}\).
Returns:

Represents the new approxiamte Hessian \(H_{k+1}\).

Return type:

numpy.ndarray, shape(n,n)

oppy.unconOpt.quasinewton.dfp(s, y, H)

DFP update formula for quasinewton().

Calculating the bfgs update during the formula:

Parameters:
  • s (numpy.ndarray, shape(n,)) – Represents \(x_{k+1} - x_k\).
  • y (numpy.ndarray, shape(n,)) – Represents \(\nabla f(x_{k+1}) - \nabla f(x_k)\).
  • H (numpy.ndarray, shape(n,n)) – Represents the old approximate Hessian \(H_k\).
Returns:

Represents the new approxiamte Hessian \(H_{k+1}\).

Return type:

numpy.ndarray, shape(n,n)

oppy.unconOpt.quasinewton.sr1(s, y, H)

SR1 update formula for quasinewton().

Calculating the bfgs update during the formula:

Parameters:
  • s (numpy.ndarray, shape(n,)) – Represents \(x_{k+1} - x_k\).
  • y (numpy.ndarray, shape(n,)) – Represents \(\nabla f(x_k+1) - \nabla f(x_k)\).
  • H (numpy.ndarray, shape(n,n)) – Represents the old approximate Hessian \(H_k\).
Returns:

Represents the new approxiamte Hessian \(H_{k+1}\).

Return type:

numpy.ndarray, shape(n,n)

oppy.unconOpt.quasinewton.broyden_fam(s, y, H, phi)

At the moment not really tested.

oppy.unconOpt.wolfe module

(Powell-)Wolfe Linesearch.

Given a starting point \(x\) and a search direction \(d\) the stepsize algorithm wolfe() computes a stepsize \(t\) which satisfies the sufficient decrease conditions:

\[(1) \, f(x+td) \leq f(x) + c_1 t d^T \nabla f(x)\]
\[(2) \, d^T \nabla f(x) \geq c_2 d^T \nabla f(x)\]
oppy.unconOpt.wolfe.wolfe(fhandle, dfhandle, x, d, options=None)

Compute a stepsize satisfying the (Powell-)Wolfe conditon.

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.
  • options (Options, optional) –

    Options for wolfe, 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,
      Weightingparameter, initialized to 10^-4.
    • c2float,
      Weightingparameter, initialized to 9/10.
    • betafloat,
      Parameter to reduce stepsize, initialized to 1/2.
    • nmaxint or float,
      Maximum number of iterations, initialized to 30.
    • inner_prodcallable(),
      Callable function that computes the inner product between a search direction \(d \in \mathbb{R}^n\) and the gradient \(\nabla f\) for a vector \(x \in \mathbb{R}^n\). The default is None and therefore initialized to the standard \(\mathbb{R}^n\) scalarproduct, i.e.
      \[\mathrm{inner\_prod} : (x,d) \mapsto x^T d\]

      Note

      Please note, that if you want to change the inner_prod function in the linesearch algorithm, that you have to change the norm function in your optimizer to

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

  • t (float) – Final steplength

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

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

References

  1. Ulbrich, S. Ulbrich. Nichtlineare Optimierung, see Ulbrich and Ulbrich [UU12] p. 39.

Examples

Minimizing problem of the negative cosine function

>>> from oppy.tests.costfunctions import negativeCosine, \
... gradient_negativeCosine, hessian_negativeCosine
>>> from oppy.unconOpt.newton import newton
>>> from oppy.unconOpt.wolfe import wolfe
>>> from oppy.options.options import Options
>>> import numpy as np

The Options for the wolfe linesearch method are set as:

>>> options_wolfe = Options()
>>> options_wolfe.initstep = 1
>>> options_wolfe.beta = 0.5
>>> options_wolfe.c1 = 1e-4
>>> options_wolfe.c2 = 2/10
>>> options_wolfe.nmax = 30
>>> def linesearch(x,d):
...     t, flag = wolfe(negativeCosine, gradient_negativeCosine,x, d,\
...                     options_wolfe)
...     return(t, flag)

Override default linesearch of the newton() method with wolfe() linesearch and set returns option to 'all'

>>> options_newton = Options()
>>> options_newton.linesearch = linesearch
>>> options_newton.results = 'all'

Set initial guess for the newton() method

>>> x0 = np.array([1.655])
>>> results_newton = newton(negativeCosine, gradient_negativeCosine,\
...                     hessian_negativeCosine, x0, options_newton)

Print optimization Results

>>> print(results_newton)
          X: array([[1.65500000e+00, -3.37913930e-01,  1.34775758e-02,
                     -8.16104275e-07]])
 additional: [{'linesearch': <function linesearch at 0x7f88c69c8940>,
               'results': 'all', 'tol_abs': 1e-06, 'tol_rel': 1e-06,
               'nmax': 100, 'disp': False, 'inner_prod': None, 'fargs': (),
               'parent_results': None, 'alpha1': 0.1, 'alpha2': 0.1}]
 iterations: 3
 linesearch: array([0., 2., 1., 1.])
        res: array([9.96456965e-01, 3.31519715e-01, 1.34771678e-02,
                    8.16104275e-07])
          x: array([-8.16104275e-07])