Log In | Get Help   
Home My Page Projects Code Snippets Project Openings Mareframe
Summary Activity Forums Tracker Lists Tasks Docs Surveys News SCM Files
[mareframe] Annotation of /trunk/gadget/simann.cc
[mareframe] / trunk / gadget / simann.cc Repository:
ViewVC logotype

Annotation of /trunk/gadget/simann.cc

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1 - (view) (download)

1 : agomez 1 /* ABSTRACT: */
2 :     /* Simulated annealing is a global optimization method that distinguishes */
3 :     /* different local optima. Starting from an initial point, the algorithm */
4 :     /* takes a step and the function is evaluated. When minimizing a function,*/
5 :     /* any downhill step is accepted and the process repeats from this new */
6 :     /* point. An uphill step may be accepted (thus, it can escape from local */
7 :     /* optima). This uphill decision is made by the Metropolis criteria. As */
8 :     /* the optimization process proceeds, the length of the steps decline and */
9 :     /* the algorithm closes in on the global optimum. Since the algorithm */
10 :     /* makes very few assumptions regarding the function to be optimized, it */
11 :     /* is quite robust with respect to non-quadratic surfaces. The degree of */
12 :     /* robustness can be adjusted by the user. In fact, simulated annealing */
13 :     /* can be used as a local optimizer for difficult functions. */
14 :     /* */
15 :     /* The author can be contacted at h2zr1001@vm.cis.smu.edu */
16 :    
17 :     /* This file is a translation of a fortran code, which is an example of the*/
18 :     /* Corana et al. simulated annealing algorithm for multimodal and robust */
19 :     /* optimization as implemented and modified by by Goffe et al. */
20 :     /* */
21 :     /* Use the sample function from Judge with the following suggestions */
22 :     /* to get a feel for how SA works. When you've done this, you should be */
23 :     /* ready to use it on most any function with a fair amount of expertise. */
24 :     /* 1. Run the program as is to make sure it runs okay. Take a look at */
25 :     /* the intermediate output and see how it optimizes as temperature */
26 :     /* (T) falls. Notice how the optimal point is reached and how */
27 :     /* falling T reduces VM. */
28 :     /* 2. Look through the documentation to SA so the following makes a */
29 :     /* bit of sense. In line with the paper, it shouldn't be that hard */
30 :     /* to figure out. The core of the algorithm is described on pp. 4-6 */
31 :     /* and on pp. 28. Also see Corana et al. pp. 264-9. */
32 :     /* 3. To see the importance of different temperatures, try starting */
33 :     /* with a very low one (say T = 10E-5). You'll see (i) it never */
34 :     /* escapes from the local optima (in annealing terminology, it */
35 :     /* quenches) & (ii) the step length (VM) will be quite small. This */
36 :     /* is a key part of the algorithm: as temperature (T) falls, step */
37 :     /* length does too. In a minor point here, note how VM is quickly */
38 :     /* reset from its initial value. Thus, the input VM is not very */
39 :     /* important. This is all the more reason to examine VM once the */
40 :     /* algorithm is underway. */
41 :     /* 4. To see the effect of different parameters and their effect on */
42 :     /* the speed of the algorithm, try RT = .95 & RT = .1. Notice the */
43 :     /* vastly different speed for optimization. Also try NT = 20. Note */
44 :     /* that this sample function is quite easy to optimize, so it will */
45 :     /* tolerate big changes in these parameters. RT and NT are the */
46 :     /* parameters one should adjust to modify the runtime of the */
47 :     /* algorithm and its robustness. */
48 :     /* 5. Try constraining the algorithm with either LB or UB. */
49 :    
50 :     /* Synopsis: */
51 :     /* This routine implements the continuous simulated annealing global */
52 :     /* optimization algorithm described in Corana et al.'s article */
53 :     /* "Minimizing Multimodal Functions of Continuous Variables with the */
54 :     /* "Simulated Annealing" Algorithm" in the September 1987 (vol. 13, */
55 :     /* no. 3, pp. 262-280) issue of the ACM Transactions on Mathematical */
56 :     /* Software. */
57 :    
58 :     /* A very quick (perhaps too quick) overview of SA: */
59 :     /* SA tries to find the global optimum of an N dimensional function. */
60 :     /* It moves both up and downhill and as the optimization process */
61 :     /* proceeds, it focuses on the most promising area. */
62 :     /* To start, it randomly chooses a trial point within the step length */
63 :     /* VM (a vector of length N) of the user selected starting point. The */
64 :     /* function is evaluated at this trial point and its value is compared */
65 :     /* to its value at the initial point. */
66 :     /* In a maximization problem, all uphill moves are accepted and the */
67 :     /* algorithm continues from that trial point. Downhill moves may be */
68 :     /* accepted; the decision is made by the Metropolis criteria. It uses T */
69 :     /* (temperature) and the size of the downhill move in a probabilistic */
70 :     /* manner. The smaller T and the size of the downhill move are, the more */
71 :     /* likely that move will be accepted. If the trial is accepted, the */
72 :     /* algorithm moves on from that point. If it is rejected, another point */
73 :     /* is chosen instead for a trial evaluation. */
74 :     /* Each element of VM periodically adjusted so that half of all */
75 :     /* function evaluations in that direction are accepted. */
76 :     /* A fall in T is imposed upon the system with the RT variable by */
77 :     /* T(i+1) = RT*T(i) where i is the ith iteration. Thus, as T declines, */
78 :     /* downhill moves are less likely to be accepted and the percentage of */
79 :     /* rejections rise. Given the scheme for the selection for VM, VM falls. */
80 :     /* Thus, as T declines, VM falls and SA focuses upon the most promising */
81 :     /* area for optimization. */
82 :    
83 :     /* The importance of the parameter T: */
84 :     /* The parameter T is crucial in using SA successfully. It influences */
85 :     /* VM, the step length over which the algorithm searches for optima. For */
86 :     /* a small intial T, the step length may be too small; thus not enough */
87 :     /* of the function might be evaluated to find the global optima. The user */
88 :     /* should carefully examine VM in the intermediate output (set IPRINT = */
89 :     /* 1) to make sure that VM is appropriate. The relationship between the */
90 :     /* initial temperature and the resulting step length is function */
91 :     /* dependent. */
92 :     /* To determine the starting temperature that is consistent with */
93 :     /* optimizing a function, it is worthwhile to run a trial run first. Set */
94 :     /* RT = 1.5 and T = 1.0. With RT > 1.0, the temperature increases and VM */
95 :     /* rises as well. Then select the T that produces a large enough VM. */
96 :    
97 :     /* For modifications to the algorithm and many details on its use, */
98 :     /* (particularly for econometric applications) see Goffe, Ferrier */
99 :     /* and Rogers, "Global Optimization of Statistical Functions with */
100 :     /* the Simulated Annealing," Journal of Econometrics (forthcoming) */
101 :     /* For a pre-publication copy, contact */
102 :     /* Bill Goffe */
103 :     /* Department of Economics */
104 :     /* Southern Methodist University */
105 :     /* Dallas, TX 75275 */
106 :     /* h2zr1001 @ smuvm1 (Bitnet) */
107 :     /* h2zr1001 @ vm.cis.smu.edu (Internet) */
108 :    
109 :     /* As far as possible, the parameters here have the same name as in */
110 :     /* the description of the algorithm on pp. 266-8 of Corana et al. */
111 :    
112 :     /* Input Parameters: */
113 :     /* Note: The suggested values generally come from Corana et al. To */
114 :     /* drastically reduce runtime, see Goffe et al., pp. 17-8 for */
115 :     /* suggestions on choosing the appropriate RT and NT. */
116 :     /* n - Number of variables in the function to be optimized. (INT) */
117 :     /* x - The starting values for the variables of the function to be */
118 :     /* optimized. (DP(N)) */
119 :     /* max - Denotes whether the function should be maximized or */
120 :     /* minimized. A true value denotes maximization while a false */
121 :     /* value denotes minimization. */
122 :     /* RT - The temperature reduction factor. The value suggested by */
123 :     /* Corana et al. is .85. See Goffe et al. for more advice. (DP) */
124 :     /* EPS - Error tolerance for termination. If the final function */
125 :     /* values from the last neps temperatures differ from the */
126 :     /* corresponding value at the current temperature by less than */
127 :     /* EPS and the final function value at the current temperature */
128 :     /* differs from the current optimal function value by less than */
129 :     /* EPS, execution terminates and IER = 0 is returned. (EP) */
130 :     /* NS - Number of cycles. After NS*N function evaluations, each element */
131 :     /* of VM is adjusted so that approximately half of all function */
132 :     /* evaluations are accepted. The suggested value is 20. (INT) */
133 :     /* nt - Number of iterations before temperature reduction. After */
134 :     /* NT*NS*N function evaluations, temperature (T) is changed */
135 :     /* by the factor RT. Value suggested by Corana et al. is */
136 :     /* MAX(100, 5*N). See Goffe et al. for further advice. (INT) */
137 :     /* NEPS - Number of final function values used to decide upon termi- */
138 :     /* nation. See EPS. Suggested value is 4. (INT) */
139 :     /* maxevl - The maximum number of function evaluations. If it is */
140 :     /* exceeded, IER = 1. (INT) */
141 :     /* lb - The lower bound for the allowable solution variables. (DP(N)) */
142 :     /* ub - The upper bound for the allowable solution variables. (DP(N)) */
143 :     /* If the algorithm chooses X(I) .LT. LB(I) or X(I) .GT. UB(I), */
144 :     /* I = 1, N, a point is from inside is randomly selected. This */
145 :     /* This focuses the algorithm on the region inside UB and LB. */
146 :     /* Unless the user wishes to concentrate the search to a par- */
147 :     /* ticular region, UB and LB should be set to very large positive */
148 :     /* and negative values, respectively. Note that the starting */
149 :     /* vector X should be inside this region. Also note that LB and */
150 :     /* UB are fixed in position, while VM is centered on the last */
151 :     /* accepted trial set of variables that optimizes the function. */
152 :     /* c - Vector that controls the step length adjustment. The suggested */
153 :     /* value for all elements is 2.0. (DP(N)) */
154 :     /* t - On input, the initial temperature. See Goffe et al. for advice. */
155 :     /* On output, the final temperature. (DP) */
156 :     /* vm - The step length vector. On input it should encompass the */
157 :     /* region of interest given the starting value X. For point */
158 :     /* X(I), the next trial point is selected is from X(I) - VM(I) */
159 :     /* to X(I) + VM(I). Since VM is adjusted so that about half */
160 :     /* of all points are accepted, the input value is not very */
161 :     /* important (i.e. is the value is off, SA adjusts VM to the */
162 :     /* correct value). (DP(N)) */
163 :    
164 :     /* Output Parameters: */
165 :     /* xopt - The variables that optimize the function. (DP(N)) */
166 :     /* fopt - The optimal value of the function. (DP) */
167 :    
168 :     /* JMB this has been modified to work with the gadget object structure */
169 :     /* This means that the function has been replaced by a call to ecosystem */
170 :     /* object, and we can use the vector objects that have been defined */
171 :    
172 :     #include "gadget.h"
173 :     #include "optinfo.h"
174 :     #include "mathfunc.h"
175 :     #include "doublevector.h"
176 :     #include "intvector.h"
177 :     #include "errorhandler.h"
178 :     #include "ecosystem.h"
179 :     #include "global.h"
180 :    
181 :     extern Ecosystem* EcoSystem;
182 :    
183 :    
184 :     void OptInfoSimann::OptimiseLikelihood() {
185 :    
186 :     //set initial values
187 :     int nacc = 0; //The number of accepted function evaluations
188 :     int nrej = 0; //The number of rejected function evaluations
189 :     int naccmet = 0; //The number of metropolis accepted function evaluations
190 :    
191 :     double tmp, p, pp, ratio, nsdiv;
192 :     double fopt, funcval, trialf;
193 :     int a, i, j, k, l, offset, quit;
194 :     int rchange, rcheck, rnumber; //Used to randomise the order of the parameters
195 :    
196 :     handle.logMessage(LOGINFO, "\nStarting Simulated Annealing optimisation algorithm\n");
197 :     int nvars = EcoSystem->numOptVariables();
198 :     DoubleVector x(nvars);
199 :     DoubleVector init(nvars);
200 :     DoubleVector trialx(nvars, 0.0);
201 :     DoubleVector bestx(nvars);
202 :     DoubleVector scalex(nvars);
203 :     DoubleVector lowerb(nvars);
204 :     DoubleVector upperb(nvars);
205 :     DoubleVector fstar(tempcheck);
206 :     DoubleVector vm(nvars, vminit);
207 :     IntVector param(nvars, 0);
208 :     IntVector nacp(nvars, 0);
209 :    
210 :     EcoSystem->resetVariables(); //JMB need to reset variables in case they have been scaled
211 :     if (scale)
212 :     EcoSystem->scaleVariables();
213 :     EcoSystem->getOptScaledValues(x);
214 :     EcoSystem->getOptLowerBounds(lowerb);
215 :     EcoSystem->getOptUpperBounds(upperb);
216 :     EcoSystem->getOptInitialValues(init);
217 :    
218 :     for (i = 0; i < nvars; i++) {
219 :     bestx[i] = x[i];
220 :     param[i] = i;
221 :     }
222 :    
223 :     if (scale) {
224 :     for (i = 0; i < nvars; i++) {
225 :     scalex[i] = x[i];
226 :     // Scaling the bounds, because the parameters are scaled
227 :     lowerb[i] = lowerb[i] / init[i];
228 :     upperb[i] = upperb[i] / init[i];
229 :     if (lowerb[i] > upperb[i]) {
230 :     tmp = lowerb[i];
231 :     lowerb[i] = upperb[i];
232 :     upperb[i] = tmp;
233 :     }
234 :     }
235 :     }
236 :    
237 :     //funcval is the function value at x
238 :     funcval = EcoSystem->SimulateAndUpdate(x);
239 :     if (funcval != funcval) { //check for NaN
240 :     handle.logMessage(LOGINFO, "Error starting Simulated Annealing optimisation with f(x) = infinity");
241 :     converge = -1;
242 :     iters = 1;
243 :     return;
244 :     }
245 :    
246 :     //the function is to be minimised so switch the sign of funcval (and trialf)
247 :     funcval = -funcval;
248 :     offset = EcoSystem->getFuncEval(); //number of function evaluations done before loop
249 :     nacc++;
250 :     cs /= lratio; //JMB save processing time
251 :     nsdiv = 1.0 / ns;
252 :     fopt = funcval;
253 :     for (i = 0; i < tempcheck; i++)
254 :     fstar[i] = funcval;
255 :    
256 :     //Start the main loop. Note that it terminates if
257 :     //(i) the algorithm succesfully optimises the function or
258 :     //(ii) there are too many function evaluations
259 :     while (1) {
260 :     for (a = 0; a < nt; a++) {
261 :     //Randomize the order of the parameters once in a while, to avoid
262 :     //the order having an influence on which changes are accepted
263 :     rchange = 0;
264 :     while (rchange < nvars) {
265 :     rnumber = rand() % nvars;
266 :     rcheck = 1;
267 :     for (i = 0; i < rchange; i++)
268 :     if (param[i] == rnumber)
269 :     rcheck = 0;
270 :     if (rcheck) {
271 :     param[rchange] = rnumber;
272 :     rchange++;
273 :     }
274 :     }
275 :    
276 :     for (j = 0; j < ns; j++) {
277 :     for (l = 0; l < nvars; l++) {
278 :     //Generate trialx, the trial value of x
279 :     for (i = 0; i < nvars; i++) {
280 :     if (i == param[l]) {
281 :     trialx[i] = x[i] + ((randomNumber() * 2.0) - 1.0) * vm[i];
282 :    
283 :     //If trialx is out of bounds, try again until we find a point that is OK
284 :     if ((trialx[i] < lowerb[i]) || (trialx[i] > upperb[i])) {
285 :     //JMB - this used to just select a random point between the bounds
286 :     k = 0;
287 :     while ((trialx[i] < lowerb[i]) || (trialx[i] > upperb[i])) {
288 :     trialx[i] = x[i] + ((randomNumber() * 2.0) - 1.0) * vm[i];
289 :     k++;
290 :     if (k > 10) //we've had 10 tries to find a point neatly, so give up
291 :     trialx[i] = lowerb[i] + (upperb[i] - lowerb[i]) * randomNumber();
292 :     }
293 :     }
294 :    
295 :     } else
296 :     trialx[i] = x[i];
297 :     }
298 :    
299 :     //Evaluate the function with the trial point trialx and return as -trialf
300 :     trialf = EcoSystem->SimulateAndUpdate(trialx);
301 :     trialf = -trialf;
302 :    
303 :     //If too many function evaluations occur, terminate the algorithm
304 :     iters = EcoSystem->getFuncEval() - offset;
305 :     if (iters > simanniter) {
306 :     handle.logMessage(LOGINFO, "\nStopping Simulated Annealing optimisation algorithm\n");
307 :     handle.logMessage(LOGINFO, "The optimisation stopped after", iters, "function evaluations");
308 :     handle.logMessage(LOGINFO, "The temperature was reduced to", t);
309 :     handle.logMessage(LOGINFO, "The optimisation stopped because the maximum number of function evaluations");
310 :     handle.logMessage(LOGINFO, "was reached and NOT because an optimum was found for this run");
311 :     handle.logMessage(LOGINFO, "Number of directly accepted points", nacc);
312 :     handle.logMessage(LOGINFO, "Number of metropolis accepted points", naccmet);
313 :     handle.logMessage(LOGINFO, "Number of rejected points", nrej);
314 :    
315 :     score = EcoSystem->SimulateAndUpdate(bestx);
316 :     handle.logMessage(LOGINFO, "\nSimulated Annealing finished with a likelihood score of", score);
317 :     return;
318 :     }
319 :    
320 :     //Accept the new point if the new function value better
321 :     if ((trialf - funcval) > verysmall) {
322 :     for (i = 0; i < nvars; i++)
323 :     x[i] = trialx[i];
324 :     funcval = trialf;
325 :     nacc++;
326 :     nacp[param[l]]++; //JMB - not sure about this ...
327 :    
328 :     } else {
329 :     //Accept according to metropolis condition
330 :     p = expRep((trialf - funcval) / t);
331 :     pp = randomNumber();
332 :     if (pp < p) {
333 :     //Accept point
334 :     for (i = 0; i < nvars; i++)
335 :     x[i] = trialx[i];
336 :     funcval = trialf;
337 :     naccmet++;
338 :     nacp[param[l]]++;
339 :     } else {
340 :     //Reject point
341 :     nrej++;
342 :     }
343 :     }
344 :    
345 :     // JMB added check for really silly values
346 :     if (isZero(trialf)) {
347 :     handle.logMessage(LOGINFO, "Error in Simulated Annealing optimisation after", iters, "function evaluations, f(x) = 0");
348 :     converge = -1;
349 :     return;
350 :     }
351 :    
352 :     //If greater than any other point, record as new optimum
353 :     if ((trialf > fopt) && (trialf == trialf)) {
354 :     for (i = 0; i < nvars; i++)
355 :     bestx[i] = trialx[i];
356 :     fopt = trialf;
357 :    
358 :     if (scale) {
359 :     for (i = 0; i < nvars; i++)
360 :     scalex[i] = bestx[i] * init[i];
361 :     EcoSystem->storeVariables(-fopt, scalex);
362 :     } else
363 :     EcoSystem->storeVariables(-fopt, bestx);
364 :    
365 :     handle.logMessage(LOGINFO, "\nNew optimum found after", iters, "function evaluations");
366 :     handle.logMessage(LOGINFO, "The likelihood score is", -fopt, "at the point");
367 :     EcoSystem->writeBestValues();
368 :     }
369 :     }
370 :     }
371 :    
372 :     //Adjust vm so that approximately half of all evaluations are accepted
373 :     for (i = 0; i < nvars; i++) {
374 :     ratio = nsdiv * nacp[i];
375 :     nacp[i] = 0;
376 :     if (ratio > uratio) {
377 :     vm[i] = vm[i] * (1.0 + cs * (ratio - uratio));
378 :     } else if (ratio < lratio) {
379 :     vm[i] = vm[i] / (1.0 + cs * (lratio - ratio));
380 :     }
381 :    
382 :     if (vm[i] < rathersmall)
383 :     vm[i] = rathersmall;
384 :     if (vm[i] > (upperb[i] - lowerb[i]))
385 :     vm[i] = upperb[i] - lowerb[i];
386 :     }
387 :     }
388 :    
389 :     //Check termination criteria
390 :     for (i = tempcheck - 1; i > 0; i--)
391 :     fstar[i] = fstar[i - 1];
392 :     fstar[0] = funcval;
393 :    
394 :     quit = 0;
395 :     if (fabs(fopt - funcval) < simanneps) {
396 :     quit = 1;
397 :     for (i = 0; i < tempcheck - 1; i++)
398 :     if (fabs(fstar[i + 1] - fstar[i]) > simanneps)
399 :     quit = 0;
400 :     }
401 :    
402 :     handle.logMessage(LOGINFO, "Checking convergence criteria after", iters, "function evaluations ...");
403 :    
404 :     //Terminate SA if appropriate
405 :     if (quit) {
406 :     handle.logMessage(LOGINFO, "\nStopping Simulated Annealing optimisation algorithm\n");
407 :     handle.logMessage(LOGINFO, "The optimisation stopped after", iters, "function evaluations");
408 :     handle.logMessage(LOGINFO, "The temperature was reduced to", t);
409 :     handle.logMessage(LOGINFO, "The optimisation stopped because an optimum was found for this run");
410 :     handle.logMessage(LOGINFO, "Number of directly accepted points", nacc);
411 :     handle.logMessage(LOGINFO, "Number of metropolis accepted points", naccmet);
412 :     handle.logMessage(LOGINFO, "Number of rejected points", nrej);
413 :    
414 :     converge = 1;
415 :     score = EcoSystem->SimulateAndUpdate(bestx);
416 :     handle.logMessage(LOGINFO, "\nSimulated Annealing finished with a likelihood score of", score);
417 :     return;
418 :     }
419 :    
420 :     //If termination criteria is not met, prepare for another loop.
421 :     t *= rt;
422 :     if (t < rathersmall)
423 :     t = rathersmall; //JMB make sure temperature doesnt get too small
424 :    
425 :     handle.logMessage(LOGINFO, "Reducing the temperature to", t);
426 :     funcval = fopt;
427 :     for (i = 0; i < nvars; i++)
428 :     x[i] = bestx[i];
429 :     }
430 :     }

root@forge.cesga.es
ViewVC Help
Powered by ViewVC 1.0.0  

Powered By FusionForge