Actual source code: morethuente.c
1: /*$Id$*/
3: #include "morethuente.h" /*I "tao_solver.h" I*/
4: #include <stdlib.h>
6: /*
7: The subroutine mcstep is taken from the work of Jorge Nocedal.
8: this is a variant of More' and Thuente's routine.
10: subroutine mcstep
12: the purpose of mcstep is to compute a safeguarded step for
13: a linesearch and to update an interval of uncertainty for
14: a minimizer of the function.
16: the parameter stx contains the step with the least function
17: value. the parameter stp contains the current step. it is
18: assumed that the derivative at stx is negative in the
19: direction of the step. if bracket is set true then a
20: minimizer has been bracketed in an interval of uncertainty
21: with endpoints stx and sty.
23: the subroutine statement is
25: subroutine mcstep(stx,fx,dx,sty,fy,dy,stp,fp,dp,bracket,
26: stpmin,stpmax,info)
28: where
30: stx, fx, and dx are variables which specify the step,
31: the function, and the derivative at the best step obtained
32: so far. The derivative must be negative in the direction
33: of the step, that is, dx and stp-stx must have opposite
34: signs. On output these parameters are updated appropriately.
36: sty, fy, and dy are variables which specify the step,
37: the function, and the derivative at the other endpoint of
38: the interval of uncertainty. On output these parameters are
39: updated appropriately.
41: stp, fp, and dp are variables which specify the step,
42: the function, and the derivative at the current step.
43: If bracket is set true then on input stp must be
44: between stx and sty. On output stp is set to the new step.
46: bracket is a logical variable which specifies if a minimizer
47: has been bracketed. If the minimizer has not been bracketed
48: then on input bracket must be set false. If the minimizer
49: is bracketed then on output bracket is set true.
51: stpmin and stpmax are input variables which specify lower
52: and upper bounds for the step.
54: info is an integer output variable set as follows:
55: if info = 1,2,3,4,5, then the step has been computed
56: according to one of the five cases below. otherwise
57: info = 0, and this indicates improper input parameters.
59: subprograms called
61: fortran-supplied ... abs,max,min,sqrt
63: argonne national laboratory. minpack project. june 1983
64: jorge j. more', david j. thuente
66: */
70: static int TaoStep_LineSearch(TAO_SOLVER tao,
71: double *stx, double *fx, double *dx,
72: double *sty, double *fy, double *dy,
73: double *stp, double *fp, double *dp)
74: {
75: TAO_LINESEARCH *neP = (TAO_LINESEARCH *) tao->linectx;
76: double gamma1, p, q, r, s, sgnd, stpc, stpf, stpq, theta;
77: int bound;
79: TaoFunctionBegin;
81: // Check the input parameters for errors
82: neP->infoc = 0;
83: if (neP->bracket && (*stp <= TaoMin(*stx,*sty) || (*stp >= TaoMax(*stx,*sty)))) SETERRQ(1,"bad stp in bracket");
84: if (*dx * (*stp-*stx) >= 0.0) SETERRQ(1,"dx * (stp-stx) >= 0.0");
85: if (neP->stepmax < neP->stepmin) SETERRQ(1,"stepmax > stepmin");
87: // Determine if the derivatives have opposite sign */
88: sgnd = *dp * (*dx / TaoAbsDouble(*dx));
90: if (*fp > *fx) {
91: // Case 1: a higher function value.
92: // The minimum is bracketed. If the cubic step is closer
93: // to stx than the quadratic step, the cubic step is taken,
94: // else the average of the cubic and quadratic steps is taken.
96: neP->infoc = 1;
97: bound = 1;
98: theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
99: s = TaoMax(TaoAbsDouble(theta),TaoAbsDouble(*dx));
100: s = TaoMax(s,TaoAbsDouble(*dp));
101: gamma1 = s*sqrt(pow(theta/s,2.0) - (*dx/s)*(*dp/s));
102: if (*stp < *stx) gamma1 = -gamma1;
103: /* Can p be 0? Check */
104: p = (gamma1 - *dx) + theta;
105: q = ((gamma1 - *dx) + gamma1) + *dp;
106: r = p/q;
107: stpc = *stx + r*(*stp - *stx);
108: stpq = *stx + ((*dx/((*fx-*fp)/(*stp-*stx)+*dx))*0.5) * (*stp - *stx);
110: if (TaoAbsDouble(stpc-*stx) < TaoAbsDouble(stpq-*stx)) {
111: stpf = stpc;
112: }
113: else {
114: stpf = stpc + 0.5*(stpq - stpc);
115: }
116: neP->bracket = 1;
117: }
118: else if (sgnd < 0.0) {
119: // Case 2: A lower function value and derivatives of
120: // opposite sign. The minimum is bracketed. If the cubic
121: // step is closer to stx than the quadratic (secant) step,
122: // the cubic step is taken, else the quadratic step is taken.
124: neP->infoc = 2;
125: bound = 0;
126: theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp;
127: s = TaoMax(TaoAbsDouble(theta),TaoAbsDouble(*dx));
128: s = TaoMax(s,TaoAbsDouble(*dp));
129: gamma1 = s*sqrt(pow(theta/s,2.0) - (*dx/s)*(*dp/s));
130: if (*stp > *stx) gamma1 = -gamma1;
131: p = (gamma1 - *dp) + theta;
132: q = ((gamma1 - *dp) + gamma1) + *dx;
133: r = p/q;
134: stpc = *stp + r*(*stx - *stp);
135: stpq = *stp + (*dp/(*dp-*dx))*(*stx - *stp);
137: if (TaoAbsDouble(stpc-*stp) > TaoAbsDouble(stpq-*stp)) {
138: stpf = stpc;
139: }
140: else {
141: stpf = stpq;
142: }
143: neP->bracket = 1;
144: }
145: else if (TaoAbsDouble(*dp) < TaoAbsDouble(*dx)) {
146: // Case 3: A lower function value, derivatives of the
147: // same sign, and the magnitude of the derivative decreases.
148: // The cubic step is only used if the cubic tends to infinity
149: // in the direction of the step or if the minimum of the cubic
150: // is beyond stp. Otherwise the cubic step is defined to be
151: // either stepmin or stepmax. The quadratic (secant) step is also
152: // computed and if the minimum is bracketed then the the step
153: // closest to stx is taken, else the step farthest away is taken.
155: neP->infoc = 3;
156: bound = 1;
157: theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp;
158: s = TaoMax(TaoAbsDouble(theta),TaoAbsDouble(*dx));
159: s = TaoMax(s,TaoAbsDouble(*dp));
161: // The case gamma1 = 0 only arises if the cubic does not tend
162: // to infinity in the direction of the step.
163: gamma1 = s*sqrt(TaoMax(0.0,pow(theta/s,2.0) - (*dx/s)*(*dp/s)));
164: if (*stp > *stx) gamma1 = -gamma1;
165: p = (gamma1 - *dp) + theta;
166: q = (gamma1 + (*dx - *dp)) + gamma1;
167: r = p/q;
168: if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r*(*stx - *stp);
169: else if (*stp > *stx) stpc = neP->stepmax;
170: else stpc = neP->stepmin;
171: stpq = *stp + (*dp/(*dp-*dx)) * (*stx - *stp);
173: if (neP->bracket) {
174: if (TaoAbsDouble(*stp-stpc) < TaoAbsDouble(*stp-stpq)) {
175: stpf = stpc;
176: }
177: else {
178: stpf = stpq;
179: }
180: }
181: else {
182: if (TaoAbsDouble(*stp-stpc) > TaoAbsDouble(*stp-stpq)) {
183: stpf = stpc;
184: }
185: else {
186: stpf = stpq;
187: }
188: }
189: }
190: else {
191: // Case 4: A lower function value, derivatives of the
192: // same sign, and the magnitude of the derivative does
193: // not decrease. If the minimum is not bracketed, the step
194: // is either stpmin or stpmax, else the cubic step is taken.
196: neP->infoc = 4;
197: bound = 0;
198: if (neP->bracket) {
199: theta = 3*(*fp - *fy)/(*sty - *stp) + *dy + *dp;
200: s = TaoMax(TaoAbsDouble(theta),TaoAbsDouble(*dy));
201: s = TaoMax(s,TaoAbsDouble(*dp));
202: gamma1 = s*sqrt(pow(theta/s,2.0) - (*dy/s)*(*dp/s));
203: if (*stp > *sty) gamma1 = -gamma1;
204: p = (gamma1 - *dp) + theta;
205: q = ((gamma1 - *dp) + gamma1) + *dy;
206: r = p/q;
207: stpc = *stp + r*(*sty - *stp);
208: stpq = *stp + (*dp/(*dp-*dx)) * (*stx - *stp);
210: stpf = stpc;
211: }
212: else if (*stp > *stx) {
213: stpf = neP->stepmax;
214: }
215: else {
216: stpf = neP->stepmin;
217: }
218: }
219:
220: // Update the interval of uncertainty. This update does not
221: // depend on the new step or the case analysis above.
223: if (*fp > *fx) {
224: *sty = *stp;
225: *fy = *fp;
226: *dy = *dp;
227: }
228: else {
229: if (sgnd < 0.0) {
230: *sty = *stx;
231: *fy = *fx;
232: *dy = *dx;
233: }
234: *stx = *stp;
235: *fx = *fp;
236: *dx = *dp;
237: }
238:
239: // Compute the new step and safeguard it.
240: stpf = TaoMin(neP->stepmax,stpf);
241: stpf = TaoMax(neP->stepmin,stpf);
242: *stp = stpf;
243: if (neP->bracket && bound) {
244: if (*sty > *stx) {
245: *stp = TaoMin(*stx+0.66*(*sty-*stx),*stp);
246: }
247: else {
248: *stp = TaoMax(*stx+0.66*(*sty-*stx),*stp);
249: }
250: }
251: TaoFunctionReturn(0);
252: }
256: /* @
257: TaoApply_LineSearch - This routine performs a line search algorithm.
259: Input Parameters:
260: + tao - TAO_SOLVER context
261: . X - current iterate (on output X contains new iterate, X + step*S)
262: . S - search direction
263: . f - objective function evaluated at X
264: . G - gradient evaluated at X
265: . W - work vector
266: - step - initial estimate of step length
268: Output parameters:
269: + f - objective function evaluated at new iterate, X + step*S
270: . G - gradient evaluated at new iterate, X + step*S
271: . X - new iterate
272: . gnorm - 2-norm of G
273: - step - final step length
276: Info is set to one of:
277: + 0 - the line search succeeds; the sufficient decrease
278: condition and the directional derivative condition hold
280: negative number if an input parameter is invalid
281: . -1 - step < 0
282: . -2 - ftol < 0
283: . -3 - rtol < 0
284: . -4 - gtol < 0
285: . -5 - stepmin < 0
286: . -6 - stepmax < stepmin
287: - -7 - maxfev < 0
289: positive number > 1 if the line search otherwise terminates
290: + 2 - Relative width of the interval of uncertainty is
291: at most rtol.
292: . 3 - Maximum number of function evaluations (maxfev) has
293: been reached.
294: . 4 - Step is at the lower bound, stepmin.
295: . 5 - Step is at the upper bound, stepmax.
296: . 6 - Rounding errors may prevent further progress.
297: There may not be a step that satisfies the
298: sufficient decrease and curvature conditions.
299: Tolerances may be too small.
300: - 7 - Search direction is not a descent direction.
302: Notes:
303: This algorithm is taken from More' and Thuente, "Line search algorithms
304: with guaranteed sufficient decrease", Argonne National Laboratory,
305: Technical Report MCS-P330-1092.
307: Notes:
308: This routine is used within the following TAO unconstrained minimization
309: solvers: Newton linesearch (tao_nls), limited memory variable metric
310: (tao_lmvm), and conjugate gradient (tao_cg).
312: Level: advanced
314: .keywords: TAO_SOLVER, linesearch
315: @ */
316: static int TaoApply_LineSearch(TAO_SOLVER tao,TaoVec* X,TaoVec* G,TaoVec* S,TaoVec* W,double *f, double *f_full,
317: double *step,int *info2,void*ctx)
318: {
319: TAO_LINESEARCH *neP = (TAO_LINESEARCH *) tao->linectx;
320: double xtrapf = 4.0;
321: double finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym;
322: double dgx, dgy, dg, fx, fy, stx, sty, dgtest, ftest1=0.0;
323: int info, i, stage1;
325: #if defined(PETSC_USE_COMPLEX)
326: PetscScalar cdginit, cdg, cstep = 0.0;
327: #endif
329: TaoFunctionBegin;
330: /* neP->stepmin - lower bound for step */
331: /* neP->stepmax - upper bound for step */
332: /* neP->rtol - relative tolerance for an acceptable step */
333: /* neP->ftol - tolerance for sufficient decrease condition */
334: /* neP->gtol - tolerance for curvature condition */
335: /* neP->nfev - number of function evaluations */
336: /* neP->maxfev - maximum number of function evaluations */
338: *info2 = 0;
340: /* Check input parameters for errors */
341: if (*step < 0.0) {
342: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: step (%g) < 0\n",*step); CHKERRQ(info);
343: *info2 = -1;
344: }
345:
346: if (neP->ftol < 0.0) {
347: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: ftol (%g) < 0\n",neP->ftol); CHKERRQ(info);
348: *info2 = -2;
349: }
350:
351: if (neP->rtol < 0.0) {
352: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: rtol (%g) < 0\n",neP->rtol); CHKERRQ(info);
353: *info2 = -3;
354: }
356: if (neP->gtol < 0.0) {
357: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: gtol (%g) < 0\n",neP->gtol); CHKERRQ(info);
358: *info2 = -4;
359: }
361: if (neP->stepmin < 0.0) {
362: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: stepmin (%g) < 0\n",neP->stepmin); CHKERRQ(info);
363: *info2 = -5;
364: }
366: if (neP->stepmax < neP->stepmin) {
367: info = PetscInfo2(tao, "TaoApply_LineSearch: Line search error: stepmax (%g) < stepmin (%g)\n", neP->stepmax,neP->stepmin); CHKERRQ(info);
368: *info2 = -6;
369: }
370:
371: if (neP->maxfev < 0) {
372: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: maxfev (%d) < 0\n",neP->maxfev); CHKERRQ(info);
373: *info2 = -7;
374: }
376: if (TaoInfOrNaN(*f)) {
377: info = PetscInfo1(tao, "TaoApply_LineSearch: Line search error: function (%g) inf or nan\n", *f); CHKERRQ(info);
378: *info2 = -8;
379: }
381: /* Check that search direction is a descent direction */
383: #if defined(PETSC_USE_COMPLEX)
384: info = G->Dot(S,&cdginit);CHKERRQ(info); dginit = TaoReal(cdginit);
385: #else
386: info = G->Dot(S,&dginit);CHKERRQ(info);
387: #endif
389: if (TaoInfOrNaN(dginit)) {
390: info = PetscInfo1(tao,"TaoApply_LineSearch: Line search error: dginit (%g) inf or nan\n", dginit); CHKERRQ(info);
391: *info2 = -9;
392: }
394: if (dginit >= 0.0) {
395: info = PetscInfo(tao,"TaoApply_LineSearch:Search direction not a descent direction\n"); CHKERRQ(info);
396: *info2 = -10;
397: }
399: if (*info2) {
400: TaoFunctionReturn(0);
401: }
403: /* Initialization */
404: neP->bracket = 0;
405: *info2 = 0;
406: stage1 = 1;
407: finit = *f;
408: dgtest = neP->ftol * dginit;
409: width = neP->stepmax - neP->stepmin;
410: width1 = width * 2.0;
411: info = W->CopyFrom(X);CHKERRQ(info);
412: /* Variable dictionary:
413: stx, fx, dgx - the step, function, and derivative at the best step
414: sty, fy, dgy - the step, function, and derivative at the other endpoint
415: of the interval of uncertainty
416: step, f, dg - the step, function, and derivative at the current step */
418: stx = 0.0;
419: fx = finit;
420: dgx = dginit;
421: sty = 0.0;
422: fy = finit;
423: dgy = dginit;
424:
425: neP->nfev = 0;
426: for (i=0; i< neP->maxfev; i++) {
427: /* Set min and max steps to correspond to the interval of uncertainty */
428: if (neP->bracket) {
429: neP->stepmin = TaoMin(stx,sty);
430: neP->stepmax = TaoMax(stx,sty);
431: }
432: else {
433: neP->stepmin = stx;
434: neP->stepmax = *step + xtrapf * (*step - stx);
435: }
437: /* Force the step to be within the bounds */
438: *step = TaoMax(*step,neP->stepmin);
439: *step = TaoMin(*step,neP->stepmax);
440:
441: /* If an unusual termination is to occur, then let step be the lowest
442: point obtained thus far */
443: if (((neP->bracket) && (*step <= neP->stepmin || *step >= neP->stepmax)) ||
444: ((neP->bracket) && (neP->stepmax - neP->stepmin <= neP->rtol * neP->stepmax)) ||
445: (neP->nfev >= neP->maxfev - 1) || (neP->infoc == 0)) {
446: *step = stx;
447: }
449: #if defined(PETSC_USE_COMPLEX)
450: cstep = *step;
451: info = X->Waxpby(cstep,S,1.0,W);CHKERRQ(info);
452: #else
453: info = X->Waxpby(*step,S,1.0,W);CHKERRQ(info); /* X = W + step*S */
454: #endif
455: info = TaoComputeMeritFunctionGradient(tao,X,f,G);CHKERRQ(info);
456: if (0 == i) {
457: *f_full = *f;
458: }
460: neP->nfev++;
461: #if defined(PETSC_USE_COMPLEX)
462: info = G->Dot(S,&cdg);CHKERRQ(info); dg = TaoReal(cdg);
463: #else
464: info = G->Dot(S,&dg);CHKERRQ(info); /* dg = G^T S */
465: #endif
467: if (TaoInfOrNaN(*f) || TaoInfOrNaN(dg)) {
468: // User provided compute function generated Not-a-Number, assume
469: // domain violation and set function value and directional
470: // derivative to infinity.
471: *f = TAO_INFINITY;
472: dg = TAO_INFINITY;
473: }
475: ftest1 = finit + *step * dgtest;
477: /* Convergence testing */
478: if ((*f <= ftest1) && (TaoAbsDouble(dg) <= neP->gtol*(-dginit))) {
479: info = PetscInfo(tao, "TaoApply_LineSearch:Line search success: Sufficient decrease and directional deriv conditions hold\n"); CHKERRQ(info);
480: *info2 = 0;
481: break;
482: }
484: /* Checks for bad cases */
485: if (((neP->bracket) && (*step <= neP->stepmin||*step >= neP->stepmax)) || (!neP->infoc)) {
486: info = PetscInfo(tao,"TaoApply_LineSearch:Rounding errors may prevent further progress. May not be a step satisfying\n"); CHKERRQ(info);
487: info = PetscInfo(tao,"TaoApply_LineSearch:sufficient decrease and curvature conditions. Tolerances may be too small.\n"); CHKERRQ(info);
488: *info2 = 6; break;
489: }
490: if ((*step == neP->stepmax) && (*f <= ftest1) && (dg <= dgtest)) {
491: info = PetscInfo1(tao,"TaoApply_LineSearch:Step is at the upper bound, stepmax (%g)\n",neP->stepmax); CHKERRQ(info);
492: *info2 = 5; break;
493: }
494: if ((*step == neP->stepmin) && (*f >= ftest1) && (dg >= dgtest)) {
495: info = PetscInfo1(tao,"TaoApply_LineSearch:Step is at the lower bound, stepmin (%g)\n",neP->stepmin); CHKERRQ(info);
496: *info2 = 4; break;
497: }
498: if (neP->nfev >= neP->maxfev) {
499: info = PetscInfo2(tao,"TaoApply_LineSearch:Number of line search function evals (%d) > maximum (%d)\n",neP->nfev,neP->maxfev); CHKERRQ(info);
500: *info2 = 3; break;
501: }
502: if ((neP->bracket) && (neP->stepmax - neP->stepmin <= neP->rtol*neP->stepmax)){
503: info = PetscInfo1(tao,"TaoApply_LineSearch:Relative width of interval of uncertainty is at most rtol (%g)\n",neP->rtol); CHKERRQ(info);
504: *info2 = 2; break;
505: }
507: /* In the first stage, we seek a step for which the modified function
508: has a nonpositive value and nonnegative derivative */
509: if ((stage1) && (*f <= ftest1) && (dg >= dginit * TaoMin(neP->ftol, neP->gtol))) {
510: stage1 = 0;
511: }
513: /* A modified function is used to predict the step only if we
514: have not obtained a step for which the modified function has a
515: nonpositive function value and nonnegative derivative, and if a
516: lower function value has been obtained but the decrease is not
517: sufficient */
519: if ((stage1) && (*f <= fx) && (*f > ftest1)) {
520: fm = *f - *step * dgtest; /* Define modified function */
521: fxm = fx - stx * dgtest; /* and derivatives */
522: fym = fy - sty * dgtest;
523: dgm = dg - dgtest;
524: dgxm = dgx - dgtest;
525: dgym = dgy - dgtest;
526:
527: /* Update the interval of uncertainty and compute the new step */
528: info = TaoStep_LineSearch(tao,&stx,&fxm,&dgxm,&sty,&fym,&dgym,step,&fm,&dgm);CHKERRQ(info);
529:
530: fx = fxm + stx * dgtest; /* Reset the function and */
531: fy = fym + sty * dgtest; /* gradient values */
532: dgx = dgxm + dgtest;
533: dgy = dgym + dgtest;
534: }
535: else {
536: /* Update the interval of uncertainty and compute the new step */
537: info = TaoStep_LineSearch(tao,&stx,&fx,&dgx,&sty,&fy,&dgy,step,f,&dg);CHKERRQ(info);
538: }
539:
540: /* Force a sufficient decrease in the interval of uncertainty */
541: if (neP->bracket) {
542: if (TaoAbsDouble(sty - stx) >= 0.66 * width1) *step = stx + 0.5*(sty - stx);
543: width1 = width;
544: width = TaoAbsDouble(sty - stx);
545: }
546: }
547:
548: /* Finish computations */
549: info = PetscInfo2(tao,"TaoApply_LineSearch:%d function evals in line search, step = %10.4f\n",neP->nfev,*step); CHKERRQ(info);
550: TaoFunctionReturn(0);
551: }
555: /* @
556: TaoApply_BoundLineSearch - This routine performs a line search algorithm.
558: Input Parameters:
559: + tao - TAO_SOLVER context
560: . X - current iterate (on output X contains new iterate, X + step*S)
561: . S - search direction
562: . f - objective function evaluated at X
563: . G - gradient evaluated at X
564: . W - work vector
565: - step - initial estimate of step length
567: Output parameters:
568: + f - objective function evaluated at new iterate, X + step*S
569: . G - gradient evaluated at new iterate, X + step*S
570: . X - new iterate
571: . gnorm - 2-norm of G
572: - step - final step length
575: Info is set to one of:
576: + 0 - the line search succeeds; the sufficient decrease
577: condition and the directional derivative condition hold
579: negative number if an input parameter is invalid
580: . -1 - step < 0
581: . -2 - ftol < 0
582: . -3 - rtol < 0
583: . -4 - gtol < 0
584: . -5 - stepmin < 0
585: . -6 - stepmax < stepmin
586: - -7 - maxfev < 0
588: positive number > 1 if the line search otherwise terminates
589: + 2 - Relative width of the interval of uncertainty is
590: at most rtol.
591: . 3 - Maximum number of function evaluations (maxfev) has
592: been reached.
593: . 4 - Step is at the lower bound, stepmin.
594: . 5 - Step is at the upper bound, stepmax.
595: . 6 - Rounding errors may prevent further progress.
596: There may not be a step that satisfies the
597: sufficient decrease and curvature conditions.
598: Tolerances may be too small.
599: - 7 - Search direction is not a descent direction.
601: Notes:
602: This algorithm is is a modification of the algorithm by More' and
603: Thuente. The modifications concern bounds. This algorithm steps
604: in the direction passed into this routine. This point get
605: projected back into the feasible set. In the context of bound
606: constrained optimization, there may not be a point in the piecewise
607: linear path that satisfies the Wolfe conditions. When the active
608: set is changing, decrease in the objective function may be
609: sufficient to terminate this line search.
611: Note: Much of this coded is identical to the More' Thuente line search.
612: Variations to the code are commented.
614: Notes:
615: This routine is used within the following TAO bound constrained minimization
616: solvers: Newton linesearch (tao_tron) and limited memory variable metric
617: (tao_blmvm).
619: Level: advanced
621: .keywords: TAO_SOLVER, linesearch
622: @ */
623: static int TaoApply_BoundLineSearch(TAO_SOLVER tao,TaoVec* X,TaoVec* G,TaoVec* S,TaoVec* W,double *f, double *f_full,
624: double *step,int *info2,void*ctx)
625: {
626: TAO_LINESEARCH *neP = (TAO_LINESEARCH *) tao->linectx;
627: TaoVec *XL,*XU;
628: double xtrapf = 4.0;
629: double finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym;
630: double dgx, dgy, dg, fx, fy, stx, sty, dgtest, ftest1=0.0, ftest2=0.0;
631: double bstepmin1, bstepmin2, bstepmax;
632: double dg1, dg2;
633: int info, i, stage1;
635: #if defined(PETSC_USE_COMPLEX)
636: PetscScalar cdginit, cdg, cstep = 0.0;
637: #endif
639: TaoFunctionBegin;
640: /* neP->stepmin - lower bound for step */
641: /* neP->stepmax - upper bound for step */
642: /* neP->rtol - relative tolerance for an acceptable step */
643: /* neP->ftol - tolerance for sufficient decrease condition */
644: /* neP->gtol - tolerance for curvature condition */
645: /* neP->nfev - number of function evaluations */
646: /* neP->maxfev - maximum number of function evaluations */
648: /* Check input parameters for errors */
649: if (*step < 0.0) {
650: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Line search error: step (%g) < 0\n",*step); CHKERRQ(info);
651: *info2 = -1; TaoFunctionReturn(0);
652: }
653: else if (neP->ftol < 0.0) {
654: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Line search error: ftol (%g) < 0\n",neP->ftol); CHKERRQ(info);
655: *info2 = -2; TaoFunctionReturn(0);
656: }
657: else if (neP->rtol < 0.0) {
658: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Line search error: rtol (%g) < 0\n",neP->rtol); CHKERRQ(info);
659: *info2 = -3; TaoFunctionReturn(0);
660: }
661: else if (neP->gtol < 0.0) {
662: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Line search error: gtol (%g) < 0\n",neP->gtol); CHKERRQ(info);
663: *info2 = -4; TaoFunctionReturn(0);
664: }
665: else if (neP->stepmin < 0.0) {
666: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Line search error: stepmin (%g) < 0\n",neP->stepmin); CHKERRQ(info);
667: *info2 = -5; TaoFunctionReturn(0);
668: }
669: else if (neP->stepmax < neP->stepmin) {
670: info = PetscInfo2(tao,"TaoApply_BoundLineSearch:Line search error: stepmax (%g) < stepmin (%g)\n", neP->stepmax,neP->stepmin); CHKERRQ(info);
671: *info2 = -6; TaoFunctionReturn(0);
672: }
673: else if (neP->maxfev < 0.0) {
674: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Line search error: maxfev (%d) < 0\n",neP->maxfev); CHKERRQ(info);
675: *info2 = -7; TaoFunctionReturn(0);
676: }
678: /* Compute step length needed to make all variables equal a bound */
679: /* Compute the smallest steplength that will make one nonbinding */
680: /* variable equal the bound */
681: info = TaoGetVariableBounds(tao,&XL,&XU); CHKERRQ(info);
682: info = S->Negate(); CHKERRQ(info);
683: info = S->BoundGradientProjection(S,XL,X,XU); CHKERRQ(info);
684: info = S->Negate(); CHKERRQ(info);
685: info = X->StepBoundInfo(XL,XU,S,&bstepmin1,&bstepmin2,&bstepmax); CHKERRQ(info);
686: neP->stepmax=TaoMin(bstepmax,1.0e+15);
688: /* Check that search direction is a descent direction */
690: #if defined(PETSC_USE_COMPLEX)
691: info = G->Dot(S,&cdginit);CHKERRQ(info); dginit = TaoReal(cdginit);
692: #else
693: info = G->Dot(S,&dginit);CHKERRQ(info);
694: #endif
696: if (dginit >= 0.0) {
697: info = PetscInfo(tao,"TaoApply_BoundLineSearch:Search direction not a descent direction\n"); CHKERRQ(info);
698: *info2 = 7; TaoFunctionReturn(0);
699: }
701: /* Initialization */
702: neP->bracket = 0;
703: *info2 = 0;
704: stage1 = 1;
705: finit = *f;
706: dgtest = neP->ftol * dginit;
707: width = neP->stepmax - neP->stepmin;
708: width1 = width * 2.0;
709: info = W->CopyFrom(X);CHKERRQ(info);
710: /* Variable dictionary:
711: stx, fx, dgx - the step, function, and derivative at the best step
712: sty, fy, dgy - the step, function, and derivative at the other endpoint
713: of the interval of uncertainty
714: step, f, dg - the step, function, and derivative at the current step */
716: stx = 0.0;
717: fx = finit;
718: dgx = dginit;
719: sty = 0.0;
720: fy = finit;
721: dgy = dginit;
722:
723: neP->nfev = 0;
724: for (i=0; i< neP->maxfev; i++) {
725: /* Set min and max steps to correspond to the interval of uncertainty */
726: if (neP->bracket) {
727: neP->stepmin = TaoMin(stx,sty);
728: neP->stepmax = TaoMax(stx,sty);
729: } else {
730: neP->stepmin = stx;
731: neP->stepmax = *step + xtrapf * (*step - stx);
732: }
734: /* Force the step to be within the bounds */
735: *step = TaoMax(*step,neP->stepmin);
736: *step = TaoMin(*step,neP->stepmax);
738: /* If an unusual termination is to occur, then let step be the lowest
739: point obtained thus far */
740: if (((neP->bracket) && (*step <= neP->stepmin || *step >= neP->stepmax)) ||
741: ((neP->bracket) && (neP->stepmax - neP->stepmin <= neP->rtol * neP->stepmax)) ||
742: (neP->nfev >= neP->maxfev - 1) || (neP->infoc == 0)) {
743: *step = stx;
744: }
745:
746: #if defined(PETSC_USE_COMPLEX)
747: cstep = *step;
748: info = X->Waxpby(cstep,S,1.0,W);CHKERRQ(info);
749: #else
750: info = X->Waxpby(*step,S,1.0,W);CHKERRQ(info); /* X = W + step*S */
751: #endif
753: info=X->Median(XL,X,XU);CHKERRQ(info);
754: info = TaoComputeMeritFunctionGradient(tao,X,f,G);CHKERRQ(info);
755: if (0 == i) {
756: *f_full = *f;
757: }
759: neP->nfev++;
760: #if defined(PETSC_USE_COMPLEX)
761: info = G->Dot(S,&cdg);CHKERRQ(info); dg = TaoReal(cdg);
762: #else
763: info = G->Dot(W,&dg1);CHKERRQ(info); /* dg = G^T S */
764: info = G->Dot(X,&dg2);CHKERRQ(info); /* dg = G^T S */
765: dg = (dg2-dg1) / (*step);
766: #endif
768: if ((*f != *f) || (dg1 != dg1) || (dg2 != dg2) || (dg != dg)) {
769: // User provided compute function generated Not-a-Number, assume
770: // domain violation and set function value and directional
771: // derivative to infinity.
772: *f = TAO_INFINITY;
773: dg = TAO_INFINITY;
774: dg1 = TAO_INFINITY;
775: dg2 = TAO_INFINITY;
776: }
778: ftest1 = finit + (*step) * dgtest;
779: ftest2 = finit + (*step) * dgtest * neP->ftol; // Armijo
781: /* Convergence testing */
782: if ((*f <= ftest1) && (TaoAbsDouble(dg) <= neP->gtol*(-dginit))) {
783: info = PetscInfo(tao, "TaoApply_BoundLineSearch:Line search success: Sufficient decrease and directional deriv conditions hold\n"); CHKERRQ(info);
784: *info2 = 0;
785: break;
786: }
788: /* Check Armijo if beyond the first breakpoint */
789: if ((*f <= ftest2) && (*step >= bstepmin2)) {
790: info = PetscInfo(tao,"TaoApply_BoundLineSearch:Line search success: Sufficient decrease\n"); CHKERRQ(info);
791: *info2 = 0;
792: break;
793: }
795: /* Checks for bad cases */
796: if (((neP->bracket) && (*step <= neP->stepmin||*step >= neP->stepmax)) || (!neP->infoc)) {
797: info = PetscInfo(tao,"TaoApply_LineSearch:Rounding errors may prevent further progress. May not be a step satisfying\n"); CHKERRQ(info);
798: info = PetscInfo(tao,"TaoApply_BoundLineSearch:sufficient decrease and curvature conditions. Tolerances may be too small.\n"); CHKERRQ(info);
799: *info2 = 6; break;
800: }
801: if ((*step == neP->stepmax) && (*f <= ftest1) && (dg <= dgtest)) {
802: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Step is at the upper bound, stepmax (%g)\n",neP->stepmax); CHKERRQ(info);
803: *info2 = 5; break;
804: }
805: if ((*step == neP->stepmin) && (*f >= ftest1) && (dg >= dgtest)) {
806: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Step is at the lower bound, stepmin (%g)\n",neP->stepmin); CHKERRQ(info);
807: *info2 = 4; break;
808: }
809: if (neP->nfev >= neP->maxfev) {
810: info = PetscInfo2(tao,"TaoApply_BoundLineSearch:Number of line search function evals (%d) > maximum (%d)\n",neP->nfev,neP->maxfev); CHKERRQ(info);
811: *info2 = 3; break;
812: }
813: if ((neP->bracket) && (neP->stepmax - neP->stepmin <= neP->rtol*neP->stepmax)){
814: info = PetscInfo1(tao,"TaoApply_BoundLineSearch:Relative width of interval of uncertainty is at most rtol (%g)\n",neP->rtol); CHKERRQ(info);
815: *info2 = 2; break;
816: }
818: /* In the first stage, we seek a step for which the modified function
819: has a nonpositive value and nonnegative derivative */
820: if ((stage1) && (*f <= ftest1) && (dg >= dginit * TaoMin(neP->ftol, neP->gtol)))
821: stage1 = 0;
823: /* A modified function is used to predict the step only if we
824: have not obtained a step for which the modified function has a
825: nonpositive function value and nonnegative derivative, and if a
826: lower function value has been obtained but the decrease is not
827: sufficient */
829: if ((stage1) && (*f <= fx) && (*f > ftest1)) {
830: fm = *f - *step * dgtest; /* Define modified function */
831: fxm = fx - stx * dgtest; /* and derivatives */
832: fym = fy - sty * dgtest;
833: dgm = dg - dgtest;
834: dgxm = dgx - dgtest;
835: dgym = dgy - dgtest;
837: /* Update the interval of uncertainty and compute the new step */
838: info = TaoStep_LineSearch(tao,&stx,&fxm,&dgxm,&sty,&fym,&dgym,step,&fm,&dgm);CHKERRQ(info);
840: fx = fxm + stx * dgtest; /* Reset the function and */
841: fy = fym + sty * dgtest; /* gradient values */
842: dgx = dgxm + dgtest;
843: dgy = dgym + dgtest;
844: } else {
845: /* Update the interval of uncertainty and compute the new step */
846: info = TaoStep_LineSearch(tao,&stx,&fx,&dgx,&sty,&fy,&dgy,step,f,&dg);CHKERRQ(info);
847: }
849: /* Force a sufficient decrease in the interval of uncertainty */
850: if (neP->bracket) {
851: if (TaoAbsDouble(sty - stx) >= 0.66 * width1) *step = stx + 0.5*(sty - stx);
852: width1 = width;
853: width = TaoAbsDouble(sty - stx);
854: }
855: }
857: /* Finish computations */
858: info = PetscInfo2(tao,"TaoApply_BoundLineSearch:%d function evals in line search, step = %10.4e\n",neP->nfev,*step); CHKERRQ(info);
859: TaoFunctionReturn(0);
860: }
864: static int TaoDestroy_LineSearch(TAO_SOLVER tao, void *ctx)
865: {
866: int info;
868: TaoFunctionBegin;
869: info = TaoFree(tao->linectx);CHKERRQ(info);
870: TaoFunctionReturn(0);
871: }
875: static int TaoSetOptions_LineSearch(TAO_SOLVER tao, void*linectx)
876: {
877: TAO_LINESEARCH *ctx = (TAO_LINESEARCH *)linectx;
878: int info;
880: TaoFunctionBegin;
881: info = TaoOptionsHead("More-Thuente line line search options for unconstrained minimization");CHKERRQ(info);
882: info = TaoOptionInt("-tao_ls_maxfev","max function evals in line search","",ctx->maxfev,&ctx->maxfev,0);CHKERRQ(info);
883: info = TaoOptionDouble("-tao_ls_ftol","tol for sufficient decrease","",ctx->ftol,&ctx->ftol,0);CHKERRQ(info);
884: info = TaoOptionDouble("-tao_ls_gtol","tol for curvature condition","",ctx->gtol,&ctx->gtol,0);CHKERRQ(info);
885: info = TaoOptionDouble("-tao_ls_rtol","relative tol for acceptable step","",ctx->rtol,&ctx->rtol,0);CHKERRQ(info);
886: info = TaoOptionDouble("-tao_ls_stepmin","lower bound for step","",ctx->stepmin,&ctx->stepmin,0);CHKERRQ(info);
887: info = TaoOptionDouble("-tao_ls_stepmax","upper bound for step","",ctx->stepmax,&ctx->stepmax,0);CHKERRQ(info);
888: info = TaoOptionsTail();CHKERRQ(info);
889: TaoFunctionReturn(0);
890: }
894: static int TaoView_LineSearch(TAO_SOLVER tao, void*ctx)
895: {
896: TAO_LINESEARCH *ls = (TAO_LINESEARCH *)ctx;
897: int info;
899: TaoFunctionBegin;
900: info=TaoPrintInt(tao," More'-Thuente Line search: maxf=%d,",ls->maxfev);CHKERRQ(info);
901: info=TaoPrintDouble(tao," ftol=%g,",ls->ftol);CHKERRQ(info);
902: info=TaoPrintDouble(tao," rtol=%g,",ls->rtol);CHKERRQ(info);
903: info=TaoPrintDouble(tao," gtol=%g\n",ls->gtol);CHKERRQ(info);
904: TaoFunctionReturn(0);
905: }
909: /*@C
910: TaoCreateMoreThuenteLineSearch - Create a line search
912: Input Parameters:
913: + tao - TAO_SOLVER context
914: . fftol - the sufficient descent parameter , greater than 0.
915: - ggtol - the curvature tolerance, greater than 0, less than 1.
918: Note:
919: If either fftol or ggtol is 0, default parameters will be used.
921: Note:
922: This algorithm is taken from More' and Thuente, "Line search algorithms
923: with guaranteed sufficient decrease", Argonne National Laboratory,
924: Technical Report MCS-P330-1092.
926: Note:
927: This line search enforces the strong Wolfe conditions for unconstrained
928: optimization. This routine is used within the following TAO unconstrained
929: minimization solvers: Newton linesearch (tao_nls), limited memory variable
930: metric (tao_lmvm), and nonlinear conjugate gradient methods.
932: Level: developer
934: .keywords: TAO_SOLVER, linesearch
935: @*/
936: int TaoCreateMoreThuenteLineSearch(TAO_SOLVER tao, double fftol, double ggtol)
937: {
938: int info;
939: TAO_LINESEARCH *neP;
941: TaoFunctionBegin;
943: info = TaoNew(TAO_LINESEARCH,&neP); CHKERRQ(info);
944: info = PetscLogObjectMemory(tao,sizeof(TAO_LINESEARCH)); CHKERRQ(info);
946: if (fftol<=0) fftol=0.0001;
947: if (ggtol<=0) ggtol=0.9;
949: neP->ftol = fftol;
950: neP->rtol = 1.0e-10;
951: neP->gtol = ggtol;
952: neP->stepmin = 1.0e-20;
953: neP->stepmax = 1.0e+20;
954: neP->nfev = 0;
955: neP->bracket = 0;
956: neP->infoc = 1;
957: neP->maxfev = 30;
959: info = TaoSetLineSearch(tao,0,
960: TaoSetOptions_LineSearch,
961: TaoApply_LineSearch,
962: TaoView_LineSearch,
963: TaoDestroy_LineSearch,
964: (void *) neP);CHKERRQ(info);
966: TaoFunctionReturn(0);
967: }
971: /*@C
972: TaoCreateMoreThuenteBoundLineSearch - Create a line search
974: Input Parameters:
975: + tao - TAO_SOLVER context
976: . fftol - the sufficient descent parameter , greater than 0.
977: - ggtol - the curvature tolerance, greater than 0, less than 1.
980: Note:
981: If either fftol or ggtol is 0, default parameters will be used.
983: Note:
984: This algorithm is is a modification of the algorithm by More' and Thuente.
985: The modifications concern bounds. This algorithm steps in the direction
986: passed into this routine. This point get projected back into the feasible set.
987: It tries to satisfy the Wolfe conditions, but in the context of bound constrained
988: optimization, there may not be a point in the piecewise linear
989: path that satisfies the Wolfe conditions. When the active set
990: is changing, decrease in the objective function may be sufficient
991: to terminate this line search.
993: Note:
994: This routine is used within the following TAO bound constrained
995: minimization solvers: Newton trust region (tao_tron) and limited memory variable
996: metric (tao_blmvm).
998: Level: developer
1000: .keywords: TAO_SOLVER, linesearch
1001: @*/
1002: int TaoCreateMoreThuenteBoundLineSearch(TAO_SOLVER tao, double fftol, double ggtol)
1003: {
1004: int info;
1005: TAO_LINESEARCH *neP;
1007: TaoFunctionBegin;
1009: info = TaoNew(TAO_LINESEARCH,&neP); CHKERRQ(info);
1010: info = PetscLogObjectMemory(tao,sizeof(TAO_LINESEARCH)); CHKERRQ(info);
1012: if (fftol<=0) fftol=0.0001;
1013: if (ggtol<=0) ggtol=0.9;
1015: neP->ftol = fftol;
1016: neP->rtol = 1.0e-10;
1017: neP->gtol = ggtol;
1018: neP->stepmin = 1.0e-20;
1019: neP->stepmax = 1.0e+20;
1020: neP->nfev = 0;
1021: neP->bracket = 0;
1022: neP->infoc = 1;
1023: neP->maxfev = 30;
1025: info = TaoSetLineSearch(tao,0,
1026: TaoSetOptions_LineSearch,
1027: TaoApply_BoundLineSearch,
1028: TaoView_LineSearch,
1029: TaoDestroy_LineSearch,
1030: (void *) neP);CHKERRQ(info);
1032: TaoFunctionReturn(0);
1033: }