File: lbfgsb.py

package info (click to toggle)
python-scipy 0.6.0-12
  • links: PTS, VCS
  • area: main
  • in suites: lenny
  • size: 32,016 kB
  • ctags: 46,675
  • sloc: cpp: 124,854; ansic: 110,614; python: 108,664; fortran: 76,260; objc: 424; makefile: 384; sh: 10
file content (271 lines) | stat: -rwxr-xr-x 9,104 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
## Automatically adapted for scipy Oct 07, 2005 by convertcode.py


## License for the Python wrapper
## ==============================

## Copyright (c) 2004 David M. Cooke <cookedm@physics.mcmaster.ca>

## Permission is hereby granted, free of charge, to any person obtaining a copy of
## this software and associated documentation files (the "Software"), to deal in
## the Software without restriction, including without limitation the rights to
## use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies
## of the Software, and to permit persons to whom the Software is furnished to do
## so, subject to the following conditions:

## The above copyright notice and this permission notice shall be included in all
## copies or substantial portions of the Software.

## THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
## IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
## FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
## AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
## LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
## OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
## SOFTWARE.

## Modifications by Travis Oliphant and Enthought, Inc.  for inclusion in SciPy

from numpy import zeros, float64, array, int32
import _lbfgsb
import optimize

approx_fprime = optimize.approx_fprime

def fmin_l_bfgs_b(func, x0, fprime=None, args=(),
                  approx_grad=0,
                  bounds=None, m=10, factr=1e7, pgtol=1e-5,
                  epsilon=1e-8,
                  iprint=-1, maxfun=15000):
    """
    Minimize a function func using the L-BFGS-B algorithm.

    Arguments:

    func    -- function to minimize. Called as func(x, *args)

    x0      -- initial guess to minimum

    fprime  -- gradient of func. If None, then func returns the function
               value and the gradient ( f, g = func(x, *args) ), unless
               approx_grad is True then func returns only f.
               Called as fprime(x, *args)

    args    -- arguments to pass to function

    approx_grad -- if true, approximate the gradient numerically and func returns
                   only function value.

    bounds  -- a list of (min, max) pairs for each element in x, defining
               the bounds on that parameter. Use None for one of min or max
               when there is no bound in that direction

    m       -- the maximum number of variable metric corrections
               used to define the limited memory matrix. (the limited memory BFGS
               method does not store the full hessian but uses this many terms in an
               approximation to it).

    factr   -- The iteration stops when
               (f^k - f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= factr*epsmch

               where epsmch is the machine precision, which is automatically
               generated by the code. Typical values for factr: 1e12 for
               low accuracy; 1e7 for moderate accuracy; 10.0 for extremely
               high accuracy.

    pgtol   -- The iteration will stop when
                  max{|proj g_i | i = 1, ..., n} <= pgtol
               where pg_i is the ith component of the projected gradient.

    epsilon -- step size used when approx_grad is true, for numerically
               calculating the gradient

    iprint  -- controls the frequency of output. <0 means no output.

    maxfun  -- maximum number of function evaluations.


    Returns:
    x, f, d = fmin_lbfgs_b(func, x0, ...)

    x -- position of the minimum
    f -- value of func at the minimum
    d -- dictionary of information from routine
        d['warnflag'] is
            0 if converged,
            1 if too many function evaluations,
            2 if stopped for another reason, given in d['task']
        d['grad'] is the gradient at the minimum (should be 0 ish)
        d['funcalls'] is the number of function calls made.


   License of L-BFGS-B (Fortran code)
   ==================================

   The version included here (in fortran code) is 2.1 (released in 1997). It was
   written by Ciyou Zhu, Richard Byrd, and Jorge Nocedal <nocedal@ece.nwu.edu>. It
   carries the following condition for use:

   This software is freely available, but we expect that all publications
   describing  work using this software , or all commercial products using it,
   quote at least one of the references given below.

   References
     * R. H. Byrd, P. Lu and J. Nocedal. A Limited Memory Algorithm for Bound
       Constrained Optimization, (1995), SIAM Journal on Scientific and
       Statistical Computing , 16, 5, pp. 1190-1208.
     * C. Zhu, R. H. Byrd and J. Nocedal. L-BFGS-B: Algorithm 778: L-BFGS-B,
       FORTRAN routines for large scale bound constrained optimization (1997),
       ACM Transactions on Mathematical Software, Vol 23, Num. 4, pp. 550 - 560.

    See also:
        scikits.openopt, which offers a unified syntax to call this and other solvers

        fmin, fmin_powell, fmin_cg,
               fmin_bfgs, fmin_ncg -- multivariate local optimizers
        leastsq -- nonlinear least squares minimizer

        fmin_l_bfgs_b, fmin_tnc,
               fmin_cobyla -- constrained multivariate optimizers

        anneal, brute -- global optimizers

        fminbound, brent, golden, bracket -- local scalar minimizers

        fsolve -- n-dimenstional root-finding

        brentq, brenth, ridder, bisect, newton -- one-dimensional root-finding

        fixed_point -- scalar fixed-point finder

    """
    n = len(x0)

    if bounds is None:
        bounds = [(None,None)] * n
    if len(bounds) != n:
        raise ValueError('length of x0 != length of bounds')

    if approx_grad:
        def func_and_grad(x):
            f = func(x, *args)
            g = approx_fprime(x, func, epsilon, *args)
            return f, g
    elif fprime is None:
        def func_and_grad(x):
            f, g = func(x, *args)
            return f, g
    else:
        def func_and_grad(x):
            f = func(x, *args)
            g = fprime(x, *args)
            return f, g

    nbd = zeros((n,), int32)
    low_bnd = zeros((n,), float64)
    upper_bnd = zeros((n,), float64)
    bounds_map = {(None, None): 0,
              (1, None) : 1,
              (1, 1) : 2,
              (None, 1) : 3}
    for i in range(0, n):
        l,u = bounds[i]
        if l is not None:
            low_bnd[i] = l
            l = 1
        if u is not None:
            upper_bnd[i] = u
            u = 1
        nbd[i] = bounds_map[l, u]

    x = array(x0, float64)
    f = array(0.0, float64)
    g = zeros((n,), float64)
    wa = zeros((2*m*n+4*n + 12*m**2 + 12*m,), float64)
    iwa = zeros((3*n,), int32)
    task = zeros(1, 'S60')
    csave = zeros(1,'S60')
    lsave = zeros((4,), int32)
    isave = zeros((44,), int32)
    dsave = zeros((29,), float64)

    task[:] = 'START'

    n_function_evals = 0
    while 1:
