1 | /*************************************************************************
|
---|
2 | Copyright (c) 2007-2008, Sergey Bochkanov (ALGLIB project).
|
---|
3 |
|
---|
4 | >>> SOURCE LICENSE >>>
|
---|
5 | This program is free software; you can redistribute it and/or modify
|
---|
6 | it under the terms of the GNU General Public License as published by
|
---|
7 | the Free Software Foundation (www.fsf.org); either version 2 of the
|
---|
8 | License, or (at your option) any later version.
|
---|
9 |
|
---|
10 | This program is distributed in the hope that it will be useful,
|
---|
11 | but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
12 | MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
13 | GNU General Public License for more details.
|
---|
14 |
|
---|
15 | A copy of the GNU General Public License is available at
|
---|
16 | http://www.fsf.org/licensing/licenses
|
---|
17 |
|
---|
18 | >>> END OF LICENSE >>>
|
---|
19 | *************************************************************************/
|
---|
20 |
|
---|
21 | using System;
|
---|
22 |
|
---|
23 | namespace alglib
|
---|
24 | {
|
---|
25 | public class minlbfgs
|
---|
26 | {
|
---|
27 | public struct minlbfgsstate
|
---|
28 | {
|
---|
29 | public int n;
|
---|
30 | public int m;
|
---|
31 | public double epsg;
|
---|
32 | public double epsf;
|
---|
33 | public double epsx;
|
---|
34 | public int maxits;
|
---|
35 | public int flags;
|
---|
36 | public bool xrep;
|
---|
37 | public double stpmax;
|
---|
38 | public int nfev;
|
---|
39 | public int mcstage;
|
---|
40 | public int k;
|
---|
41 | public int q;
|
---|
42 | public int p;
|
---|
43 | public double[] rho;
|
---|
44 | public double[,] y;
|
---|
45 | public double[,] s;
|
---|
46 | public double[] theta;
|
---|
47 | public double[] d;
|
---|
48 | public double stp;
|
---|
49 | public double[] work;
|
---|
50 | public double fold;
|
---|
51 | public double gammak;
|
---|
52 | public double[] x;
|
---|
53 | public double f;
|
---|
54 | public double[] g;
|
---|
55 | public bool needfg;
|
---|
56 | public bool xupdated;
|
---|
57 | public AP.rcommstate rstate;
|
---|
58 | public int repiterationscount;
|
---|
59 | public int repnfev;
|
---|
60 | public int repterminationtype;
|
---|
61 | public linmin.linminstate lstate;
|
---|
62 | };
|
---|
63 |
|
---|
64 |
|
---|
65 | public struct minlbfgsreport
|
---|
66 | {
|
---|
67 | public int iterationscount;
|
---|
68 | public int nfev;
|
---|
69 | public int terminationtype;
|
---|
70 | };
|
---|
71 |
|
---|
72 |
|
---|
73 |
|
---|
74 |
|
---|
75 | /*************************************************************************
|
---|
76 | LIMITED MEMORY BFGS METHOD FOR LARGE SCALE OPTIMIZATION
|
---|
77 |
|
---|
78 | The subroutine minimizes function F(x) of N arguments by using a quasi-
|
---|
79 | Newton method (LBFGS scheme) which is optimized to use a minimum amount
|
---|
80 | of memory.
|
---|
81 |
|
---|
82 | The subroutine generates the approximation of an inverse Hessian matrix by
|
---|
83 | using information about the last M steps of the algorithm (instead of N).
|
---|
84 | It lessens a required amount of memory from a value of order N^2 to a
|
---|
85 | value of order 2*N*M.
|
---|
86 |
|
---|
87 | INPUT PARAMETERS:
|
---|
88 | N - problem dimension. N>0
|
---|
89 | M - number of corrections in the BFGS scheme of Hessian
|
---|
90 | approximation update. Recommended value: 3<=M<=7. The smaller
|
---|
91 | value causes worse convergence, the bigger will not cause a
|
---|
92 | considerably better convergence, but will cause a fall in the
|
---|
93 | performance. M<=N.
|
---|
94 | X - initial solution approximation, array[0..N-1].
|
---|
95 |
|
---|
96 | OUTPUT PARAMETERS:
|
---|
97 | State - structure used for reverse communication.
|
---|
98 |
|
---|
99 | This function initializes State structure with default optimization
|
---|
100 | parameters (stopping conditions, step size, etc.). Use MinLBFGSSet??????()
|
---|
101 | functions to tune optimization parameters.
|
---|
102 |
|
---|
103 | After all optimization parameters are tuned, you should use
|
---|
104 | MinLBFGSIteration() function to advance algorithm iterations.
|
---|
105 |
|
---|
106 | NOTES:
|
---|
107 |
|
---|
108 | 1. you may tune stopping conditions with MinLBFGSSetCond() function
|
---|
109 | 2. if target function contains exp() or other fast growing functions, and
|
---|
110 | optimization algorithm makes too large steps which leads to overflow,
|
---|
111 | use MinLBFGSSetStpMax() function to bound algorithm's steps. However,
|
---|
112 | L-BFGS rarely needs such a tuning.
|
---|
113 |
|
---|
114 |
|
---|
115 | -- ALGLIB --
|
---|
116 | Copyright 02.04.2010 by Bochkanov Sergey
|
---|
117 | *************************************************************************/
|
---|
118 | public static void minlbfgscreate(int n,
|
---|
119 | int m,
|
---|
120 | ref double[] x,
|
---|
121 | ref minlbfgsstate state)
|
---|
122 | {
|
---|
123 | minlbfgscreatex(n, m, ref x, 0, ref state);
|
---|
124 | }
|
---|
125 |
|
---|
126 |
|
---|
127 | /*************************************************************************
|
---|
128 | This function sets stopping conditions for L-BFGS optimization algorithm.
|
---|
129 |
|
---|
130 | INPUT PARAMETERS:
|
---|
131 | State - structure which stores algorithm state between calls and
|
---|
132 | which is used for reverse communication. Must be initialized
|
---|
133 | with MinLBFGSCreate()
|
---|
134 | EpsG - >=0
|
---|
135 | The subroutine finishes its work if the condition
|
---|
136 | ||G||<EpsG is satisfied, where ||.|| means Euclidian norm,
|
---|
137 | G - gradient.
|
---|
138 | EpsF - >=0
|
---|
139 | The subroutine finishes its work if on k+1-th iteration
|
---|
140 | the condition |F(k+1)-F(k)|<=EpsF*max{|F(k)|,|F(k+1)|,1}
|
---|
141 | is satisfied.
|
---|
142 | EpsX - >=0
|
---|
143 | The subroutine finishes its work if on k+1-th iteration
|
---|
144 | the condition |X(k+1)-X(k)| <= EpsX is fulfilled.
|
---|
145 | MaxIts - maximum number of iterations. If MaxIts=0, the number of
|
---|
146 | iterations is unlimited.
|
---|
147 |
|
---|
148 | Passing EpsG=0, EpsF=0, EpsX=0 and MaxIts=0 (simultaneously) will lead to
|
---|
149 | automatic stopping criterion selection (small EpsX).
|
---|
150 |
|
---|
151 | -- ALGLIB --
|
---|
152 | Copyright 02.04.2010 by Bochkanov Sergey
|
---|
153 | *************************************************************************/
|
---|
154 | public static void minlbfgssetcond(ref minlbfgsstate state,
|
---|
155 | double epsg,
|
---|
156 | double epsf,
|
---|
157 | double epsx,
|
---|
158 | int maxits)
|
---|
159 | {
|
---|
160 | System.Diagnostics.Debug.Assert((double)(epsg)>=(double)(0), "MinLBFGSSetCond: negative EpsG!");
|
---|
161 | System.Diagnostics.Debug.Assert((double)(epsf)>=(double)(0), "MinLBFGSSetCond: negative EpsF!");
|
---|
162 | System.Diagnostics.Debug.Assert((double)(epsx)>=(double)(0), "MinLBFGSSetCond: negative EpsX!");
|
---|
163 | System.Diagnostics.Debug.Assert(maxits>=0, "MinLBFGSSetCond: negative MaxIts!");
|
---|
164 | if( (double)(epsg)==(double)(0) & (double)(epsf)==(double)(0) & (double)(epsx)==(double)(0) & maxits==0 )
|
---|
165 | {
|
---|
166 | epsx = 1.0E-6;
|
---|
167 | }
|
---|
168 | state.epsg = epsg;
|
---|
169 | state.epsf = epsf;
|
---|
170 | state.epsx = epsx;
|
---|
171 | state.maxits = maxits;
|
---|
172 | }
|
---|
173 |
|
---|
174 |
|
---|
175 | /*************************************************************************
|
---|
176 | This function turns on/off reporting.
|
---|
177 |
|
---|
178 | INPUT PARAMETERS:
|
---|
179 | State - structure which stores algorithm state between calls and
|
---|
180 | which is used for reverse communication. Must be
|
---|
181 | initialized with MinLBFGSCreate()
|
---|
182 | NeedXRep- whether iteration reports are needed or not
|
---|
183 |
|
---|
184 | Usually algorithm returns from MinLBFGSIteration() only when it needs
|
---|
185 | function/gradient/ (which is indicated by NeedFG field. However, with this
|
---|
186 | function we can let it stop after each iteration (one iteration may
|
---|
187 | include more than one function evaluation), which is indicated by XUpdated
|
---|
188 | field.
|
---|
189 |
|
---|
190 |
|
---|
191 | -- ALGLIB --
|
---|
192 | Copyright 02.04.2010 by Bochkanov Sergey
|
---|
193 | *************************************************************************/
|
---|
194 | public static void minlbfgssetxrep(ref minlbfgsstate state,
|
---|
195 | bool needxrep)
|
---|
196 | {
|
---|
197 | state.xrep = needxrep;
|
---|
198 | }
|
---|
199 |
|
---|
200 |
|
---|
201 | /*************************************************************************
|
---|
202 | This function sets maximum step length
|
---|
203 |
|
---|
204 | INPUT PARAMETERS:
|
---|
205 | State - structure which stores algorithm state between calls and
|
---|
206 | which is used for reverse communication. Must be
|
---|
207 | initialized with MinLBFGSCreate()
|
---|
208 | StpMax - maximum step length, >=0. Set StpMax to 0.0, if you don't
|
---|
209 | want to limit step length.
|
---|
210 |
|
---|
211 | Use this subroutine when you optimize target function which contains exp()
|
---|
212 | or other fast growing functions, and optimization algorithm makes too
|
---|
213 | large steps which leads to overflow. This function allows us to reject
|
---|
214 | steps that are too large (and therefore expose us to the possible
|
---|
215 | overflow) without actually calculating function value at the x+stp*d.
|
---|
216 |
|
---|
217 | -- ALGLIB --
|
---|
218 | Copyright 02.04.2010 by Bochkanov Sergey
|
---|
219 | *************************************************************************/
|
---|
220 | public static void minlbfgssetstpmax(ref minlbfgsstate state,
|
---|
221 | double stpmax)
|
---|
222 | {
|
---|
223 | System.Diagnostics.Debug.Assert((double)(stpmax)>=(double)(0), "MinLBFGSSetStpMax: StpMax<0!");
|
---|
224 | state.stpmax = stpmax;
|
---|
225 | }
|
---|
226 |
|
---|
227 |
|
---|
228 | /*************************************************************************
|
---|
229 | Extended subroutine for internal use only.
|
---|
230 |
|
---|
231 | Accepts additional parameters:
|
---|
232 |
|
---|
233 | Flags - additional settings:
|
---|
234 | * Flags = 0 means no additional settings
|
---|
235 | * Flags = 1 "do not allocate memory". used when solving
|
---|
236 | a many subsequent tasks with same N/M values.
|
---|
237 | First call MUST be without this flag bit set,
|
---|
238 | subsequent calls of MinLBFGS with same
|
---|
239 | MinLBFGSState structure can set Flags to 1.
|
---|
240 |
|
---|
241 | -- ALGLIB --
|
---|
242 | Copyright 02.04.2010 by Bochkanov Sergey
|
---|
243 | *************************************************************************/
|
---|
244 | public static void minlbfgscreatex(int n,
|
---|
245 | int m,
|
---|
246 | ref double[] x,
|
---|
247 | int flags,
|
---|
248 | ref minlbfgsstate state)
|
---|
249 | {
|
---|
250 | bool allocatemem = new bool();
|
---|
251 | int i_ = 0;
|
---|
252 |
|
---|
253 | System.Diagnostics.Debug.Assert(n>=1, "MinLBFGS: N too small!");
|
---|
254 | System.Diagnostics.Debug.Assert(m>=1, "MinLBFGS: M too small!");
|
---|
255 | System.Diagnostics.Debug.Assert(m<=n, "MinLBFGS: M too large!");
|
---|
256 |
|
---|
257 | //
|
---|
258 | // Initialize
|
---|
259 | //
|
---|
260 | state.n = n;
|
---|
261 | state.m = m;
|
---|
262 | state.flags = flags;
|
---|
263 | allocatemem = flags%2==0;
|
---|
264 | flags = flags/2;
|
---|
265 | if( allocatemem )
|
---|
266 | {
|
---|
267 | state.rho = new double[m-1+1];
|
---|
268 | state.theta = new double[m-1+1];
|
---|
269 | state.y = new double[m-1+1, n-1+1];
|
---|
270 | state.s = new double[m-1+1, n-1+1];
|
---|
271 | state.d = new double[n-1+1];
|
---|
272 | state.x = new double[n-1+1];
|
---|
273 | state.g = new double[n-1+1];
|
---|
274 | state.work = new double[n-1+1];
|
---|
275 | }
|
---|
276 | minlbfgssetcond(ref state, 0, 0, 0, 0);
|
---|
277 | minlbfgssetxrep(ref state, false);
|
---|
278 | minlbfgssetstpmax(ref state, 0);
|
---|
279 |
|
---|
280 | //
|
---|
281 | // Prepare first run
|
---|
282 | //
|
---|
283 | state.k = 0;
|
---|
284 | for(i_=0; i_<=n-1;i_++)
|
---|
285 | {
|
---|
286 | state.x[i_] = x[i_];
|
---|
287 | }
|
---|
288 | state.rstate.ia = new int[6+1];
|
---|
289 | state.rstate.ra = new double[4+1];
|
---|
290 | state.rstate.stage = -1;
|
---|
291 | }
|
---|
292 |
|
---|
293 |
|
---|
294 | /*************************************************************************
|
---|
295 | L-BFGS iterations
|
---|
296 |
|
---|
297 | Called after initialization with MinLBFGSCreate() function.
|
---|
298 |
|
---|
299 | INPUT PARAMETERS:
|
---|
300 | State - structure which stores algorithm state between calls and
|
---|
301 | which is used for reverse communication. Must be initialized
|
---|
302 | with MinLBFGSCreate()
|
---|
303 |
|
---|
304 | RESULT:
|
---|
305 | * if function returned False, iterative proces has converged.
|
---|
306 | Use MinLBFGSResults() to obtain optimization results.
|
---|
307 | * if subroutine returned True, then, depending on structure fields, we
|
---|
308 | have one of the following situations
|
---|
309 |
|
---|
310 |
|
---|
311 | === FUNC/GRAD REQUEST ===
|
---|
312 | State.NeedFG is True => function value/gradient are needed.
|
---|
313 | Caller should calculate function value State.F and gradient
|
---|
314 | State.G[0..N-1] at State.X[0..N-1] and call MinLBFGSIteration() again.
|
---|
315 |
|
---|
316 | === NEW INTERATION IS REPORTED ===
|
---|
317 | State.XUpdated is True => one more iteration was made.
|
---|
318 | State.X contains current position, State.F contains function value at X.
|
---|
319 | You can read info from these fields, but never modify them because they
|
---|
320 | contain the only copy of optimization algorithm state.
|
---|
321 |
|
---|
322 |
|
---|
323 | One and only one of these fields (NeedFG, XUpdated) is true on return. New
|
---|
324 | iterations are reported only when reports are explicitly turned on by
|
---|
325 | MinLBFGSSetXRep() function, so if you never called it, you can expect that
|
---|
326 | NeedFG is always True.
|
---|
327 |
|
---|
328 |
|
---|
329 | -- ALGLIB --
|
---|
330 | Copyright 20.03.2009 by Bochkanov Sergey
|
---|
331 | *************************************************************************/
|
---|
332 | public static bool minlbfgsiteration(ref minlbfgsstate state)
|
---|
333 | {
|
---|
334 | bool result = new bool();
|
---|
335 | int n = 0;
|
---|
336 | int m = 0;
|
---|
337 | int maxits = 0;
|
---|
338 | double epsf = 0;
|
---|
339 | double epsg = 0;
|
---|
340 | double epsx = 0;
|
---|
341 | int i = 0;
|
---|
342 | int j = 0;
|
---|
343 | int ic = 0;
|
---|
344 | int mcinfo = 0;
|
---|
345 | double v = 0;
|
---|
346 | double vv = 0;
|
---|
347 | int i_ = 0;
|
---|
348 |
|
---|
349 |
|
---|
350 | //
|
---|
351 | // Reverse communication preparations
|
---|
352 | // I know it looks ugly, but it works the same way
|
---|
353 | // anywhere from C++ to Python.
|
---|
354 | //
|
---|
355 | // This code initializes locals by:
|
---|
356 | // * random values determined during code
|
---|
357 | // generation - on first subroutine call
|
---|
358 | // * values from previous call - on subsequent calls
|
---|
359 | //
|
---|
360 | if( state.rstate.stage>=0 )
|
---|
361 | {
|
---|
362 | n = state.rstate.ia[0];
|
---|
363 | m = state.rstate.ia[1];
|
---|
364 | maxits = state.rstate.ia[2];
|
---|
365 | i = state.rstate.ia[3];
|
---|
366 | j = state.rstate.ia[4];
|
---|
367 | ic = state.rstate.ia[5];
|
---|
368 | mcinfo = state.rstate.ia[6];
|
---|
369 | epsf = state.rstate.ra[0];
|
---|
370 | epsg = state.rstate.ra[1];
|
---|
371 | epsx = state.rstate.ra[2];
|
---|
372 | v = state.rstate.ra[3];
|
---|
373 | vv = state.rstate.ra[4];
|
---|
374 | }
|
---|
375 | else
|
---|
376 | {
|
---|
377 | n = -983;
|
---|
378 | m = -989;
|
---|
379 | maxits = -834;
|
---|
380 | i = 900;
|
---|
381 | j = -287;
|
---|
382 | ic = 364;
|
---|
383 | mcinfo = 214;
|
---|
384 | epsf = -338;
|
---|
385 | epsg = -686;
|
---|
386 | epsx = 912;
|
---|
387 | v = 585;
|
---|
388 | vv = 497;
|
---|
389 | }
|
---|
390 | if( state.rstate.stage==0 )
|
---|
391 | {
|
---|
392 | goto lbl_0;
|
---|
393 | }
|
---|
394 | if( state.rstate.stage==1 )
|
---|
395 | {
|
---|
396 | goto lbl_1;
|
---|
397 | }
|
---|
398 | if( state.rstate.stage==2 )
|
---|
399 | {
|
---|
400 | goto lbl_2;
|
---|
401 | }
|
---|
402 | if( state.rstate.stage==3 )
|
---|
403 | {
|
---|
404 | goto lbl_3;
|
---|
405 | }
|
---|
406 |
|
---|
407 | //
|
---|
408 | // Routine body
|
---|
409 | //
|
---|
410 |
|
---|
411 | //
|
---|
412 | // Unload frequently used variables from State structure
|
---|
413 | // (just for typing convinience)
|
---|
414 | //
|
---|
415 | n = state.n;
|
---|
416 | m = state.m;
|
---|
417 | epsg = state.epsg;
|
---|
418 | epsf = state.epsf;
|
---|
419 | epsx = state.epsx;
|
---|
420 | maxits = state.maxits;
|
---|
421 | state.repterminationtype = 0;
|
---|
422 | state.repiterationscount = 0;
|
---|
423 | state.repnfev = 0;
|
---|
424 |
|
---|
425 | //
|
---|
426 | // Calculate F/G at the initial point
|
---|
427 | //
|
---|
428 | clearrequestfields(ref state);
|
---|
429 | state.needfg = true;
|
---|
430 | state.rstate.stage = 0;
|
---|
431 | goto lbl_rcomm;
|
---|
432 | lbl_0:
|
---|
433 | if( ! state.xrep )
|
---|
434 | {
|
---|
435 | goto lbl_4;
|
---|
436 | }
|
---|
437 | clearrequestfields(ref state);
|
---|
438 | state.xupdated = true;
|
---|
439 | state.rstate.stage = 1;
|
---|
440 | goto lbl_rcomm;
|
---|
441 | lbl_1:
|
---|
442 | lbl_4:
|
---|
443 | state.repnfev = 1;
|
---|
444 | state.fold = state.f;
|
---|
445 | v = 0.0;
|
---|
446 | for(i_=0; i_<=n-1;i_++)
|
---|
447 | {
|
---|
448 | v += state.g[i_]*state.g[i_];
|
---|
449 | }
|
---|
450 | v = Math.Sqrt(v);
|
---|
451 | if( (double)(v)<=(double)(epsg) )
|
---|
452 | {
|
---|
453 | state.repterminationtype = 4;
|
---|
454 | result = false;
|
---|
455 | return result;
|
---|
456 | }
|
---|
457 |
|
---|
458 | //
|
---|
459 | // Choose initial step
|
---|
460 | //
|
---|
461 | if( (double)(state.stpmax)==(double)(0) )
|
---|
462 | {
|
---|
463 | state.stp = Math.Min(1.0/v, 1);
|
---|
464 | }
|
---|
465 | else
|
---|
466 | {
|
---|
467 | state.stp = Math.Min(1.0/v, state.stpmax);
|
---|
468 | }
|
---|
469 | for(i_=0; i_<=n-1;i_++)
|
---|
470 | {
|
---|
471 | state.d[i_] = -state.g[i_];
|
---|
472 | }
|
---|
473 |
|
---|
474 | //
|
---|
475 | // Main cycle
|
---|
476 | //
|
---|
477 | lbl_6:
|
---|
478 | if( false )
|
---|
479 | {
|
---|
480 | goto lbl_7;
|
---|
481 | }
|
---|
482 |
|
---|
483 | //
|
---|
484 | // Main cycle: prepare to 1-D line search
|
---|
485 | //
|
---|
486 | state.p = state.k%m;
|
---|
487 | state.q = Math.Min(state.k, m-1);
|
---|
488 |
|
---|
489 | //
|
---|
490 | // Store X[k], G[k]
|
---|
491 | //
|
---|
492 | for(i_=0; i_<=n-1;i_++)
|
---|
493 | {
|
---|
494 | state.s[state.p,i_] = -state.x[i_];
|
---|
495 | }
|
---|
496 | for(i_=0; i_<=n-1;i_++)
|
---|
497 | {
|
---|
498 | state.y[state.p,i_] = -state.g[i_];
|
---|
499 | }
|
---|
500 |
|
---|
501 | //
|
---|
502 | // Minimize F(x+alpha*d)
|
---|
503 | // Calculate S[k], Y[k]
|
---|
504 | //
|
---|
505 | state.mcstage = 0;
|
---|
506 | if( state.k!=0 )
|
---|
507 | {
|
---|
508 | state.stp = 1.0;
|
---|
509 | }
|
---|
510 | linmin.linminnormalized(ref state.d, ref state.stp, n);
|
---|
511 | linmin.mcsrch(n, ref state.x, ref state.f, ref state.g, ref state.d, ref state.stp, state.stpmax, ref mcinfo, ref state.nfev, ref state.work, ref state.lstate, ref state.mcstage);
|
---|
512 | lbl_8:
|
---|
513 | if( state.mcstage==0 )
|
---|
514 | {
|
---|
515 | goto lbl_9;
|
---|
516 | }
|
---|
517 | clearrequestfields(ref state);
|
---|
518 | state.needfg = true;
|
---|
519 | state.rstate.stage = 2;
|
---|
520 | goto lbl_rcomm;
|
---|
521 | lbl_2:
|
---|
522 | linmin.mcsrch(n, ref state.x, ref state.f, ref state.g, ref state.d, ref state.stp, state.stpmax, ref mcinfo, ref state.nfev, ref state.work, ref state.lstate, ref state.mcstage);
|
---|
523 | goto lbl_8;
|
---|
524 | lbl_9:
|
---|
525 | if( ! state.xrep )
|
---|
526 | {
|
---|
527 | goto lbl_10;
|
---|
528 | }
|
---|
529 |
|
---|
530 | //
|
---|
531 | // report
|
---|
532 | //
|
---|
533 | clearrequestfields(ref state);
|
---|
534 | state.xupdated = true;
|
---|
535 | state.rstate.stage = 3;
|
---|
536 | goto lbl_rcomm;
|
---|
537 | lbl_3:
|
---|
538 | lbl_10:
|
---|
539 | state.repnfev = state.repnfev+state.nfev;
|
---|
540 | state.repiterationscount = state.repiterationscount+1;
|
---|
541 | for(i_=0; i_<=n-1;i_++)
|
---|
542 | {
|
---|
543 | state.s[state.p,i_] = state.s[state.p,i_] + state.x[i_];
|
---|
544 | }
|
---|
545 | for(i_=0; i_<=n-1;i_++)
|
---|
546 | {
|
---|
547 | state.y[state.p,i_] = state.y[state.p,i_] + state.g[i_];
|
---|
548 | }
|
---|
549 |
|
---|
550 | //
|
---|
551 | // Stopping conditions
|
---|
552 | //
|
---|
553 | if( state.repiterationscount>=maxits & maxits>0 )
|
---|
554 | {
|
---|
555 |
|
---|
556 | //
|
---|
557 | // Too many iterations
|
---|
558 | //
|
---|
559 | state.repterminationtype = 5;
|
---|
560 | result = false;
|
---|
561 | return result;
|
---|
562 | }
|
---|
563 | v = 0.0;
|
---|
564 | for(i_=0; i_<=n-1;i_++)
|
---|
565 | {
|
---|
566 | v += state.g[i_]*state.g[i_];
|
---|
567 | }
|
---|
568 | if( (double)(Math.Sqrt(v))<=(double)(epsg) )
|
---|
569 | {
|
---|
570 |
|
---|
571 | //
|
---|
572 | // Gradient is small enough
|
---|
573 | //
|
---|
574 | state.repterminationtype = 4;
|
---|
575 | result = false;
|
---|
576 | return result;
|
---|
577 | }
|
---|
578 | if( (double)(state.fold-state.f)<=(double)(epsf*Math.Max(Math.Abs(state.fold), Math.Max(Math.Abs(state.f), 1.0))) )
|
---|
579 | {
|
---|
580 |
|
---|
581 | //
|
---|
582 | // F(k+1)-F(k) is small enough
|
---|
583 | //
|
---|
584 | state.repterminationtype = 1;
|
---|
585 | result = false;
|
---|
586 | return result;
|
---|
587 | }
|
---|
588 | v = 0.0;
|
---|
589 | for(i_=0; i_<=n-1;i_++)
|
---|
590 | {
|
---|
591 | v += state.s[state.p,i_]*state.s[state.p,i_];
|
---|
592 | }
|
---|
593 | if( (double)(Math.Sqrt(v))<=(double)(epsx) )
|
---|
594 | {
|
---|
595 |
|
---|
596 | //
|
---|
597 | // X(k+1)-X(k) is small enough
|
---|
598 | //
|
---|
599 | state.repterminationtype = 2;
|
---|
600 | result = false;
|
---|
601 | return result;
|
---|
602 | }
|
---|
603 |
|
---|
604 | //
|
---|
605 | // If Wolfe conditions are satisfied, we can update
|
---|
606 | // limited memory model.
|
---|
607 | //
|
---|
608 | // However, if conditions are not satisfied (NFEV limit is met,
|
---|
609 | // function is too wild, ...), we'll skip L-BFGS update
|
---|
610 | //
|
---|
611 | if( mcinfo!=1 )
|
---|
612 | {
|
---|
613 |
|
---|
614 | //
|
---|
615 | // Skip update.
|
---|
616 | //
|
---|
617 | // In such cases we'll initialize search direction by
|
---|
618 | // antigradient vector, because it leads to more
|
---|
619 | // transparent code with less number of special cases
|
---|
620 | //
|
---|
621 | state.fold = state.f;
|
---|
622 | for(i_=0; i_<=n-1;i_++)
|
---|
623 | {
|
---|
624 | state.d[i_] = -state.g[i_];
|
---|
625 | }
|
---|
626 | }
|
---|
627 | else
|
---|
628 | {
|
---|
629 |
|
---|
630 | //
|
---|
631 | // Calculate Rho[k], GammaK
|
---|
632 | //
|
---|
633 | v = 0.0;
|
---|
634 | for(i_=0; i_<=n-1;i_++)
|
---|
635 | {
|
---|
636 | v += state.y[state.p,i_]*state.s[state.p,i_];
|
---|
637 | }
|
---|
638 | vv = 0.0;
|
---|
639 | for(i_=0; i_<=n-1;i_++)
|
---|
640 | {
|
---|
641 | vv += state.y[state.p,i_]*state.y[state.p,i_];
|
---|
642 | }
|
---|
643 | if( (double)(v)==(double)(0) | (double)(vv)==(double)(0) )
|
---|
644 | {
|
---|
645 |
|
---|
646 | //
|
---|
647 | // Rounding errors make further iterations impossible.
|
---|
648 | //
|
---|
649 | state.repterminationtype = -2;
|
---|
650 | result = false;
|
---|
651 | return result;
|
---|
652 | }
|
---|
653 | state.rho[state.p] = 1/v;
|
---|
654 | state.gammak = v/vv;
|
---|
655 |
|
---|
656 | //
|
---|
657 | // Calculate d(k+1) = -H(k+1)*g(k+1)
|
---|
658 | //
|
---|
659 | // for I:=K downto K-Q do
|
---|
660 | // V = s(i)^T * work(iteration:I)
|
---|
661 | // theta(i) = V
|
---|
662 | // work(iteration:I+1) = work(iteration:I) - V*Rho(i)*y(i)
|
---|
663 | // work(last iteration) = H0*work(last iteration)
|
---|
664 | // for I:=K-Q to K do
|
---|
665 | // V = y(i)^T*work(iteration:I)
|
---|
666 | // work(iteration:I+1) = work(iteration:I) +(-V+theta(i))*Rho(i)*s(i)
|
---|
667 | //
|
---|
668 | // NOW WORK CONTAINS d(k+1)
|
---|
669 | //
|
---|
670 | for(i_=0; i_<=n-1;i_++)
|
---|
671 | {
|
---|
672 | state.work[i_] = state.g[i_];
|
---|
673 | }
|
---|
674 | for(i=state.k; i>=state.k-state.q; i--)
|
---|
675 | {
|
---|
676 | ic = i%m;
|
---|
677 | v = 0.0;
|
---|
678 | for(i_=0; i_<=n-1;i_++)
|
---|
679 | {
|
---|
680 | v += state.s[ic,i_]*state.work[i_];
|
---|
681 | }
|
---|
682 | state.theta[ic] = v;
|
---|
683 | vv = v*state.rho[ic];
|
---|
684 | for(i_=0; i_<=n-1;i_++)
|
---|
685 | {
|
---|
686 | state.work[i_] = state.work[i_] - vv*state.y[ic,i_];
|
---|
687 | }
|
---|
688 | }
|
---|
689 | v = state.gammak;
|
---|
690 | for(i_=0; i_<=n-1;i_++)
|
---|
691 | {
|
---|
692 | state.work[i_] = v*state.work[i_];
|
---|
693 | }
|
---|
694 | for(i=state.k-state.q; i<=state.k; i++)
|
---|
695 | {
|
---|
696 | ic = i%m;
|
---|
697 | v = 0.0;
|
---|
698 | for(i_=0; i_<=n-1;i_++)
|
---|
699 | {
|
---|
700 | v += state.y[ic,i_]*state.work[i_];
|
---|
701 | }
|
---|
702 | vv = state.rho[ic]*(-v+state.theta[ic]);
|
---|
703 | for(i_=0; i_<=n-1;i_++)
|
---|
704 | {
|
---|
705 | state.work[i_] = state.work[i_] + vv*state.s[ic,i_];
|
---|
706 | }
|
---|
707 | }
|
---|
708 | for(i_=0; i_<=n-1;i_++)
|
---|
709 | {
|
---|
710 | state.d[i_] = -state.work[i_];
|
---|
711 | }
|
---|
712 |
|
---|
713 | //
|
---|
714 | // Next step
|
---|
715 | //
|
---|
716 | state.fold = state.f;
|
---|
717 | state.k = state.k+1;
|
---|
718 | }
|
---|
719 | goto lbl_6;
|
---|
720 | lbl_7:
|
---|
721 | result = false;
|
---|
722 | return result;
|
---|
723 |
|
---|
724 | //
|
---|
725 | // Saving state
|
---|
726 | //
|
---|
727 | lbl_rcomm:
|
---|
728 | result = true;
|
---|
729 | state.rstate.ia[0] = n;
|
---|
730 | state.rstate.ia[1] = m;
|
---|
731 | state.rstate.ia[2] = maxits;
|
---|
732 | state.rstate.ia[3] = i;
|
---|
733 | state.rstate.ia[4] = j;
|
---|
734 | state.rstate.ia[5] = ic;
|
---|
735 | state.rstate.ia[6] = mcinfo;
|
---|
736 | state.rstate.ra[0] = epsf;
|
---|
737 | state.rstate.ra[1] = epsg;
|
---|
738 | state.rstate.ra[2] = epsx;
|
---|
739 | state.rstate.ra[3] = v;
|
---|
740 | state.rstate.ra[4] = vv;
|
---|
741 | return result;
|
---|
742 | }
|
---|
743 |
|
---|
744 |
|
---|
745 | /*************************************************************************
|
---|
746 | L-BFGS algorithm results
|
---|
747 |
|
---|
748 | Called after MinLBFGSIteration() returned False.
|
---|
749 |
|
---|
750 | INPUT PARAMETERS:
|
---|
751 | State - algorithm state (used by MinLBFGSIteration).
|
---|
752 |
|
---|
753 | OUTPUT PARAMETERS:
|
---|
754 | X - array[0..N-1], solution
|
---|
755 | Rep - optimization report:
|
---|
756 | * Rep.TerminationType completetion code:
|
---|
757 | * -2 rounding errors prevent further improvement.
|
---|
758 | X contains best point found.
|
---|
759 | * -1 incorrect parameters were specified
|
---|
760 | * 1 relative function improvement is no more than
|
---|
761 | EpsF.
|
---|
762 | * 2 relative step is no more than EpsX.
|
---|
763 | * 4 gradient norm is no more than EpsG
|
---|
764 | * 5 MaxIts steps was taken
|
---|
765 | * 7 stopping conditions are too stringent,
|
---|
766 | further improvement is impossible
|
---|
767 | * Rep.IterationsCount contains iterations count
|
---|
768 | * NFEV countains number of function calculations
|
---|
769 |
|
---|
770 | -- ALGLIB --
|
---|
771 | Copyright 02.04.2010 by Bochkanov Sergey
|
---|
772 | *************************************************************************/
|
---|
773 | public static void minlbfgsresults(ref minlbfgsstate state,
|
---|
774 | ref double[] x,
|
---|
775 | ref minlbfgsreport rep)
|
---|
776 | {
|
---|
777 | int i_ = 0;
|
---|
778 |
|
---|
779 | x = new double[state.n-1+1];
|
---|
780 | for(i_=0; i_<=state.n-1;i_++)
|
---|
781 | {
|
---|
782 | x[i_] = state.x[i_];
|
---|
783 | }
|
---|
784 | rep.iterationscount = state.repiterationscount;
|
---|
785 | rep.nfev = state.repnfev;
|
---|
786 | rep.terminationtype = state.repterminationtype;
|
---|
787 | }
|
---|
788 |
|
---|
789 |
|
---|
790 | /*************************************************************************
|
---|
791 | Clears request fileds (to be sure that we don't forgot to clear something)
|
---|
792 | *************************************************************************/
|
---|
793 | private static void clearrequestfields(ref minlbfgsstate state)
|
---|
794 | {
|
---|
795 | state.needfg = false;
|
---|
796 | state.xupdated = false;
|
---|
797 | }
|
---|
798 | }
|
---|
799 | }
|
---|