xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision b2ccae6bdc8edea944f1c160ca3b2eb32c69ecb2)
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   PetscBool          isbt;
23 
24   PetscFunctionBegin;
25   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
26   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, SNESLINESEARCHBT, &isbt));
27   if (isbt) bt->alpha = alpha;
28   PetscFunctionReturn(PETSC_SUCCESS);
29 }
30 
31 /*@
32   SNESLineSearchBTGetAlpha - Gets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` variant that was set with `SNESLineSearchBTSetAlpha()`
33 
34   Input Parameter:
35 . linesearch - linesearch context
36 
37   Output Parameter:
38 . alpha - The descent parameter
39 
40   Level: intermediate
41 
42 .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchGetLambda()`, `SNESLineSearchGetTolerances()`, `SNESLINESEARCHBT`, `SNESLineSearchBTSetAlpha()`
43 @*/
44 PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
45 {
46   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
47   PetscBool          isbt;
48 
49   PetscFunctionBegin;
50   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
51   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, SNESLINESEARCHBT, &isbt));
52   PetscCheck(isbt, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Not for type %s", ((PetscObject)linesearch)->type_name);
53   *alpha = bt->alpha;
54   PetscFunctionReturn(PETSC_SUCCESS);
55 }
56 
57 static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch)
58 {
59   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
60   PetscBool          changed_y, changed_w;
61   Vec                X, F, Y, W, G;
62   SNES               snes;
63   PetscReal          fnorm, xnorm, ynorm, gnorm;
64   PetscReal          lambda, lambdatemp, lambdaprev, minlambda, initslope, alpha, stol;
65   PetscReal          t1, t2, a, b, d;
66   PetscReal          f;
67   PetscReal          g, gprev;
68   PetscViewer        monitor;
69   PetscInt           max_it, count;
70   Mat                jac;
71   SNESObjectiveFn   *objective;
72   const char *const  ordStr[] = {"Linear", "Quadratic", "Cubic"};
73 
74   PetscFunctionBegin;
75   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G));
76   PetscCall(SNESLineSearchGetNorms(linesearch, NULL, &fnorm, NULL));
77   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
78   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
79   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
80   PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, NULL, NULL, NULL, NULL, &max_it));
81   PetscCall(SNESGetTolerances(snes, NULL, NULL, &stol, NULL, NULL));
82   PetscCall(SNESGetObjective(snes, &objective, NULL));
83   alpha = bt->alpha;
84 
85   PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
86   PetscCheck(jac || objective, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
87 
88   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
89   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
90 
91   PetscCall(VecNormBegin(Y, NORM_2, &ynorm));
92   PetscCall(VecNormBegin(X, NORM_2, &xnorm));
93   PetscCall(VecNormEnd(Y, NORM_2, &ynorm));
94   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
95 
96   if (ynorm == 0.0) {
97     if (monitor) {
98       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
99       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Initial direction and size is 0\n"));
100       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
101     }
102     PetscCall(VecCopy(X, W));
103     PetscCall(VecCopy(F, G));
104     PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
105     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
106     PetscFunctionReturn(PETSC_SUCCESS);
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 = 0.5 * PetscSqr(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(PETSC_SUCCESS);
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 = 0.5 * 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;
160   }
161 
162   if (!objective) PetscCall(PetscInfo(snes, "Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
163   if (g <= 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: old obj %14.12e new obj %14.12e\n", (double)f, (double)g));
170       }
171       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
172     }
173   } else {
174     if (stol * xnorm > ynorm) {
175       /* Since the full step didn't give sufficient decrease and the step is tiny, exit */
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(PETSC_SUCCESS);
184     }
185     /* Here to avoid -Wmaybe-uninitiliazed warnings */
186     lambdaprev = lambda;
187     gprev      = g;
188     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR) {
189       /* Fit points with quadratic */
190       lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
191       lambda     = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
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(PETSC_SUCCESS);
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 = 0.5 * 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(PETSC_SUCCESS);
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     }
228     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR && g <= f + lambda * alpha * initslope) { /* sufficient reduction */
229       if (monitor) {
230         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
231         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, lambda=%18.16e\n", (double)lambda));
232         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
233       }
234     } else {
235       for (count = 0; count < max_it; 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(PETSC_SUCCESS);
249         }
250         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
251           /* Fit points with cubic */
252           t1 = g - f - lambda * initslope;
253           t2 = gprev - f - lambdaprev * initslope;
254           a  = (t1 / (lambda * lambda) - t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
255           b  = (-lambdaprev * t1 / (lambda * lambda) + lambda * t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
256           d  = b * b - 3 * a * initslope;
257           if (d < 0.0) d = 0.0;
258           if (a == 0.0) lambdatemp = -initslope / (2.0 * b);
259           else lambdatemp = (-b + PetscSqrtReal(d)) / (3.0 * a);
260         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
261           lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
262         } else if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { /* Just backtrack */
263           lambdatemp = .5 * lambda;
264         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "Line search order %" PetscInt_FMT " for type bt", linesearch->order);
265         lambdaprev = lambda;
266         gprev      = g;
267 
268         lambda = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
269         PetscCall(VecWAXPY(W, -lambda, Y, X));
270         if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
271         if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
272           PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while looking for good step length! %" PetscInt_FMT " \n", count));
273           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));
274           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
275           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
276           PetscFunctionReturn(PETSC_SUCCESS);
277         }
278         if (objective) {
279           PetscCall(SNESComputeObjective(snes, W, &g));
280         } else {
281           PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
282           if (linesearch->ops->vinorm) {
283             gnorm = fnorm;
284             PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
285           } else {
286             PetscCall(VecNorm(G, NORM_2, &gnorm));
287           }
288           g = 0.5 * PetscSqr(gnorm);
289         }
290         if (PetscIsInfOrNanReal(g)) {
291           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
292           PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
293           PetscFunctionReturn(PETSC_SUCCESS);
294         }
295         if (g <= f + lambda * alpha * initslope) { /* is reduction enough? */
296           if (monitor) {
297             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
298             if (!objective) {
299               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: %s step, current gnorm %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)gnorm, (double)lambda));
300               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
301             } else {
302               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: %s step, obj %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)g, (double)lambda));
303               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
304             }
305           }
306           break;
307         } else if (monitor) {
308           PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
309           if (!objective) {
310             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));
311             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
312           } else {
313             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));
314             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
315           }
316         }
317       }
318     }
319   }
320 
321   /* postcheck */
322   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
323   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
324   if (changed_y) {
325     if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
326     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
327   }
328   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
329     PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
330     if (linesearch->ops->vinorm) {
331       gnorm = fnorm;
332       PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
333     } else {
334       PetscCall(VecNorm(G, NORM_2, &gnorm));
335     }
336     PetscCall(VecNorm(Y, NORM_2, &ynorm));
337     if (PetscIsInfOrNanReal(gnorm)) {
338       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
339       PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
340       PetscFunctionReturn(PETSC_SUCCESS);
341     }
342   }
343 
344   /* copy the solution over */
345   PetscCall(VecCopy(W, X));
346   PetscCall(VecCopy(G, F));
347   PetscCall(VecNorm(X, NORM_2, &xnorm));
348   PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm));
349   PetscFunctionReturn(PETSC_SUCCESS);
350 }
351 
352 static PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
353 {
354   PetscBool          isascii;
355   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
356 
357   PetscFunctionBegin;
358   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
359   if (isascii) {
360     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
361       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n"));
362     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
363       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n"));
364     }
365     PetscCall(PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha));
366   }
367   PetscFunctionReturn(PETSC_SUCCESS);
368 }
369 
370 static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
371 {
372   PetscFunctionBegin;
373   PetscCall(PetscFree(linesearch->data));
374   PetscFunctionReturn(PETSC_SUCCESS);
375 }
376 
377 static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch, PetscOptionItems PetscOptionsObject)
378 {
379   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
380 
381   PetscFunctionBegin;
382   PetscOptionsHeadBegin(PetscOptionsObject, "SNESLineSearch BT options");
383   PetscCall(PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL));
384   PetscOptionsHeadEnd();
385   PetscFunctionReturn(PETSC_SUCCESS);
386 }
387 
388 /*MC
389    SNESLINESEARCHBT - Backtracking line search {cite}`dennis:83`.
390 
391    This line search finds the minimum of a polynomial fitting either $1/2 ||F(x_k + \lambda Y_k)||_2^2$,
392    or the objective function $G(x_k + \lambda Y_k)$ if it is provided with `SNESSetObjective()`.
393    If this fit does not satisfy the sufficient decrease conditions, the interval shrinks
394    and the fit is reattempted at most `max_it` times or until $\lambda$ is below `minlambda`.
395 
396    Options Database Keys:
397 +  -snes_linesearch_alpha <1e\-4>      - slope descent parameter
398 .  -snes_linesearch_damping <1.0>      - initial `lambda` on entry to the line search
399 .  -snes_linesearch_max_it <40>        - maximum number of shrinking iterations in the line search
400 .  -snes_linesearch_minlambda <1e\-12> - minimum `lambda` (scaling of solution update) allowed
401 -  -snes_linesearch_order <3>          - order of the polynomial fit, must be 1, 2, or 3. With order 1, it performs a simple backtracking without any curve fitting
402 
403    Level: advanced
404 
405    Note:
406    This line search will always produce a step that is less than or equal to, in length, the full step size.
407 
408 .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
409 M*/
410 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
411 {
412   SNESLineSearch_BT *bt;
413 
414   PetscFunctionBegin;
415   linesearch->ops->apply          = SNESLineSearchApply_BT;
416   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
417   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
418   linesearch->ops->reset          = NULL;
419   linesearch->ops->view           = SNESLineSearchView_BT;
420   linesearch->ops->setup          = NULL;
421 
422   PetscCall(PetscNew(&bt));
423 
424   linesearch->data   = (void *)bt;
425   linesearch->max_it = 40;
426   linesearch->order  = SNES_LINESEARCH_ORDER_CUBIC;
427   bt->alpha          = 1e-4;
428   PetscFunctionReturn(PETSC_SUCCESS);
429 }
430