#        x, f, g, wa, iwa, task, csave, lsave, isave, dsave = \
        _lbfgsb.setulb(m, x, low_bnd, upper_bnd, nbd, f, g, factr,
                       pgtol, wa, iwa, task, iprint, csave, lsave,
                       isave, dsave)
        task_str = task.tostring()
        if task_str.startswith('FG'):
            # minimization routine wants f and g at the current x
            n_function_evals += 1
            # Overwrite f and g:
            f, g = func_and_grad(x)
        elif task_str.startswith('NEW_X'):
            # new iteration
            if n_function_evals > maxfun:
                task[:] = 'STOP: TOTAL NO. of f AND g EVALUATIONS EXCEEDS LIMIT'
        else:
            break

    task_str = task.tostring().strip('\x00').strip()
    if task_str.startswith('CONV'):
        warnflag = 0
    elif n_function_evals > maxfun:
        warnflag = 1
    else:
        warnflag = 2


    d = {'grad' : g,
         'task' : task_str,
         'funcalls' : n_function_evals,
         'warnflag' : warnflag
        }
    return x, f, d

if __name__ == '__main__':
    def func(x):
        f = 0.25*(x[0]-1)**2
        for i in range(1, x.shape[0]):
            f += (x[i] - x[i-1]**2)**2
        f *= 4
        return f
    def grad(x):
        g = zeros(x.shape, float64)
        t1 = x[1] - x[0]**2
        g[0] = 2*(x[0]-1) - 16*x[0]*t1
        for i in range(1, g.shape[0]-1):
            t2 = t1
            t1 = x[i+1] - x[i]**2
            g[i] = 8*t2 - 16*x[i]*t1
        g[-1] = 8*t1
        return g

    factr = 1e7
    pgtol = 1e-5

    n=25
    m=10

    bounds = [(None,None)] * n
    for i in range(0, n, 2):
        bounds[i] = (1.0, 100)
    for i in range(1, n, 2):
        bounds[i] = (-100, 100)

    x0 = zeros((n,), float64)
    x0[:] = 3

    x, f, d = fmin_l_bfgs_b(func, x0, fprime=grad, m=m,
                            factr=factr, pgtol=pgtol)
    print x
    print f
    print d
    x, f, d = fmin_l_bfgs_b(func, x0, approx_grad=1,
                            m=m, factr=factr, pgtol=pgtol)
    print x
    print f
    print d