#include <Uncmin.h>
Public Types | |
typedef FUNC | function_type |
Type definition of function object to minimize. | |
Public Member Functions | |
Uncmin (FUNC *f) | |
Constructor that sets dimension and function to minimize. | |
Uncmin (int d=1) | |
Constructor that sets dimension, but not function to minimize. | |
~Uncmin () | |
Destructor of uncmin. | |
int | Minimize (V &start, V &xpls, double &fpls, V &gpls, M &hpls) |
Method to compute values that minimize the function. | |
int | Minimize (V &start, V &xpls, double &fpls, V &gpls) |
Version of Minimize in which the hessian is allocated locally rather than passed as an argument to the function. | |
FUNC * | GetFunction () |
Returns pointer to object containing function to minimize, gradient, and hessian. | |
void | SetFunction (FUNC *f) |
Assigns pointer to object containing function to minimize, gradient, and hessian. | |
void | SetFuncExpensive (int iexp) |
Assigns flag indicating whether function is expensive to evaluate. | |
int | SetScaleArg (V typsiz) |
Assigns typical magnitude of each argument to function. | |
void | SetScaleFunc (double fscale) |
Assigns typical magnitude of function near minimum. | |
void | SetNumDigits (double ndigit) |
Assigns the number of reliable digits returned by f_to_minimize. | |
void | SetMaxIter (int maxiter) |
Assigns the maximum number of iterations. | |
void | SetTolerances (double step, double grad, double macheps=-1.0) |
Assigns step tolerance, gradient tolerance, and the machine epsilon. | |
void | SetMaxStep (double stepmx) |
Assigns maximum allowable scaled step length at any iteration. | |
void | SetTrustRegionRadius (double dlt) |
Assigns initial trust region radius. | |
void | SetChecks (int check_gradient, int check_hessian=0) |
Assigns flags to check analytic gradient and hessian. | |
void | SetPrint (std::FILE *file, int print_results, int print_iterations=0) |
Assigns flags to indicate whether results are printed and specify which file is receive the print output. | |
void | SetMethod (int method) |
Assigns method to use to solve minimization problem. | |
void | GetStoppingCriteria (double &gradtol, double &steptol, int &iterations) |
Returns stopping criteria computed at the end of the last call to Minimize. | |
int | GetMessage () |
Returns message index generated in the last call to Minimize. |
Definition at line 201 of file Uncmin.h.
typedef FUNC Uncmin< V, M, FUNC >::function_type |
Constructor that sets dimension and function to minimize.
The minclass function object must define:
1. a method, f_to_minimize, to minimize. f_to_minimize must have the form public static double f_to_minimize(double x[]) where x is the vector of arguments to the function and the return value is the value of the function evaluated at x.
2. a method, gradient, that has the form public static void gradient(double x[], double g[]) where g is the gradient of f evaluated at x. This method will have an empty body if the user does not wish to provide an analytic estimate of the gradient.
3. a method, hessian, that has the form public static void hessian(double x[],double h[][]) where h is the Hessian of f evaluated at x. This method will have an empty body if the user does not wish to provide an analytic estimate of the Hessian. If the user wants Uncmin++ to check the Hessian, then the hessian method should only fill the lower triangle (and diagonal)of h.
Minimization parameters are set at default values.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | *f | Pointer to function object f. |
Definition at line 448 of file Uncmin.h.
00449 : 00450 minclass(f) 00451 { 00452 if (minclass) 00453 n = minclass->dim(); 00454 else 00455 n = 0; 00456 00457 mfile = 0; 00458 00459 mLastGradCrit = 0.0; 00460 mLastStepCrit = 0.0; 00461 mLastIteration = 0; 00462 mLastMessage = 0; 00463 00464 /* Set default parameter values */ 00465 dfault(); }
Constructor that sets dimension, but not function to minimize.
Note: Function to minimize must be set by calling SetFunction. Minimization parameters are set to default values.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | d | Dimension d (number of function arguments). |
Definition at line 484 of file Uncmin.h.
00485 : 00486 minclass(0), typsiz(d) 00487 { 00488 00489 n = d; 00490 00491 mfile = 0; 00492 00493 mLastGradCrit = 0.0; 00494 mLastStepCrit = 0.0; 00495 mLastIteration = 0; 00496 mLastMessage = 0; 00497 00498 /* Set default parameter values */ 00499 dfault(); 00500 }
int Uncmin< V, M, FUNC >::Minimize | ( | V & | start, | |
V & | xpls, | |||
double & | fpls, | |||
V & | gpls, | |||
M & | hpls | |||
) | [inline] |
Method to compute values that minimize the function.
Translated by Brad Hanson from optdrv_f77 by Steve Verrill.
Returns 0 if mLastMessage is 0, 1, 2, or 3 (indicating optimal solution probably found), otherwise returns non-zero.
On exit data member mLastMessage is set to indicate status of solution.
Possible values of mLastMessage are: 0 = Optimal solution found, terminated with gradient small. 1 = Terminated with gradient small, xpls is probably optimal. 2 = Terminated with stepsize small, xpls is probably optimal. 3 = Lower point cannot be found, xpls is probably optimal. 4 = Iteration limit (default 150) exceeded. 5 = Too many large steps, function may be unbounded. -1 = Analytic gradient check requested, but no analytic gradient supplied. -2 = Analytic hessian check requested, but no analytic hessian supplied. -3 = Illegal dimension. -4 = Illegal tolerance. -5 = Illegal iteration limit. -6 = Minimization function has no good digits. -7 = Iteration limit exceeded in lnsrch. -20 = Function not defined at starting value. -21 = Check of analytic gradient failed. -22 = Check of analytic hessian failed.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | &start | Address of vector containing initial estimate of location of minimum. |
[out] | &xpls | Address of vector containing local minimum (size must be n on input). |
[out] | &fpls | Address of function value at local minimum. |
[out] | &gpls | Address of vector containing gradient at local minimum (size must be n on input). |
[out] | &hpls | Address of matrix containing hessian at local minimum (size must be n X n on input). Only the lower half should be input; the upper off-diagonal entries can be set to zero. On return, the lower half of hpls contains the Cholesky factor of the Hessian, evaluated at the solution (i.e., H = LL'). |
Definition at line 1000 of file Uncmin.h.
Referenced by Uncmin< V, M, FUNC >::Minimize().
01003 { 01004 01005 bool noupdt; 01006 bool mxtake; 01007 01008 int i, j, num5, remain, ilow, ihigh; 01009 int icscmx; 01010 int iretcd; 01011 int itncnt; 01012 01013 /* 01014 Return code. 01015 itrmcd = 0: Optimal solution found 01016 itrmcd = 1: Terminated with gradient small, X is probably optimal 01017 itrmcd = 2: Terminated with stepsize small, X is probably optimal 01018 itrmcd = 3: Lower point cannot be found, X is probably optimal 01019 itrmcd = 4: Iteration limit (150) exceeded 01020 itrmcd = 5: Too many large steps, function may be unbounded. 01021 */ 01022 int itrmcd; 01023 01024 double rnf, analtl, dltsav, amusav, dlpsav, phisav, phpsav; 01025 01026 double f; 01027 01028 double amu = 0.0; 01029 double dltp = 0.0; 01030 double phi = 0.0; 01031 double phip0 = 0.0; 01032 01033 int method = algorithm; 01034 01035 mLastGradCrit = 0.0; 01036 mLastStepCrit = 0.0; 01037 mLastIteration = 0; 01038 mLastMessage = -9999; 01039 01040 iagflg = minclass->HasAnalyticGradient(); 01041 iahflg = minclass->HasAnalyticHessian(); 01042 01043 /* 01044 msg 01045 Flag to set various options. 01046 Bits which indicate flags which control checks of gradient and hessian and what output is printed. 01047 msg should be zero or a sum of the following values which indicate which flags to set or not set. 01048 01049 1 = Turn off check of analytic gradient. 01050 2 = Turn off check of analytic hessian. 01051 4 = Turn off printing of results. 01052 8 = Turn on printing of results at each iteration. 01053 */ 01054 int msg = (fCheckGradient == 0) * 1 + (fCheckHessian == 0) * 2 + (fPrintResults == 0) * 4 01055 + (fPrintIterationResults != 0) * 8; 01056 01057 /* Check that dimensions match */ 01058 if (start.size() != n || xpls.size() != n || gpls.size() != n || hpls.num_rows() != n 01059 || hpls.num_columns() != n) // Refactored function "num_columns", ww, 12-22-2007 01060 { 01061 mLastMessage = -3; 01062 return -1; 01063 } 01064 01065 /* Check that starting values are valid */ 01066 if (!(minclass->ValidParameters(start))) 01067 { 01068 mLastMessage = -20; 01069 return -1; 01070 } 01071 01072 /* Workspace */ 01073 M &a = hpls; 01074 V udiag(n); //workspace (for diagonal of Hessian) 01075 V g(n); // workspace (for gradient at current iterate) 01076 V p(n); // workspace for step 01077 V sx(n); // workspace (for scaling vector) 01078 V wrk0(n), wrk1(n), wrk2(n), wrk3(n); 01079 01080 V x(start); 01081 01082 double dlt = mdlt; 01083 01084 dltsav = amusav = dlpsav = phisav = phpsav = 0.0; 01085 01086 // INITIALIZATION 01087 01088 for (i = 1; i <= n; i++) 01089 { 01090 p(i) = 0.0; 01091 } 01092 01093 itncnt = 0; 01094 iretcd = -1; 01095 01096 /* Check for valid options */ 01097 if (optchk(x, sx, msg, method)) 01098 { 01099 mLastMessage = msg; 01100 return -1; 01101 } 01102 01103 rnf = max(std::pow(10.0, -ndigit), epsm); 01104 analtl = max(.01, std::sqrt(rnf)); 01105 01106 if (!(msg & 4)) 01107 { 01108 num5 = n/5; 01109 remain = n%5; 01110 01111 std::fprintf(mfile, "\n\nOPTDRV Typical x\n\n"); 01112 01113 ilow = -4; 01114 ihigh = 0; 01115 01116 for (i = 1; i <= num5; i++) 01117 { 01118 ilow += 5; 01119 ihigh += 5; 01120 01121 std::fprintf(mfile, "%d--%d ", ilow, ihigh); 01122 01123 for (j = 1; j <= 5; j++) 01124 { 01125 std::fprintf(mfile, "%lf ", (double) typsiz(ilow+j-1)); 01126 } 01127 std::fprintf(mfile, "\n"); 01128 } 01129 01130 ilow += 5; 01131 ihigh = ilow + remain - 1; 01132 01133 std::fprintf(mfile, "%d--%d ", ilow, ihigh); 01134 01135 for (j = 1; j <= remain; j++) 01136 { 01137 std::fprintf(mfile, "%lf ", (double) typsiz(ilow+j-1)); 01138 } 01139 01140 std::fprintf(mfile, "\n"); 01141 01142 std::fprintf(mfile, "\n\nOPTDRV Scaling vector for x\n\n"); 01143 01144 ilow = -4; 01145 ihigh = 0; 01146 01147 for (i = 1; i <= num5; i++) 01148 { 01149 ilow += 5; 01150 ihigh += 5; 01151 01152 std::fprintf(mfile, "%d--%d ", ilow, ihigh); 01153 01154 for (j = 1; j <= 5; j++) 01155 { 01156 std::fprintf(mfile, "%lf ", (double) sx(ilow+j-1)); 01157 } 01158 std::fprintf(mfile, "\n"); 01159 } 01160 01161 ilow += 5; 01162 ihigh = ilow + remain - 1; 01163 01164 std::fprintf(mfile, "%d--%d ", ilow, ihigh); 01165 01166 for (j = 1; j <= remain; j++) 01167 { 01168 std::fprintf(mfile, "%lf ", (double) sx(ilow+j-1)); 01169 } 01170 01171 std::fprintf(mfile, "\n"); 01172 01173 std::fprintf(mfile, "\n\nOPTDRV Typical f = %f\n", fscale); 01174 std::fprintf(mfile, "OPTDRV Number of good digits in"); 01175 std::fprintf(mfile, " f_to_minimize = %d\n", ndigit); 01176 std::fprintf(mfile, "OPTDRV Gradient flag"); 01177 std::fprintf(mfile, " = %d\n", iagflg); 01178 std::fprintf(mfile, "OPTDRV Hessian flag"); 01179 std::fprintf(mfile, " = %d\n", iahflg); 01180 std::fprintf(mfile, "OPTDRV Expensive function calculation flag"); 01181 std::fprintf(mfile, " = %d\n", iexp); 01182 std::fprintf(mfile, "OPTDRV Method to use"); 01183 std::fprintf(mfile, " = %d\n", method); 01184 std::fprintf(mfile, "OPTDRV Iteration limit"); 01185 std::fprintf(mfile, " = %d\n", itnlim); 01186 std::fprintf(mfile, "OPTDRV Machine epsilon"); 01187 std::fprintf(mfile, " = %e\n", epsm); 01188 std::fprintf(mfile, "OPTDRV Maximum step size"); 01189 std::fprintf(mfile, " = %f\n", stepmx); 01190 std::fprintf(mfile, "OPTDRV Step tolerance"); 01191 std::fprintf(mfile, " = %e\n", steptl); 01192 std::fprintf(mfile, "OPTDRV Gradient tolerance"); 01193 std::fprintf(mfile, " = %e\n", gradtl); 01194 std::fprintf(mfile, "OPTDRV Trust region radius"); 01195 std::fprintf(mfile, " = %f\n", dlt); 01196 std::fprintf(mfile, "OPTDRV Relative noise in"); 01197 std::fprintf(mfile, " f_to_minimize = %e\n", rnf); 01198 std::fprintf(mfile, "OPTDRV Analytical fd tolerance"); 01199 std::fprintf(mfile, " = %e\n", analtl); 01200 } 01201 01202 // EVALUATE FCN(X) 01203 f = minclass->f_to_minimize(x); 01204 01205 // EVALUATE ANALYTIC OR FINITE DIFFERENCE GRADIENT AND CHECK ANALYTIC 01206 // GRADIENT, IF REQUESTED. 01207 if (iagflg == 0) 01208 { 01209 fstofd(x, f, g, sx, rnf); 01210 } 01211 else 01212 { 01213 minclass->gradient(x, g); 01214 01215 if (!(msg & 1)) 01216 { 01217 if (grdchk(x, f, g, typsiz, sx, fscale, rnf, analtl, wrk1, msg)) 01218 { 01219 mLastMessage = msg; 01220 return -1; // msg = -21 01221 } 01222 } 01223 } 01224 01225 optstp(x, f, g, wrk1, itncnt, icscmx, itrmcd, sx, fscale, itnlim, iretcd, mxtake, msg); 01226 01227 if (itrmcd != 0) 01228 { 01229 fpls = f; 01230 01231 for (i = 1; i <= n; i++) 01232 { 01233 xpls(i) = x(i); 01234 gpls(i) = g(i); 01235 } 01236 } 01237 else 01238 { 01239 if (iexp == 1) 01240 { 01241 // IF OPTIMIZATION FUNCTION EXPENSIVE TO EVALUATE (IEXP=1), THEN 01242 // HESSIAN WILL BE OBTAINED BY SECANT (BFGS) UPDATES. GET INITIAL HESSIAN. 01243 01244 hsnint(a, sx, method); 01245 } 01246 else 01247 { 01248 // EVALUATE ANALYTIC OR FINITE DIFFERENCE HESSIAN AND CHECK ANALYTIC 01249 // HESSIAN IF REQUESTED (ONLY IF USER-SUPPLIED ANALYTIC HESSIAN 01250 // ROUTINE minclass->hessian FILLS ONLY LOWER TRIANGULAR PART AND DIAGONAL OF A). 01251 01252 if (iahflg == 0) 01253 { 01254 if (iagflg == 1) 01255 { 01256 fstofd(x, g, a, sx, rnf, wrk1); 01257 } 01258 else 01259 { 01260 sndofd(x, f, a, sx, rnf, wrk1, wrk2); 01261 } 01262 } 01263 else 01264 { 01265 if (msg & 2) 01266 { 01267 minclass->hessian(x, a); 01268 } 01269 else 01270 { 01271 if (heschk(x, f, g, a, typsiz, sx, rnf, analtl, iagflg, udiag, wrk1, wrk2, msg)) 01272 { 01273 mLastMessage = msg; 01274 return -1; // msg = -22 01275 } 01276 // HESCHK EVALUATES minclass->hessian AND CHECKS IT AGAINST THE FINITE 01277 // DIFFERENCE HESSIAN WHICH IT CALCULATES BY CALLING FSTOFD 01278 // (IF IAGFLG .EQ. 1) OR SNDOFD (OTHERWISE). 01279 } 01280 } 01281 } 01282 01283 if (!(msg & 4)) 01284 result(x, f, g, a, p, itncnt, 1); 01285 01286 // ITERATION 01287 while (itrmcd == 0) 01288 { 01289 itncnt++; 01290 mLastIteration = itncnt; // BAH 01291 01292 // FIND PERTURBED LOCAL MODEL HESSIAN AND ITS LL+ DECOMPOSITION 01293 // (SKIP THIS STEP IF LINE SEARCH OR DOGSTEP TECHNIQUES BEING USED WITH 01294 // SECANT UPDATES. CHOLESKY DECOMPOSITION L ALREADY OBTAINED FROM 01295 // SECFAC.) 01296 if (iexp != 1 || method == 3) 01297 { 01298 chlhsn(a, sx, udiag); 01299 } 01300 01301 // SOLVE FOR NEWTON STEP: AP=-G 01302 for (i = 1; i <= n; i++) 01303 { 01304 wrk1(i) = -g(i); 01305 } 01306 01307 lltslv(a, p, wrk1); 01308 01309 // DECIDE WHETHER TO ACCEPT NEWTON STEP XPLS=X + P 01310 // OR TO CHOOSE XPLS BY A GLOBAL STRATEGY. 01311 if (iagflg == 0 && method != 1) 01312 { 01313 dltsav = dlt; 01314 01315 if (method != 2) 01316 { 01317 amusav = amu; 01318 dlpsav = dltp; 01319 phisav = phi; 01320 phpsav = phip0; 01321 } 01322 } 01323 01324 if (method == 1) 01325 { 01326 lnsrch(x, f, g, p, xpls, fpls, mxtake, iretcd, sx); 01327 if (iretcd == -7) // BAH 01328 { 01329 mLastMessage = iretcd; 01330 return -1; 01331 } 01332 } 01333 else if (method == 2) 01334 { 01335 dogdrv(x, f, g, a, p, xpls, fpls, sx, dlt, iretcd, mxtake, wrk0, wrk1, wrk2, wrk3); 01336 } 01337 else 01338 { 01339 hookdr(x, f, g, a, udiag, p, xpls, fpls, sx, dlt, iretcd, mxtake, amu, dltp, phi, phip0, 01340 wrk0, wrk1, wrk2, itncnt); 01341 } 01342 01343 // IF COULD NOT FIND SATISFACTORY STEP AND FORWARD DIFFERENCE 01344 // GRADIENT WAS USED, RETRY USING CENTRAL DIFFERENCE GRADIENT. 01345 if (iretcd == 1 && iagflg == 0) 01346 { 01347 // SET IAGFLG FOR CENTRAL DIFFERENCES 01348 iagflg = -1; 01349 01350 if (!(msg & 4)) 01351 { 01352 std::fprintf(mfile, "\nOPTDRV Shift from forward to central"); 01353 std::fprintf(mfile, " differences"); 01354 std::fprintf(mfile, "\nOPTDRV in iteration %d", itncnt); 01355 std::fprintf(mfile, "\n"); 01356 } 01357 01358 fstocd(x, sx, rnf, g); 01359 01360 if (method == 1) 01361 { 01362 // SOLVE FOR NEWTON STEP: AP=-G 01363 for (i = 1; i <= n; i++) 01364 { 01365 wrk1(i) = -g(i); 01366 } 01367 01368 lltslv(a, p, wrk1); 01369 01370 lnsrch(x, f, g, p, xpls, fpls, mxtake, iretcd, sx); 01371 if (iretcd == -7) // BAH 01372 { 01373 mLastMessage = iretcd; 01374 return -1; 01375 } 01376 } 01377 else 01378 { 01379 dlt = dltsav; 01380 if (method == 2) 01381 { 01382 // SOLVE FOR NEWTON STEP: AP=-G 01383 for (i = 1; i <= n; i++) 01384 { 01385 wrk1(i) = -g(i); 01386 } 01387 01388 lltslv(a, p, wrk1); 01389 01390 dogdrv(x, f, g, a, p, xpls, fpls, sx, dlt, iretcd, mxtake, wrk0, wrk1, wrk2, wrk3); 01391 } 01392 else 01393 { 01394 amu = amusav; 01395 dltp = dlpsav; 01396 phi = phisav; 01397 phip0 = phpsav; 01398 01399 chlhsn(a, sx, udiag); 01400 01401 // SOLVE FOR NEWTON STEP: AP=-G 01402 for (i = 1; i <= n; i++) 01403 { 01404 wrk1(i) = -g(i); 01405 } 01406 01407 lltslv(a, p, wrk1); 01408 01409 hookdr(x, f, g, a, udiag, p, xpls, fpls, sx, dlt, iretcd, mxtake, amu, dltp, phi, 01410 phip0, wrk0, wrk1, wrk2, itncnt); 01411 } 01412 } 01413 } 01414 // CALCULATE STEP FOR OUTPUT 01415 for (i = 1; i <= n; i++) 01416 { 01417 p(i) = xpls(i) - x(i); 01418 } 01419 01420 // CALCULATE GRADIENT AT XPLS 01421 if (iagflg == -1) 01422 { 01423 // CENTRAL DIFFERENCE GRADIENT 01424 fstocd(xpls, sx, rnf, gpls); 01425 } 01426 else if (iagflg == 0) 01427 { 01428 // FORWARD DIFFERENCE GRADIENT 01429 fstofd(xpls, fpls, gpls, sx, rnf); 01430 } 01431 else 01432 { 01433 // ANALYTIC GRADIENT 01434 minclass->gradient(xpls, gpls); 01435 } 01436 01437 // CHECK WHETHER STOPPING CRITERIA SATISFIED 01438 optstp(xpls, fpls, gpls, x, itncnt, icscmx, itrmcd, sx, fscale, itnlim, iretcd, mxtake, msg); 01439 01440 if (itrmcd == 0) 01441 { 01442 // EVALUATE HESSIAN AT XPLS 01443 if (iexp != 0) 01444 { 01445 if (method == 3) 01446 { 01447 secunf(x, g, a, udiag, xpls, gpls, itncnt, rnf, iagflg, noupdt, wrk1, wrk2, wrk3); 01448 } 01449 else 01450 { 01451 secfac(x, g, a, xpls, gpls, itncnt, rnf, iagflg, noupdt, wrk0, wrk1, wrk2, wrk3); 01452 } 01453 } 01454 else 01455 { 01456 if (iahflg == 1) 01457 { 01458 minclass->hessian(xpls, a); 01459 } 01460 else 01461 { 01462 if (iagflg == 1) 01463 { 01464 fstofd(xpls, gpls, a, sx, rnf, wrk1); 01465 } 01466 else 01467 { 01468 sndofd(xpls, fpls, a, sx, rnf, wrk1, wrk2); 01469 } 01470 } 01471 } 01472 01473 if (msg & 8) 01474 { 01475 result(xpls, fpls, gpls, a, p, itncnt, 1); 01476 } 01477 01478 // X <-- XPLS AND G <-- GPLS AND F <-- FPLS 01479 f = fpls; 01480 01481 for (i = 1; i <= n; i++) 01482 { 01483 x(i) = xpls(i); 01484 g(i) = gpls(i); 01485 } 01486 } 01487 } 01488 01489 // TERMINATION 01490 // ----------- 01491 // RESET XPLS,FPLS,GPLS, IF PREVIOUS ITERATE SOLUTION 01492 // 01493 if (itrmcd == 3) 01494 { 01495 fpls = f; 01496 01497 for (i = 1; i <= n; i++) 01498 { 01499 xpls(i) = x(i); 01500 gpls(i) = g(i); 01501 } 01502 } 01503 } 01504 01505 // PRINT RESULTS 01506 if (!(msg & 4)) 01507 { 01508 result(xpls, fpls, gpls, a, p, itncnt, 0); 01509 } 01510 mLastMessage = itrmcd; 01511 return (itrmcd >= 0 && itrmcd <= 3) ? 0 : 1;
int Uncmin< V, M, FUNC >::Minimize | ( | V & | start, | |
V & | xpls, | |||
double & | fpls, | |||
V & | gpls | |||
) | [inline] |
Version of Minimize in which the hessian is allocated locally rather than passed as an argument to the function.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | &start | Address of vector containing initial estimate of location of minimum. |
[out] | &xpls | Address of vector containing local minimum (size must be n on input). |
[out] | &fpls | Address of function value at local minimum. |
[out] | &gpls | Address of vector containing gradient at local minimum (size must be n on input). |
Definition at line 1532 of file Uncmin.h.
References Uncmin< V, M, FUNC >::Minimize().
01535 { 01536 01537 M *a = new M(n,n); 01538 01539 int result = Minimize(start, xpls, fpls, gpls, *a); 01540 01541 delete a; 01542 01543 return result;
FUNC * Uncmin< V, M, FUNC >::GetFunction | ( | ) | [inline] |
Returns pointer to object containing function to minimize, gradient, and hessian.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
Definition at line 766 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetFunction | ( | FUNC * | f | ) | [inline] |
Assigns pointer to object containing function to minimize, gradient, and hessian.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | *f | Pointer to function object f. |
Definition at line 786 of file Uncmin.h.
00788 { 00789 minclass = f; 00790 00791 /* If dimension does not match change n and re-intialize typsiz */ 00792 if (n != f->dim()) 00793 { 00794 n = f->dim(); 00795 V defaultTypSiz(n, 1.0); 00796 typsiz = defaultTypSiz; 00797 }
void Uncmin< V, M, FUNC >::SetFuncExpensive | ( | int | iexp | ) | [inline] |
Assigns flag indicating whether function is expensive to evaluate.
Set iexp =1 if function is expensive to evaluate, = 0 otherwise. If iexp = 1, then the Hessian will be evaluated by secant (BFGS) update rather than analytically or by finite differences.
Default value is 0 if this function is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | iexp | Set iexp =1 if function is expensive to evaluate, iexp = 0 otherwise. |
Definition at line 530 of file Uncmin.h.
int Uncmin< V, M, FUNC >::SetScaleArg | ( | V | typsiz | ) | [inline] |
Assigns typical magnitude of each argument to function.
Returns zero if dimension matches, otherwise returns 1.
The default is to set all elements of typsiz to 1 if this function is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | typsiz | Vector describing the magnitude of all function arguments. |
Definition at line 553 of file Uncmin.h.
00555 { 00556 if (typsiz.size() != n) 00557 return 1; 00558 00559 this->typsiz = typsiz; 00560 00561 return 0;
void Uncmin< V, M, FUNC >::SetScaleFunc | ( | double | fscale | ) | [inline] |
Assigns typical magnitude of function near minimum.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | fscale | Magnitude of function values near minimum. A default value of 1 is used for fscale if this function is not called. |
Definition at line 579 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetNumDigits | ( | double | ndigit | ) | [inline] |
Assigns the number of reliable digits returned by f_to_minimize.
If the argument to the function is -1 this means that f_to_minimize is expected to provide within one or two of the full number of significant digits available.
The default value of -1 is used for ndigit when this function is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | ndigit | Number of digits precision required for the solution. |
Definition at line 606 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetMaxIter | ( | int | maxiter | ) | [inline] |
Assigns the maximum number of iterations.
The default value of 150 is used if this function is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | maxiter | Maximum number of iterations. |
Definition at line 627 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetTolerances | ( | double | step, | |
double | grad, | |||
double | macheps = -1.0 | |||
) | [inline] |
Assigns step tolerance, gradient tolerance, and the machine epsilon.
Note: For any of these values <= zero, default values are used.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | step | Tolerance at which scaled distance between two successive iterates is considered close enough to zero to terminate algorithm. If the value of step is <= zero a default value is used. |
[in] | grad | Tolerance at which gradient is considered close enough to zero to terminate algorithm. If the value of grad is <= zero a default value is used. |
[in] | macheps | Largest relative spacing. Machine epsilon == 2**(-T+1) where (1 + 2**-T) == 1. This argument has a default value of -1 if not included in function call. |
Definition at line 656 of file Uncmin.h.
00659 { 00660 if (macheps <= 0.0) // use default value 00661 { 00662 00663 #ifdef BOOST_NO_LIMITS 00664 epsm = DBL_EPSILON; // from float.h 00665 #else 00666 // use epsilon from numeric_limits in standard C++ library 00667 epsm = std::numeric_limits<double>::epsilon(); 00668 #endif 00669 00670 } 00671 else 00672 { 00673 epsm = macheps; 00674 } 00675 00676 if (grad <= 0.0) 00677 { 00678 gradtl = std::pow(epsm, 1.0/3.0); // default value 00679 } 00680 else 00681 { 00682 gradtl = grad; 00683 } 00684 00685 if (step <= 0.0) 00686 { 00687 steptl = std::sqrt(epsm); // default value 00688 } 00689 else 00690 { 00691 steptl = step; 00692 } 00693
void Uncmin< V, M, FUNC >::SetMaxStep | ( | double | stepmx | ) | [inline] |
Assigns maximum allowable scaled step length at any iteration.
The maximum step length is used to prevent steps that would cause the optimization algorithm to overflow or leave the domain of iterest, as well as to detect divergence. It should be chosen small enough to prevent the first two of these occurrences, but larger than any anticipated reasonable step size. The algorithm will halt if it takes steps larger than the maximum allowable step length on 5 consecutive iterations.
If the argument to the function is <= 0 then a default value of stepmx will be computed in optchk.
A default value of 0 is used for stepmx when this function is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | stepmx | Maximum (scaled) step size. |
Definition at line 723 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetTrustRegionRadius | ( | double | dlt | ) | [inline] |
Assigns initial trust region radius.
Used when method = 2 or 3, ignored when method = 1. The value should be what the user considers a reasonable scaled step length for the first iteration, and should be less than the maximum step length.
If dlt = -1, then the length of the initial scaled gradient is used instead.
A default value of -1 is used if this function is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | dlt | Maximum (scaled) step size of first iteration. |
Definition at line 751 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetChecks | ( | int | check_gradient, | |
int | check_hessian = 0 | |||
) | [inline] |
Assigns flags to check analytic gradient and hessian.
The default is to not check the gradient and hessian, if this function is not called (fCheckGradient = 0 and fCheckHessian = 0).
Note that check_hessian can be left off function call in which case default value of 0 is used.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | check_gradient | Flag to check analytic gradient of function: 1 = yes, 0 = no. |
[in] | check_hessian | Flag to check analytic hessian of function: 1 = yes, 0 = no. |
Definition at line 820 of file Uncmin.h.
void Uncmin< V, M, FUNC >::SetPrint | ( | std::FILE * | file, | |
int | print_results, | |||
int | print_iterations = 0 | |||
) | [inline] |
Assigns flags to indicate whether results are printed and specify which file is receive the print output.
Printing is done to open file indicated by first argument.
The default action is not to print any output.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. | |
FILE | Type of file for writing print output to. |
[in] | *file | Pointer to the open file that will receive the printed results. |
[in] | print_results | Flag (if non-zero) to print standard results. |
[in] | print_iterations | Flag (if non-zero) to print intermediate results at each iteration (this argument can be left off function call in which case the default value of 0 is used). |
Definition at line 851 of file Uncmin.h.
00854 { 00855 00856 mfile = file; 00857 00858 if (!file) 00859 { 00860 fPrintResults = 0; 00861 fPrintIterationResults = 0; 00862 } 00863 else 00864 { 00865 fPrintResults = print_results; 00866 fPrintIterationResults = print_iterations; 00867 } 00868
void Uncmin< V, M, FUNC >::SetMethod | ( | int | method | ) | [inline] |
Assigns method to use to solve minimization problem.
A default value is 1 (line search) if SetMethod is not called.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[in] | method | Method indicator: 1 = line search, 2 = double dogleg, 3 = More-Hebdon. |
Definition at line 886 of file Uncmin.h.
00888 { 00889 if (method == 2 || method == 3) 00890 algorithm = method; 00891 else 00892 algorithm = 1; // use default value of 1 for any values of argument other than 2 or 3
void Uncmin< V, M, FUNC >::GetStoppingCriteria | ( | double & | gradtol, | |
double & | steptol, | |||
int & | iterations | |||
) | [inline] |
Returns stopping criteria computed at the end of the last call to Minimize.
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
[out] | &gradtol | Address of gradient tolerance, i.e., how close the gradient is to zero. |
[out] | &steptol | Address of maximum scaled difference between function parameters in consecutive iterations. |
[out] | &iterations | Address of number of iterations used. |
Definition at line 911 of file Uncmin.h.
00914 { 00915 gradtol = mLastGradCrit; 00916 steptol = mLastStepCrit; 00917 iterations = mLastIteration;
int Uncmin< V, M, FUNC >::GetMessage | ( | ) | [inline] |
Returns message index generated in the last call to Minimize.
Returns message generated in the last call to Minimize.
Possible values of message are:
0 = Optimal solution found, terminated with gradient small. 1 = Terminated with gradient small, xpls is probably optimal. 2 = Terminated with stepsize small, xpls is probably optimal. 3 = Lower point cannot be found, xpls is probably optimal. 4 = Iteration limit (default 150) exceeded. 5 = Too many large steps, function may be unbounded. -1 = Analytic gradient check requested, but no analytic gradient supplied -2 = Analytic hessian check requested, but no analytic hessian supplied -3 = Illegal dimension -4 = Illegal tolerance -5 = Illegal iteration limit -6 = Minimization function has no good digits -7 = Iteration limit exceeded in lnsrch -20 = Function not defined at starting value -21 = Check of analytic gradient failed -22 = Check of analytic hessian failed
M | SCPPNT Matrix type. | |
V | SCPPNT Vector type. | |
FUNC | Type of function to minimize. |
Definition at line 948 of file Uncmin.h.