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