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/optinfo.h
[mareframe] / trunk / gadget / optinfo.h Repository:
ViewVC logotype

Annotation of /trunk/gadget/optinfo.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1 - (view) (download)

1 : agomez 1 #ifndef optinfo_h
2 :     #define optinfo_h
3 :    
4 :     #include "maininfo.h"
5 :     #include "doublematrix.h"
6 :     #include "doublevector.h"
7 :     #include "intvector.h"
8 :    
9 :     enum OptType { OPTHOOKE = 1, OPTSIMANN, OPTBFGS };
10 :    
11 :     /**
12 :     * \class OptInfo
13 :     * \brief This is the base class used to perform the optimisation calculation for the model
14 :     *
15 :     * \note This will always be overridden by the derived classes that actually perform the optimisation calculation
16 :     */
17 :     class OptInfo {
18 :     public:
19 :     /**
20 :     * \brief This is the default OptInfo constructor
21 :     */
22 :     OptInfo() { converge = 0; iters = 0; score = 0.0; };
23 :     /**
24 :     * \brief This is the default OptInfo destructor
25 :     */
26 :     ~OptInfo() {};
27 :     /**
28 :     * \brief This is the function used to read in the optimisation parameters
29 :     * \param infile is the CommentStream to read the optimisation parameters from
30 :     * \param text is a text string used to compare parameter names
31 :     */
32 :     virtual void read(CommentStream& infile, char* text) {};
33 :     /**
34 :     * \brief This function will print information from the optimisation algorithm
35 :     * \param outfile is the ofstream that the optimisation information gets sent to
36 :     * \param prec is the precision to use in the output file
37 :     */
38 :     virtual void Print(ofstream& outfile, int prec) {};
39 :     /**
40 :     * \brief This is the function used to call the optimisation algorithms
41 :     */
42 :     virtual void OptimiseLikelihood() {};
43 :     /**
44 :     * \brief This will return the type of optimisation class
45 :     * \return type
46 :     */
47 :     OptType getType() const { return type; };
48 :     protected:
49 :     /**
50 :     * \brief This is the flag used to denote whether the optimisation converged or not
51 :     */
52 :     int converge;
53 :     /**
54 :     * \brief This is the number of iterations that took place during the optimisation
55 :     */
56 :     int iters;
57 :     /**
58 :     * \brief This is the value of the best likelihood score from the optimisation
59 :     */
60 :     double score;
61 :     /**
62 :     * \brief This denotes what type of optimisation class has been created
63 :     */
64 :     OptType type;
65 :     };
66 :    
67 :     /**
68 :     * \class OptInfoHooke
69 :     * \brief This is the class used for the Hooke & Jeeves optimisation
70 :     *
71 :     * The Hooke & Jeeves optimisation is the default optimisation, and is a simple and fast optimising method, but somewhat unreliable, which is often described as a "hill climbing" technique. From the initial starting point the algorithm takes a step in various directions, and conducts a new model run. If the new likelihood score is better than the old one then the algorithm uses the new point as it's best guess. If it is worse then the algorithm retains the old point. The search proceeds in series of these steps, each step slightly smaller than the previous one. When the algorithm finds a point which it cannot improve on with a small step in any direction then it accepts this point as being the "solution", and exits. It is recommended that you re-run the optimisation, using the final point of one run as the start of the next.
72 :     *
73 :     * The Hooke & Jeeves algorithm used in Gadget is derived from that originally presented by R. Hooke and T. A. Jeeves, ''Direct Search Solution of Numerical and Statistical Problems'' in the April 1961 (Vol. 8, pp. 212-229) issue of the Journal of the ACM, with improvements presented by Arthur F Kaupe Jr., ''Algorithm 178: Direct Search'' in the June 1963 (Vol 6, pp.313-314) issue of the Communications of the ACM.
74 :    
75 :     */
76 :     class OptInfoHooke : public OptInfo {
77 :     public:
78 :     /**
79 :     * \brief This is the default OptInfoHooke constructor
80 :     */
81 :     OptInfoHooke();
82 :     /**
83 :     * \brief This is the default OptInfoHooke destructor
84 :     */
85 :     virtual ~OptInfoHooke() {};
86 :     /**
87 :     * \brief This is the function used to read in the Hooke & Jeeves parameters
88 :     * \param infile is the CommentStream to read the optimisation parameters from
89 :     * \param text is a text string used to compare parameter names
90 :     */
91 :     virtual void read(CommentStream& infile, char* text);
92 :     /**
93 :     * \brief This function will print information from the optimisation algorithm
94 :     * \param outfile is the ofstream that the optimisation information gets sent to
95 :     * \param prec is the precision to use in the output file
96 :     */
97 :     virtual void Print(ofstream& outfile, int prec);
98 :     /**
99 :     * \brief This is the function that will calculate the likelihood score using the Hooke & Jeeves optimiser
100 :     */
101 :     virtual void OptimiseLikelihood();
102 :     private:
103 :     /**
104 :     * \brief This function will calculate the best point that can be found close to the current point
105 :     * \param delta is the DoubleVector of the steps to take when looking for the best point
106 :     * \param point is the DoubleVector that will contain the parameters corresponding to the best function value found from the search
107 :     * \param prevbest is the current best point value
108 :     * \param param is the IntVector containing the order that the parameters should be searched in
109 :     * \return the best function value found from the search
110 :     */
111 :     double bestNearby(DoubleVector& delta, DoubleVector& point, double prevbest, IntVector& param);
112 :     /**
113 :     * \brief This is the maximum number of iterations for the Hooke & Jeeves optimisation
114 :     */
115 :     int hookeiter;
116 :     /**
117 :     * \brief This is the reduction factor for the step length
118 :     */
119 :     double rho;
120 :     /**
121 :     * \brief This is the initial step length
122 :     */
123 :     double lambda;
124 :     /**
125 :     * \brief This is the minimum step length, use as the halt criteria for the optimisation process
126 :     */
127 :     double hookeeps;
128 :     /**
129 :     * \brief This is the limit when checking if a parameter is stuck on the bound
130 :     */
131 :     double bndcheck;
132 :     };
133 :    
134 :     /**
135 :     * \class OptInfoSimann
136 :     * \brief This is the class used for the Simualted Annealing optimisation
137 :     *
138 :     * Simulated Annealing is a global optimisation method that distinguishes different local optima. Starting from an initial point, the algorithm takes a step and the function is evaluated. When minimizing a function, any downhill step is accepted and the process repeats from this new point. An uphill step may be accepted (thus, it can escape from local optima). This uphill decision is made by the Metropolis criteria. It uses a parameter known as "temperature" and the size of the uphill step in a probabilistic manner, and varying the temperature will affect the number of the uphill moves that are accepted. As the optimisation process proceeds, the length of the steps decline and the algorithm closes in on the global optimum.
139 :     *
140 :     * The Simulated Annealing algorithm used in Gadget is derived from that presented by Corana et al, ''Minimising Multimodal Functions of Continuous Variables with the 'Simulated Annealing' Algorithm'' in the September 1987 (Vol. 13, pp. 262-280) issue of the ACM Transactions on Mathematical Software and Goffe et al, ''Global Optimisation of Statistical Functions with Simulated Annealing'' in the January/February 1994 (Vol. 60, pp. 65-100) issue of the Journal of Econometrics.
141 :     */
142 :     class OptInfoSimann : public OptInfo {
143 :     public:
144 :     /**
145 :     * \brief This is the default OptInfoSimann constructor
146 :     */
147 :     OptInfoSimann();
148 :     /**
149 :     * \brief This is the default OptInfoSimann destructor
150 :     */
151 :     virtual ~OptInfoSimann() {};
152 :     /**
153 :     * \brief This is the function used to read in the Simulated Annealing parameters
154 :     * \param infile is the CommentStream to read the optimisation parameters from
155 :     * \param text is a text string used to compare parameter names
156 :     */
157 :     virtual void read(CommentStream& infile, char* text);
158 :     /**
159 :     * \brief This function will print information from the optimisation algorithm
160 :     * \param outfile is the ofstream that the optimisation information gets sent to
161 :     * \param prec is the precision to use in the output file
162 :     */
163 :     virtual void Print(ofstream& outfile, int prec);
164 :     /**
165 :     * \brief This is the function that will calculate the likelihood score using the Simulated Annealing optimiser
166 :     */
167 :     virtual void OptimiseLikelihood();
168 :     private:
169 :     /**
170 :     * \brief This is the temperature reduction factor
171 :     */
172 :     double rt;
173 :     /**
174 :     * \brief This is the halt criteria for the Simulated Annealing algorithm
175 :     */
176 :     double simanneps;
177 :     /**
178 :     * \brief This is the number of loops before the step length is adjusted
179 :     */
180 :     int ns;
181 :     /**
182 :     * \brief This is the number of loops before the temperature is adjusted
183 :     */
184 :     int nt;
185 :     /**
186 :     * \brief This is the "temperature" used for the Simulated Annealing algorithm
187 :     */
188 :     double t;
189 :     /**
190 :     * \brief This is the factor used to adjust the step length
191 :     */
192 :     double cs;
193 :     /**
194 :     * \brief This is the initial value for the maximum step length
195 :     */
196 :     double vminit;
197 :     /**
198 :     * \brief This is the maximum number of function evaluations for the Simulated Annealing optimiation
199 :     */
200 :     int simanniter;
201 :     /**
202 :     * \brief This is the upper bound when adjusting the step length
203 :     */
204 :     double uratio;
205 :     /**
206 :     * \brief This is the lower bound when adjusting the step length
207 :     */
208 :     double lratio;
209 :     /**
210 :     * \brief This is the number of temperature loops to check when testing for convergence
211 :     */
212 :     int tempcheck;
213 :     /**
214 :     * \brief This is the flag to denote whether the parameters should be scaled or not (default 0, not scale)
215 :     */
216 :     int scale;
217 :     };
218 :    
219 :     /**
220 :     * \class OptInfoBFGS
221 :     * \brief This is the class used for the BFGS optimisation
222 :     *
223 :     * BFGS is a quasi-Newton global optimisation method that uses information about the gradient of the function at the current point to calculate the best direction to look in to find a better point. Using this information, the BFGS algorithm can iteratively calculate a better approximation to the inverse Hessian matrix, which will lead to a better approximation of the minimum value. From an initial starting point, the gradient of the function is calculated and then the algorithm uses this information to calculate the best direction to perform a linesearch for a point that is ''sufficiently better''. The linesearch that is used in Gadget to look for a better point in this direction is the ''Armijo'' linesearch. The algorithm will then adjust the current estimate of the inverse Hessian matrix, and restart from this new point. If a better point cannot be found, then the inverse Hessian matrix is reset and the algorithm restarts from the last accepted point.
224 :     *
225 :     * The BFGS algorithm used in Gadget is derived from that presented by Dimitri P Bertsekas, ''Nonlinear Programming'' (2nd edition, pp22-61) published by Athena Scientific.
226 :     */
227 :     class OptInfoBFGS : public OptInfo {
228 :     public:
229 :     /**
230 :     * \brief This is the default OptInfoBFGS constructor
231 :     */
232 :     OptInfoBFGS();
233 :     /**
234 :     * \brief This is the default OptInfoBFGS destructor
235 :     */
236 :     ~OptInfoBFGS() {};
237 :     /**
238 :     * \brief This is the function used to read in the BFGS parameters
239 :     * \param infile is the CommentStream to read the optimisation parameters from
240 :     * \param text is a text string used to compare parameter names
241 :     */
242 :     virtual void read(CommentStream& infile, char* text);
243 :     /**
244 :     * \brief This function will print information from the optimisation algorithm
245 :     * \param outfile is the ofstream that the optimisation information gets sent to
246 :     * \param prec is the precision to use in the output file
247 :     */
248 :     virtual void Print(ofstream& outfile, int prec);
249 :     /**
250 :     * \brief This is the function that will calculate the likelihood score using the BFGS optimiser
251 :     */
252 :     virtual void OptimiseLikelihood();
253 :     private:
254 :     /**
255 :     * \brief This function will numerically calculate the gradient of the function at the current point
256 :     * \param point is the DoubleVector that contains the parameters corresponding to the current function value
257 :     * \param pointvalue is the current function value
258 :     * \param newgrad is the DoubleVector that will contain the gradient vector for the current point
259 :     */
260 :     void gradient(DoubleVector& point, double pointvalue, DoubleVector& newgrad);
261 :     /**
262 :     * \brief This function will calculate the smallest eigenvalue of the inverse Hessian matrix
263 :     * \param M is the DoubleMatrix containing the inverse Hessian matrix
264 :     * \return the smallest eigen value of the matrix
265 :     */
266 :     double getSmallestEigenValue(DoubleMatrix M);
267 :     /**
268 :     * \brief This is the maximum number of function evaluations for the BFGS optimiation
269 :     */
270 :     int bfgsiter;
271 :     /**
272 :     * \brief This is the halt criteria for the BFGS algorithm
273 :     */
274 :     double bfgseps;
275 :     /**
276 :     * \brief This is the adjustment factor in the Armijo linesearch
277 :     */
278 :     double beta;
279 :     /**
280 :     * \brief This is the halt criteria for the Armijo linesearch
281 :     */
282 :     double sigma;
283 :     /**
284 :     * \brief This is the initial step size for the Armijo linesearch
285 :     */
286 :     double step;
287 :     /**
288 :     * \brief This is the accuracy term used when calculating the gradient
289 :     */
290 :     double gradacc;
291 :     /**
292 :     * \brief This is the factor used to adjust the gradient accuracy term
293 :     */
294 :     double gradstep;
295 :     /**
296 :     * \brief This is the halt criteria for the gradient accuracy term
297 :     */
298 :     double gradeps;
299 :     };
300 :    
301 :     #endif

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

Powered By FusionForge