Minimize a function using a nonlinear conjugate gradient algorithm.
Parameters : | f : callable, f(x, *args)
x0 : ndarray
fprime : callable, fprime(x, *args), optional
args : tuple, optional
gtol : float, optional
norm : float, optional
epsilon : float or ndarray, optional
maxiter : int, optional
full_output : bool, optional
disp : bool, optional
retall : bool, optional
callback : callable, optional
|
---|---|
Returns : | xopt : ndarray
fopt : float, optional
func_calls : int, optional
grad_calls : int, optional
warnflag : int, optional
allvecs : list of ndarray, optional
|
See also
Notes
This conjugate gradient algorithm is based on that of Polak and Ribiere [R78].
Conjugate gradient methods tend to work better when:
References
[R78] | (1, 2) Wright & Nocedal, “Numerical Optimization”, 1999, pp. 120-122. |
Examples
Example 1: seek the minimum value of the expression a*u**2 + b*u*v + c*v**2 + d*u + e*v + f for given values of the parameters and an initial guess (u, v) = (0, 0).
>>> args = (2, 3, 7, 8, 9, 10) # parameter values
>>> def f(x, *args):
... u, v = x
... a, b, c, d, e, f = args
... return a*u**2 + b*u*v + c*v**2 + d*u + e*v + f
>>> def gradf(x, *args):
... u, v = x
... a, b, c, d, e, f = args
... gu = 2*a*u + b*v + d # u-component of the gradient
... gv = b*u + 2*c*v + e # v-component of the gradient
... return np.asarray((gu, gv))
>>> x0 = np.asarray((0, 0)) # Initial guess.
>>> from scipy import optimize
>>> res1 = optimize.fmin_cg(f, x0, fprime=gradf, args=args)
>>> print 'res1 = ', res1
Optimization terminated successfully.
Current function value: 1.617021
Iterations: 2
Function evaluations: 5
Gradient evaluations: 5
res1 = [-1.80851064 -0.25531915]
Example 2: solve the same problem using the minimize function. (This myopts dictionary shows all of the available options, although in practice only non-default values would be needed. The returned value will be a dictionary.)
>>> opts = {'maxiter' : None, # default value.
... 'disp' : True, # non-default value.
... 'gtol' : 1e-5, # default value.
... 'norm' : np.inf, # default value.
... 'eps' : 1.4901161193847656e-08} # default value.
>>> res2 = optimize.minimize(f, x0, jac=gradf, args=args,
... method='CG', options=opts)
Optimization terminated successfully.
Current function value: 1.617021
Iterations: 2
Function evaluations: 5
Gradient evaluations: 5
>>> res2.x # minimum found
array([-1.80851064 -0.25531915])