xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision 6a4a1270853b5b1995c94de98f44d28bec57470a)
1 #include <petsc/private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
2 #include <petsc/private/snesimpl.h>
3 
4 typedef struct {
5   PetscReal alpha; /* sufficient decrease parameter */
6 } SNESLineSearch_BT;
7 
8 /*@
9    SNESLineSearchBTSetAlpha - Sets the descent parameter, alpha, in the `SNESLINESEARCHBT` variant.
10 
11    Input Parameters:
12 +  linesearch - linesearch context
13 -  alpha - The descent parameter
14 
15    Level: intermediate
16 
17 .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetTolerances()`, `SNESLINESEARCHBT`, `SNESLineSearchBTGetAlpha()`
18 @*/
19 PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha) {
20   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
21 
22   PetscFunctionBegin;
23   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
24   bt->alpha = alpha;
25   PetscFunctionReturn(0);
26 }
27 
28 /*@
29    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the `SNESLINESEARCHBT` variant.
30 
31    Input Parameter:
32 .  linesearch - linesearch context
33 
34    Output Parameter:
35 .  alpha - The descent parameter
36 
37    Level: intermediate
38 
39 .seealso: `SNESLineSearchGetLambda()`, `SNESLineSearchGetTolerances()` `SNESLINESEARCHBT`, `SNESLineSearchBTGetAlpha()`
40 @*/
41 PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha) {
42   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
43 
44   PetscFunctionBegin;
45   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
46   *alpha = bt->alpha;
47   PetscFunctionReturn(0);
48 }
49 
50 static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) {
51   PetscBool          changed_y, changed_w;
52   Vec                X, F, Y, W, G;
53   SNES               snes;
54   PetscReal          fnorm, xnorm, ynorm, gnorm;
55   PetscReal          lambda, lambdatemp, lambdaprev, minlambda, maxstep, initslope, alpha, stol;
56   PetscReal          t1, t2, a, b, d;
57   PetscReal          f;
58   PetscReal          g, gprev;
59   PetscViewer        monitor;
60   PetscInt           max_its, count;
61   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
62   Mat                jac;
63   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
64 
65   PetscFunctionBegin;
66   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G));
67   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
68   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
69   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
70   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
71   PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, &maxstep, NULL, NULL, NULL, &max_its));
72   PetscCall(SNESGetTolerances(snes, NULL, NULL, &stol, NULL, NULL));
73   PetscCall(SNESGetObjective(snes, &objective, NULL));
74   alpha = bt->alpha;
75 
76   PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
77   PetscCheck(jac || objective, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
78 
79   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
80   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
81 
82   PetscCall(VecNormBegin(Y, NORM_2, &ynorm));
83   PetscCall(VecNormBegin(X, NORM_2, &xnorm));
84   PetscCall(VecNormEnd(Y, NORM_2, &ynorm));
85   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
86 
87   if (ynorm == 0.0) {
88     if (monitor) {
89       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
90       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Initial direction and size is 0\n"));
91       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
92     }
93     PetscCall(VecCopy(X, W));
94     PetscCall(VecCopy(F, G));
95     PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
96     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
97     PetscFunctionReturn(0);
98   }
99   if (ynorm > maxstep) { /* Step too big, so scale back */
100     if (monitor) {
101       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
102       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep / ynorm), (double)ynorm));
103       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
104     }
105     PetscCall(VecScale(Y, maxstep / (ynorm)));
106     ynorm = maxstep;
107   }
108 
109   /* if the SNES has an objective set, use that instead of the function value */
110   if (objective) {
111     PetscCall(SNESComputeObjective(snes, X, &f));
112   } else {
113     f = fnorm * fnorm;
114   }
115 
116   /* compute the initial slope */
117   if (objective) {
118     /* slope comes from the function (assumed to be the gradient of the objective */
119     PetscCall(VecDotRealPart(Y, F, &initslope));
120   } else {
121     /* slope comes from the normal equations */
122     PetscCall(MatMult(jac, Y, W));
123     PetscCall(VecDotRealPart(F, W, &initslope));
124     if (initslope > 0.0) initslope = -initslope;
125     if (initslope == 0.0) initslope = -1.0;
126   }
127 
128   while (PETSC_TRUE) {
129     PetscCall(VecWAXPY(W, -lambda, Y, X));
130     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
131     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
132       PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while checking full step length!\n"));
133       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
134       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
135       PetscFunctionReturn(0);
136     }
137 
138     if (objective) {
139       PetscCall(SNESComputeObjective(snes, W, &g));
140     } else {
141       PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
142       if (linesearch->ops->vinorm) {
143         gnorm = fnorm;
144         PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
145       } else {
146         PetscCall(VecNorm(G, NORM_2, &gnorm));
147       }
148       g = PetscSqr(gnorm);
149     }
150     PetscCall(SNESLineSearchMonitor(linesearch));
151 
152     if (!PetscIsInfOrNanReal(g)) break;
153     if (monitor) {
154       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
155       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: objective function at lambdas = %g is Inf or Nan, cutting lambda\n", (double)lambda));
156       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
157     }
158     if (lambda <= minlambda) SNESCheckFunctionNorm(snes, g);
159     lambda = .5 * lambda;
160   }
161 
162   if (!objective) PetscCall(PetscInfo(snes, "Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
163   if (.5 * g <= .5 * f + lambda * alpha * initslope) { /* Sufficient reduction or step tolerance convergence */
164     if (monitor) {
165       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
166       if (!objective) {
167         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
168       } else {
169         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: obj %14.12e obj %14.12e\n", (double)f, (double)g));
170       }
171       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
172     }
173   } else {
174     /* Since the full step didn't work and the step is tiny, quit */
175     if (stol * xnorm > ynorm) {
176       PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
177       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
178       if (monitor) {
179         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
180         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n", (double)ynorm, (double)(stol * xnorm)));
181         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
182       }
183       PetscFunctionReturn(0);
184     }
185     /* Fit points with quadratic */
186     lambdatemp = -initslope / (g - f - 2.0 * lambda * initslope);
187     lambdaprev = lambda;
188     gprev      = g;
189     if (lambdatemp > .5 * lambda) lambdatemp = .5 * lambda;
190     if (lambdatemp <= .1 * lambda) lambda = .1 * lambda;
191     else lambda = lambdatemp;
192 
193     PetscCall(VecWAXPY(W, -lambda, Y, X));
194     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
195     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
196       PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while attempting quadratic backtracking! %" PetscInt_FMT " \n", snes->nfuncs));
197       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
198       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
199       PetscFunctionReturn(0);
200     }
201     if (objective) {
202       PetscCall(SNESComputeObjective(snes, W, &g));
203     } else {
204       PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
205       if (linesearch->ops->vinorm) {
206         gnorm = fnorm;
207         PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
208       } else {
209         PetscCall(VecNorm(G, NORM_2, &gnorm));
210       }
211       g = PetscSqr(gnorm);
212     }
213     if (PetscIsInfOrNanReal(g)) {
214       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
215       PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
216       PetscFunctionReturn(0);
217     }
218     if (monitor) {
219       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
220       if (!objective) {
221         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: gnorm after quadratic fit %14.12e\n", (double)gnorm));
222       } else {
223         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: obj after quadratic fit %14.12e\n", (double)g));
224       }
225       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
226     }
227     if (.5 * g < .5 * f + lambda * alpha * initslope) { /* sufficient reduction */
228       if (monitor) {
229         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
230         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, lambda=%18.16e\n", (double)lambda));
231         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
232       }
233     } else {
234       /* Fit points with cubic */
235       for (count = 0; count < max_its; count++) {
236         if (lambda <= minlambda) {
237           if (monitor) {
238             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
239             PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: unable to find good step length! After %" PetscInt_FMT " tries \n", count));
240             if (!objective) {
241               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope));
242             } else {
243               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope));
244             }
245             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
246           }
247           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
248           PetscFunctionReturn(0);
249         }
250         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
251           t1 = .5 * (g - f) - lambda * initslope;
252           t2 = .5 * (gprev - f) - lambdaprev * initslope;
253           a  = (t1 / (lambda * lambda) - t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
254           b  = (-lambdaprev * t1 / (lambda * lambda) + lambda * t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
255           d  = b * b - 3 * a * initslope;
256           if (d < 0.0) d = 0.0;
257           if (a == 0.0) lambdatemp = -initslope / (2.0 * b);
258           else lambdatemp = (-b + PetscSqrtReal(d)) / (3.0 * a);
259 
260         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
261           lambdatemp = -initslope / (g - f - 2.0 * initslope);
262         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt");
263         lambdaprev = lambda;
264         gprev      = g;
265         if (lambdatemp > .5 * lambda) lambdatemp = .5 * lambda;
266         if (lambdatemp <= .1 * lambda) lambda = .1 * lambda;
267         else lambda = lambdatemp;
268         PetscCall(VecWAXPY(W, -lambda, Y, X));
269         if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
270         if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
271           PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while looking for good step length! %" PetscInt_FMT " \n", count));
272           if (!objective) PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)fnorm, (double)gnorm, (double)ynorm, (double)lambda, (double)initslope));
273           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
274           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
275           PetscFunctionReturn(0);
276         }
277         if (objective) {
278           PetscCall(SNESComputeObjective(snes, W, &g));
279         } else {
280           PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
281           if (linesearch->ops->vinorm) {
282             gnorm = fnorm;
283             PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
284           } else {
285             PetscCall(VecNorm(G, NORM_2, &gnorm));
286           }
287           g = PetscSqr(gnorm);
288         }
289         if (PetscIsInfOrNanReal(g)) {
290           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
291           PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
292           PetscFunctionReturn(0);
293         }
294         if (.5 * g < .5 * f + lambda * alpha * initslope) { /* is reduction enough? */
295           if (monitor) {
296             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
297             if (!objective) {
298               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
299                 PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n", (double)gnorm, (double)lambda));
300               } else {
301                 PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n", (double)gnorm, (double)lambda));
302               }
303               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
304             } else {
305               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
306                 PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n", (double)g, (double)lambda));
307               } else {
308                 PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n", (double)g, (double)lambda));
309               }
310               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
311             }
312           }
313           break;
314         } else if (monitor) {
315           PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
316           if (!objective) {
317             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
318               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n", (double)gnorm, (double)lambda));
319             } else {
320               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n", (double)gnorm, (double)lambda));
321             }
322             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
323           } else {
324             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
325               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Cubic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n", (double)g, (double)lambda));
326             } else {
327               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n", (double)g, (double)lambda));
328             }
329             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
330           }
331         }
332       }
333     }
334   }
335 
336   /* postcheck */
337   /* update Y to lambda*Y so that W is consistent with  X - lambda*Y */
338   PetscCall(VecScale(Y, lambda));
339   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
340   if (changed_y) {
341     PetscCall(VecWAXPY(W, -1.0, Y, X));
342     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
343   }
344   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
345     PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
346     if (linesearch->ops->vinorm) {
347       gnorm = fnorm;
348       PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
349     } else {
350       PetscCall(VecNorm(G, NORM_2, &gnorm));
351     }
352     PetscCall(VecNorm(Y, NORM_2, &ynorm));
353     if (PetscIsInfOrNanReal(gnorm)) {
354       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
355       PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
356       PetscFunctionReturn(0);
357     }
358   }
359 
360   /* copy the solution over */
361   PetscCall(VecCopy(W, X));
362   PetscCall(VecCopy(G, F));
363   PetscCall(VecNorm(X, NORM_2, &xnorm));
364   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
365   PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm));
366   PetscFunctionReturn(0);
367 }
368 
369 PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) {
370   PetscBool          iascii;
371   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
372 
373   PetscFunctionBegin;
374   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
375   if (iascii) {
376     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
377       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n"));
378     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
379       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n"));
380     }
381     PetscCall(PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha));
382   }
383   PetscFunctionReturn(0);
384 }
385 
386 static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) {
387   PetscFunctionBegin;
388   PetscCall(PetscFree(linesearch->data));
389   PetscFunctionReturn(0);
390 }
391 
392 static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch, PetscOptionItems *PetscOptionsObject) {
393   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
394 
395   PetscFunctionBegin;
396   PetscOptionsHeadBegin(PetscOptionsObject, "SNESLineSearch BT options");
397   PetscCall(PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL));
398   PetscOptionsHeadEnd();
399   PetscFunctionReturn(0);
400 }
401 
402 /*MC
403    SNESLINESEARCHBT - Backtracking line search.
404 
405    This line search finds the minimum of a polynomial fitting of the L2 norm of the
406    function or the objective function if it is provided with `SNESSetObjective()`. If this fit does not satisfy the conditions for progress, the interval shrinks
407    and the fit is reattempted at most max_it times or until lambda is below minlambda.
408 
409    Options Database Keys:
410 +  -snes_linesearch_alpha <1e\-4> - slope descent parameter
411 .  -snes_linesearch_damping <1.0> - initial step length
412 .  -snes_linesearch_maxstep <length> - if the length the initial step is larger than this then the
413                                        step is scaled back to be of this length at the beginning of the line search
414 .  -snes_linesearch_max_it <40> - maximum number of shrinking step
415 .  -snes_linesearch_minlambda <1e\-12> - minimum step length allowed
416 -  -snes_linesearch_order <cubic,quadratic> - order of the approximation
417 
418    Level: advanced
419 
420    Note:
421    This line search will always produce a step that is less than or equal to, in length, the full step size.
422 
423    Reference:
424 .  - * - "Numerical Methods for Unconstrained Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
425 
426 .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
427 M*/
428 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) {
429   SNESLineSearch_BT *bt;
430 
431   PetscFunctionBegin;
432   linesearch->ops->apply          = SNESLineSearchApply_BT;
433   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
434   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
435   linesearch->ops->reset          = NULL;
436   linesearch->ops->view           = SNESLineSearchView_BT;
437   linesearch->ops->setup          = NULL;
438 
439   PetscCall(PetscNewLog(linesearch, &bt));
440 
441   linesearch->data    = (void *)bt;
442   linesearch->max_its = 40;
443   linesearch->order   = SNES_LINESEARCH_ORDER_CUBIC;
444   bt->alpha           = 1e-4;
445   PetscFunctionReturn(0);
446 }
447