xref: /petsc/src/tao/linesearch/impls/morethuente/morethuente.c (revision 3fc23ce8277e00a65c4b1be6e2b4535b56ba7978)
1 #include <petsc/private/taolinesearchimpl.h>
2 #include <../src/tao/linesearch/impls/morethuente/morethuente.h>
3 
4 /*
5    This algorithm is taken from More' and Thuente, "Line search algorithms
6    with guaranteed sufficient decrease", Argonne National Laboratory,
7    Technical Report MCS-P330-1092.
8 */
9 
10 static PetscErrorCode Tao_mcstep(TaoLineSearch ls,PetscReal *stx,PetscReal *fx,PetscReal *dx,PetscReal *sty,PetscReal *fy,PetscReal *dy,PetscReal *stp,PetscReal *fp,PetscReal *dp);
11 
12 static PetscErrorCode TaoLineSearchDestroy_MT(TaoLineSearch ls)
13 {
14   TaoLineSearch_MT *mt = (TaoLineSearch_MT*)(ls->data);
15 
16   PetscFunctionBegin;
17   PetscCall(PetscObjectDereference((PetscObject)mt->x));
18   PetscCall(VecDestroy(&mt->work));
19   PetscCall(PetscFree(ls->data));
20   PetscFunctionReturn(0);
21 }
22 
23 static PetscErrorCode TaoLineSearchSetFromOptions_MT(PetscOptionItems *PetscOptionsObject,TaoLineSearch ls)
24 {
25   PetscFunctionBegin;
26   PetscFunctionReturn(0);
27 }
28 
29 static PetscErrorCode TaoLineSearchMonitor_MT(TaoLineSearch ls)
30 {
31   TaoLineSearch_MT *mt = (TaoLineSearch_MT*)ls->data;
32 
33   PetscFunctionBegin;
34   PetscCall(PetscViewerASCIIPrintf(ls->viewer, "stx: %g, fx: %g, dgx: %g\n", (double)mt->stx, (double)mt->fx, (double)mt->dgx));
35   PetscCall(PetscViewerASCIIPrintf(ls->viewer, "sty: %g, fy: %g, dgy: %g\n", (double)mt->sty, (double)mt->fy, (double)mt->dgy));
36   PetscFunctionReturn(0);
37 }
38 
39 static PetscErrorCode TaoLineSearchApply_MT(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s)
40 {
41   TaoLineSearch_MT *mt = (TaoLineSearch_MT*)(ls->data);
42   PetscReal        xtrapf = 4.0;
43   PetscReal        finit, width, width1, dginit, fm, fxm, fym, dgm, dgxm, dgym;
44   PetscReal        dgx, dgy, dg, dg2, fx, fy, stx, sty, dgtest;
45   PetscReal        ftest1=0.0, ftest2=0.0;
46   PetscInt         i, stage1,n1,n2,nn1,nn2;
47   PetscReal        bstepmin1, bstepmin2, bstepmax;
48   PetscBool        g_computed = PETSC_FALSE; /* to prevent extra gradient computation */
49 
50   PetscFunctionBegin;
51   ls->reason = TAOLINESEARCH_CONTINUE_ITERATING;
52   PetscCall(TaoLineSearchMonitor(ls, 0, *f, 0.0));
53   /* Check work vector */
54   if (!mt->work) {
55     PetscCall(VecDuplicate(x,&mt->work));
56     mt->x = x;
57     PetscCall(PetscObjectReference((PetscObject)mt->x));
58   } else if (x != mt->x) {
59     PetscCall(VecDestroy(&mt->work));
60     PetscCall(VecDuplicate(x,&mt->work));
61     PetscCall(PetscObjectDereference((PetscObject)mt->x));
62     mt->x = x;
63     PetscCall(PetscObjectReference((PetscObject)mt->x));
64   }
65 
66   if (ls->bounded) {
67     /* Compute step length needed to make all variables equal a bound */
68     /* Compute the smallest steplength that will make one nonbinding variable
69      equal the bound */
70     PetscCall(VecGetLocalSize(ls->upper,&n1));
71     PetscCall(VecGetLocalSize(mt->x, &n2));
72     PetscCall(VecGetSize(ls->upper,&nn1));
73     PetscCall(VecGetSize(mt->x,&nn2));
74     PetscCheck(n1 == n2 && nn1 == nn2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Variable vector not compatible with bounds vector");
75     PetscCall(VecScale(s,-1.0));
76     PetscCall(VecBoundGradientProjection(s,x,ls->lower,ls->upper,s));
77     PetscCall(VecScale(s,-1.0));
78     PetscCall(VecStepBoundInfo(x,s,ls->lower,ls->upper,&bstepmin1,&bstepmin2,&bstepmax));
79     ls->stepmax = PetscMin(bstepmax,1.0e15);
80   }
81 
82   PetscCall(VecDot(g,s,&dginit));
83   if (PetscIsInfOrNanReal(dginit)) {
84     PetscCall(PetscInfo(ls,"Initial Line Search step * g is Inf or Nan (%g)\n",(double)dginit));
85     ls->reason = TAOLINESEARCH_FAILED_INFORNAN;
86     PetscFunctionReturn(0);
87   }
88   if (dginit >= 0.0) {
89     PetscCall(PetscInfo(ls,"Initial Line Search step * g is not descent direction (%g)\n",(double)dginit));
90     ls->reason = TAOLINESEARCH_FAILED_ASCENT;
91     PetscFunctionReturn(0);
92   }
93 
94   /* Initialization */
95   mt->bracket = 0;
96   stage1 = 1;
97   finit = *f;
98   dgtest = ls->ftol * dginit;
99   width = ls->stepmax - ls->stepmin;
100   width1 = width * 2.0;
101   PetscCall(VecCopy(x,mt->work));
102   /* Variable dictionary:
103    stx, fx, dgx - the step, function, and derivative at the best step
104    sty, fy, dgy - the step, function, and derivative at the other endpoint
105    of the interval of uncertainty
106    step, f, dg - the step, function, and derivative at the current step */
107 
108   stx = 0.0;
109   fx  = finit;
110   dgx = dginit;
111   sty = 0.0;
112   fy  = finit;
113   dgy = dginit;
114 
115   ls->step = ls->initstep;
116   for (i=0; i<ls->max_funcs; i++) {
117     /* Set min and max steps to correspond to the interval of uncertainty */
118     if (mt->bracket) {
119       ls->stepmin = PetscMin(stx,sty);
120       ls->stepmax = PetscMax(stx,sty);
121     } else {
122       ls->stepmin = stx;
123       ls->stepmax = ls->step + xtrapf * (ls->step - stx);
124     }
125 
126     /* Force the step to be within the bounds */
127     ls->step = PetscMax(ls->step,ls->stepmin);
128     ls->step = PetscMin(ls->step,ls->stepmax);
129 
130     /* If an unusual termination is to occur, then let step be the lowest
131      point obtained thus far */
132     if (stx != 0 && ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol * ls->stepmax)) ||
133                      (ls->nfeval + ls->nfgeval >= ls->max_funcs - 1) || mt->infoc == 0)) ls->step = stx;
134 
135     PetscCall(VecWAXPY(mt->work,ls->step,s,x));   /* W = X + step*S */
136 
137     if (ls->bounded) {
138       PetscCall(VecMedian(ls->lower, mt->work, ls->upper, mt->work));
139     }
140     if (ls->usegts) {
141       PetscCall(TaoLineSearchComputeObjectiveAndGTS(ls,mt->work,f,&dg));
142       g_computed = PETSC_FALSE;
143     } else {
144       PetscCall(TaoLineSearchComputeObjectiveAndGradient(ls,mt->work,f,g));
145       g_computed = PETSC_TRUE;
146       if (ls->bounded) {
147         PetscCall(VecDot(g,x,&dg));
148         PetscCall(VecDot(g,mt->work,&dg2));
149         dg = (dg2 - dg)/ls->step;
150       } else {
151         PetscCall(VecDot(g,s,&dg));
152       }
153     }
154 
155     /* update bracketing parameters in the MT context for printouts in monitor */
156     mt->stx = stx;
157     mt->fx = fx;
158     mt->dgx = dgx;
159     mt->sty = sty;
160     mt->fy = fy;
161     mt->dgy = dgy;
162     PetscCall(TaoLineSearchMonitor(ls, i+1, *f, ls->step));
163 
164     if (i == 0) ls->f_fullstep = *f;
165 
166     if (PetscIsInfOrNanReal(*f) || PetscIsInfOrNanReal(dg)) {
167       /* User provided compute function generated Not-a-Number, assume
168        domain violation and set function value and directional
169        derivative to infinity. */
170       *f = PETSC_INFINITY;
171       dg = PETSC_INFINITY;
172     }
173 
174     ftest1 = finit + ls->step * dgtest;
175     if (ls->bounded) ftest2 = finit + ls->step * dgtest * ls->ftol;
176 
177     /* Convergence testing */
178     if ((*f - ftest1 <= PETSC_SMALL * PetscAbsReal(finit)) &&  (PetscAbsReal(dg) + ls->gtol*dginit <= 0.0)) {
179       PetscCall(PetscInfo(ls, "Line search success: Sufficient decrease and directional deriv conditions hold\n"));
180       ls->reason = TAOLINESEARCH_SUCCESS;
181       break;
182     }
183 
184     /* Check Armijo if beyond the first breakpoint */
185     if (ls->bounded && *f <= ftest2 && ls->step >= bstepmin2) {
186       PetscCall(PetscInfo(ls,"Line search success: Sufficient decrease.\n"));
187       ls->reason = TAOLINESEARCH_SUCCESS;
188       break;
189     }
190 
191     /* Checks for bad cases */
192     if ((mt->bracket && (ls->step <= ls->stepmin || ls->step >= ls->stepmax)) || !mt->infoc) {
193       PetscCall(PetscInfo(ls,"Rounding errors may prevent further progress. May not be a step satisfying\nsufficient decrease and curvature conditions. Tolerances may be too small.\n"));
194       ls->reason = TAOLINESEARCH_HALTED_OTHER;
195       break;
196     }
197     if (ls->step == ls->stepmax && *f <= ftest1 && dg <= dgtest) {
198       PetscCall(PetscInfo(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax));
199       ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND;
200       break;
201     }
202     if (ls->step == ls->stepmin && *f >= ftest1 && dg >= dgtest) {
203       PetscCall(PetscInfo(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin));
204       ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND;
205       break;
206     }
207     if (mt->bracket && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)) {
208       PetscCall(PetscInfo(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol));
209       ls->reason = TAOLINESEARCH_HALTED_RTOL;
210       break;
211     }
212 
213     /* In the first stage, we seek a step for which the modified function
214      has a nonpositive value and nonnegative derivative */
215     if (stage1 && *f <= ftest1 && dg >= dginit * PetscMin(ls->ftol, ls->gtol)) stage1 = 0;
216 
217     /* A modified function is used to predict the step only if we
218      have not obtained a step for which the modified function has a
219      nonpositive function value and nonnegative derivative, and if a
220      lower function value has been obtained but the decrease is not
221      sufficient */
222 
223     if (stage1 && *f <= fx && *f > ftest1) {
224       fm   = *f - ls->step * dgtest;    /* Define modified function */
225       fxm  = fx - stx * dgtest;         /* and derivatives */
226       fym  = fy - sty * dgtest;
227       dgm  = dg - dgtest;
228       dgxm = dgx - dgtest;
229       dgym = dgy - dgtest;
230 
231       /* if (dgxm * (ls->step - stx) >= 0.0) */
232       /* Update the interval of uncertainty and compute the new step */
233       PetscCall(Tao_mcstep(ls,&stx,&fxm,&dgxm,&sty,&fym,&dgym,&ls->step,&fm,&dgm));
234 
235       fx  = fxm + stx * dgtest; /* Reset the function and */
236       fy  = fym + sty * dgtest; /* gradient values */
237       dgx = dgxm + dgtest;
238       dgy = dgym + dgtest;
239     } else {
240       /* Update the interval of uncertainty and compute the new step */
241       PetscCall(Tao_mcstep(ls,&stx,&fx,&dgx,&sty,&fy,&dgy,&ls->step,f,&dg));
242     }
243 
244     /* Force a sufficient decrease in the interval of uncertainty */
245     if (mt->bracket) {
246       if (PetscAbsReal(sty - stx) >= 0.66 * width1) ls->step = stx + 0.5*(sty - stx);
247       width1 = width;
248       width = PetscAbsReal(sty - stx);
249     }
250   }
251   if (ls->nfeval + ls->nfgeval > ls->max_funcs) {
252     PetscCall(PetscInfo(ls,"Number of line search function evals (%" PetscInt_FMT ") > maximum (%" PetscInt_FMT ")\n",ls->nfeval+ls->nfgeval,ls->max_funcs));
253     ls->reason = TAOLINESEARCH_HALTED_MAXFCN;
254   }
255 
256   /* Finish computations */
257   PetscCall(PetscInfo(ls,"%" PetscInt_FMT " function evals in line search, step = %g\n",ls->nfeval+ls->nfgeval,(double)ls->step));
258 
259   /* Set new solution vector and compute gradient if needed */
260   PetscCall(VecCopy(mt->work,x));
261   if (!g_computed) {
262     PetscCall(TaoLineSearchComputeGradient(ls,mt->work,g));
263   }
264   PetscFunctionReturn(0);
265 }
266 
267 /*MC
268    TAOLINESEARCHMT - Line-search type with cubic interpolation that satisfies both the sufficient decrease and
269    curvature conditions. This method can take step lengths greater than 1.
270 
271    More-Thuente line-search can be selected with "-tao_ls_type more-thuente".
272 
273    References:
274 .  * - JORGE J. MORE AND DAVID J. THUENTE, LINE SEARCH ALGORITHMS WITH GUARANTEED SUFFICIENT DECREASE.
275           ACM Trans. Math. Software 20, no. 3 (1994): 286-307.
276 
277    Level: developer
278 
279 .seealso: TaoLineSearchCreate(), TaoLineSearchSetType(), TaoLineSearchApply()
280 
281 .keywords: Tao, linesearch
282 M*/
283 PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_MT(TaoLineSearch ls)
284 {
285   TaoLineSearch_MT *ctx;
286 
287   PetscFunctionBegin;
288   PetscValidHeaderSpecific(ls,TAOLINESEARCH_CLASSID,1);
289   PetscCall(PetscNewLog(ls,&ctx));
290   ctx->bracket = 0;
291   ctx->infoc = 1;
292   ls->data = (void*)ctx;
293   ls->initstep = 1.0;
294   ls->ops->setup = NULL;
295   ls->ops->reset = NULL;
296   ls->ops->apply = TaoLineSearchApply_MT;
297   ls->ops->destroy = TaoLineSearchDestroy_MT;
298   ls->ops->setfromoptions = TaoLineSearchSetFromOptions_MT;
299   ls->ops->monitor = TaoLineSearchMonitor_MT;
300   PetscFunctionReturn(0);
301 }
302 
303 /*
304      The subroutine mcstep is taken from the work of Jorge Nocedal.
305      this is a variant of More' and Thuente's routine.
306 
307      subroutine mcstep
308 
309      the purpose of mcstep is to compute a safeguarded step for
310      a linesearch and to update an interval of uncertainty for
311      a minimizer of the function.
312 
313      the parameter stx contains the step with the least function
314      value. the parameter stp contains the current step. it is
315      assumed that the derivative at stx is negative in the
316      direction of the step. if bracket is set true then a
317      minimizer has been bracketed in an interval of uncertainty
318      with endpoints stx and sty.
319 
320      the subroutine statement is
321 
322      subroutine mcstep(stx,fx,dx,sty,fy,dy,stp,fp,dp,bracket,
323                        stpmin,stpmax,info)
324 
325      where
326 
327        stx, fx, and dx are variables which specify the step,
328          the function, and the derivative at the best step obtained
329          so far. The derivative must be negative in the direction
330          of the step, that is, dx and stp-stx must have opposite
331          signs. On output these parameters are updated appropriately.
332 
333        sty, fy, and dy are variables which specify the step,
334          the function, and the derivative at the other endpoint of
335          the interval of uncertainty. On output these parameters are
336          updated appropriately.
337 
338        stp, fp, and dp are variables which specify the step,
339          the function, and the derivative at the current step.
340          If bracket is set true then on input stp must be
341          between stx and sty. On output stp is set to the new step.
342 
343        bracket is a logical variable which specifies if a minimizer
344          has been bracketed.  If the minimizer has not been bracketed
345          then on input bracket must be set false.  If the minimizer
346          is bracketed then on output bracket is set true.
347 
348        stpmin and stpmax are input variables which specify lower
349          and upper bounds for the step.
350 
351        info is an integer output variable set as follows:
352          if info = 1,2,3,4,5, then the step has been computed
353          according to one of the five cases below. otherwise
354          info = 0, and this indicates improper input parameters.
355 
356      subprograms called
357 
358        fortran-supplied ... abs,max,min,sqrt
359 
360      argonne national laboratory. minpack project. june 1983
361      jorge j. more', david j. thuente
362 
363 */
364 
365 static PetscErrorCode Tao_mcstep(TaoLineSearch ls,PetscReal *stx,PetscReal *fx,PetscReal *dx,PetscReal *sty,PetscReal *fy,PetscReal *dy,PetscReal *stp,PetscReal *fp,PetscReal *dp)
366 {
367   TaoLineSearch_MT *mtP = (TaoLineSearch_MT *) ls->data;
368   PetscReal        gamma1, p, q, r, s, sgnd, stpc, stpf, stpq, theta;
369   PetscInt         bound;
370 
371   PetscFunctionBegin;
372   /* Check the input parameters for errors */
373   mtP->infoc = 0;
374   PetscCheck(!mtP->bracket || (*stp > PetscMin(*stx,*sty) && *stp < PetscMax(*stx,*sty)),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"bad stp in bracket");
375   PetscCheck(*dx * (*stp-*stx) < 0.0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"dx * (stp-stx) >= 0.0");
376   PetscCheck(ls->stepmax >= ls->stepmin,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"stepmax > stepmin");
377 
378   /* Determine if the derivatives have opposite sign */
379   sgnd = *dp * (*dx / PetscAbsReal(*dx));
380 
381   if (*fp > *fx) {
382     /* Case 1: a higher function value.
383      The minimum is bracketed. If the cubic step is closer
384      to stx than the quadratic step, the cubic step is taken,
385      else the average of the cubic and quadratic steps is taken. */
386 
387     mtP->infoc = 1;
388     bound = 1;
389     theta = 3 * (*fx - *fp) / (*stp - *stx) + *dx + *dp;
390     s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx));
391     s = PetscMax(s,PetscAbsReal(*dp));
392     gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s));
393     if (*stp < *stx) gamma1 = -gamma1;
394     /* Can p be 0?  Check */
395     p = (gamma1 - *dx) + theta;
396     q = ((gamma1 - *dx) + gamma1) + *dp;
397     r = p/q;
398     stpc = *stx + r*(*stp - *stx);
399     stpq = *stx + ((*dx/((*fx-*fp)/(*stp-*stx)+*dx))*0.5) * (*stp - *stx);
400 
401     if (PetscAbsReal(stpc-*stx) < PetscAbsReal(stpq-*stx)) stpf = stpc;
402     else stpf = stpc + 0.5*(stpq - stpc);
403     mtP->bracket = 1;
404   } else if (sgnd < 0.0) {
405     /* Case 2: A lower function value and derivatives of
406      opposite sign. The minimum is bracketed. If the cubic
407      step is closer to stx than the quadratic (secant) step,
408      the cubic step is taken, else the quadratic step is taken. */
409 
410     mtP->infoc = 2;
411     bound = 0;
412     theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp;
413     s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx));
414     s = PetscMax(s,PetscAbsReal(*dp));
415     gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s));
416     if (*stp > *stx) gamma1 = -gamma1;
417     p = (gamma1 - *dp) + theta;
418     q = ((gamma1 - *dp) + gamma1) + *dx;
419     r = p/q;
420     stpc = *stp + r*(*stx - *stp);
421     stpq = *stp + (*dp/(*dp-*dx))*(*stx - *stp);
422 
423     if (PetscAbsReal(stpc-*stp) > PetscAbsReal(stpq-*stp)) stpf = stpc;
424     else stpf = stpq;
425     mtP->bracket = 1;
426   } else if (PetscAbsReal(*dp) < PetscAbsReal(*dx)) {
427     /* Case 3: A lower function value, derivatives of the
428      same sign, and the magnitude of the derivative decreases.
429      The cubic step is only used if the cubic tends to infinity
430      in the direction of the step or if the minimum of the cubic
431      is beyond stp. Otherwise the cubic step is defined to be
432      either stepmin or stepmax. The quadratic (secant) step is also
433      computed and if the minimum is bracketed then the step
434      closest to stx is taken, else the step farthest away is taken. */
435 
436     mtP->infoc = 3;
437     bound = 1;
438     theta = 3*(*fx - *fp)/(*stp - *stx) + *dx + *dp;
439     s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dx));
440     s = PetscMax(s,PetscAbsReal(*dp));
441 
442     /* The case gamma1 = 0 only arises if the cubic does not tend
443        to infinity in the direction of the step. */
444     gamma1 = s*PetscSqrtScalar(PetscMax(0.0,PetscPowScalar(theta/s,2.0) - (*dx/s)*(*dp/s)));
445     if (*stp > *stx) gamma1 = -gamma1;
446     p = (gamma1 - *dp) + theta;
447     q = (gamma1 + (*dx - *dp)) + gamma1;
448     r = p/q;
449     if (r < 0.0 && gamma1 != 0.0) stpc = *stp + r*(*stx - *stp);
450     else if (*stp > *stx)         stpc = ls->stepmax;
451     else                          stpc = ls->stepmin;
452     stpq = *stp + (*dp/(*dp-*dx)) * (*stx - *stp);
453 
454     if (mtP->bracket) {
455       if (PetscAbsReal(*stp-stpc) < PetscAbsReal(*stp-stpq)) stpf = stpc;
456       else stpf = stpq;
457     } else {
458       if (PetscAbsReal(*stp-stpc) > PetscAbsReal(*stp-stpq)) stpf = stpc;
459       else stpf = stpq;
460     }
461   } else {
462     /* Case 4: A lower function value, derivatives of the
463        same sign, and the magnitude of the derivative does
464        not decrease. If the minimum is not bracketed, the step
465        is either stpmin or stpmax, else the cubic step is taken. */
466 
467     mtP->infoc = 4;
468     bound = 0;
469     if (mtP->bracket) {
470       theta = 3*(*fp - *fy)/(*sty - *stp) + *dy + *dp;
471       s = PetscMax(PetscAbsReal(theta),PetscAbsReal(*dy));
472       s = PetscMax(s,PetscAbsReal(*dp));
473       gamma1 = s*PetscSqrtScalar(PetscPowScalar(theta/s,2.0) - (*dy/s)*(*dp/s));
474       if (*stp > *sty) gamma1 = -gamma1;
475       p = (gamma1 - *dp) + theta;
476       q = ((gamma1 - *dp) + gamma1) + *dy;
477       r = p/q;
478       stpc = *stp + r*(*sty - *stp);
479       stpf = stpc;
480     } else if (*stp > *stx) {
481       stpf = ls->stepmax;
482     } else {
483       stpf = ls->stepmin;
484     }
485   }
486 
487   /* Update the interval of uncertainty.  This update does not
488      depend on the new step or the case analysis above. */
489 
490   if (*fp > *fx) {
491     *sty = *stp;
492     *fy = *fp;
493     *dy = *dp;
494   } else {
495     if (sgnd < 0.0) {
496       *sty = *stx;
497       *fy = *fx;
498       *dy = *dx;
499     }
500     *stx = *stp;
501     *fx = *fp;
502     *dx = *dp;
503   }
504 
505   /* Compute the new step and safeguard it. */
506   stpf = PetscMin(ls->stepmax,stpf);
507   stpf = PetscMax(ls->stepmin,stpf);
508   *stp = stpf;
509   if (mtP->bracket && bound) {
510     if (*sty > *stx) *stp = PetscMin(*stx+0.66*(*sty-*stx),*stp);
511     else *stp = PetscMax(*stx+0.66*(*sty-*stx),*stp);
512   }
513   PetscFunctionReturn(0);
514 }
515