xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision a38eeca9ffd3023706cf931bf4b0fca9541fca6b)
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, initslope, alpha, stol;
60   PetscReal          t1, t2, a, b, d;
61   PetscReal          f;
62   PetscReal          g, gprev;
63   PetscViewer        monitor;
64   PetscInt           max_it, 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, NULL, NULL, NULL, NULL, &max_it));
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 
104   /* if the SNES has an objective set, use that instead of the function value */
105   if (objective) {
106     PetscCall(SNESComputeObjective(snes, X, &f));
107   } else {
108     f = 0.5 * PetscSqr(fnorm);
109   }
110 
111   /* compute the initial slope */
112   if (objective) {
113     /* slope comes from the function (assumed to be the gradient of the objective) */
114     PetscCall(VecDotRealPart(Y, F, &initslope));
115   } else {
116     /* slope comes from the normal equations */
117     PetscCall(MatMult(jac, Y, W));
118     PetscCall(VecDotRealPart(F, W, &initslope));
119     if (initslope > 0.0) initslope = -initslope;
120     if (initslope == 0.0) initslope = -1.0;
121   }
122 
123   while (PETSC_TRUE) {
124     PetscCall(VecWAXPY(W, -lambda, Y, X));
125     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
126     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
127       PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while checking full step length!\n"));
128       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
129       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
130       PetscFunctionReturn(PETSC_SUCCESS);
131     }
132 
133     if (objective) {
134       PetscCall(SNESComputeObjective(snes, W, &g));
135     } else {
136       PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
137       if (linesearch->ops->vinorm) {
138         gnorm = fnorm;
139         PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
140       } else {
141         PetscCall(VecNorm(G, NORM_2, &gnorm));
142       }
143       g = 0.5 * PetscSqr(gnorm);
144     }
145     PetscCall(SNESLineSearchMonitor(linesearch));
146 
147     if (!PetscIsInfOrNanReal(g)) break;
148     if (monitor) {
149       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
150       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: objective function at lambdas = %g is Inf or Nan, cutting lambda\n", (double)lambda));
151       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
152     }
153     if (lambda <= minlambda) SNESCheckFunctionNorm(snes, g);
154     lambda *= .5;
155   }
156 
157   if (!objective) PetscCall(PetscInfo(snes, "Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
158   if (g <= f + lambda * alpha * initslope) { /* Sufficient reduction or step tolerance convergence */
159     if (monitor) {
160       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
161       if (!objective) {
162         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
163       } else {
164         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: old obj %14.12e new obj %14.12e\n", (double)f, (double)g));
165       }
166       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
167     }
168   } else {
169     if (stol * xnorm > ynorm) {
170       /* Since the full step didn't give sufficient decrease and the step is tiny, exit */
171       PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
172       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
173       if (monitor) {
174         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
175         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n", (double)ynorm, (double)(stol * xnorm)));
176         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
177       }
178       PetscFunctionReturn(PETSC_SUCCESS);
179     }
180     /* Here to avoid -Wmaybe-uninitiliazed warnings */
181     lambdaprev = lambda;
182     gprev      = g;
183     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR) {
184       /* Fit points with quadratic */
185       lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
186       lambda     = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
187 
188       PetscCall(VecWAXPY(W, -lambda, Y, X));
189       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
190       if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
191         PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while attempting quadratic backtracking! %" PetscInt_FMT " \n", snes->nfuncs));
192         snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
193         PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
194         PetscFunctionReturn(PETSC_SUCCESS);
195       }
196       if (objective) {
197         PetscCall(SNESComputeObjective(snes, W, &g));
198       } else {
199         PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
200         if (linesearch->ops->vinorm) {
201           gnorm = fnorm;
202           PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
203         } else {
204           PetscCall(VecNorm(G, NORM_2, &gnorm));
205         }
206         g = 0.5 * PetscSqr(gnorm);
207       }
208       if (PetscIsInfOrNanReal(g)) {
209         PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
210         PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
211         PetscFunctionReturn(PETSC_SUCCESS);
212       }
213       if (monitor) {
214         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
215         if (!objective) {
216           PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: gnorm after quadratic fit %14.12e\n", (double)gnorm));
217         } else {
218           PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: obj after quadratic fit %14.12e\n", (double)g));
219         }
220         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
221       }
222     }
223     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR && g <= f + lambda * alpha * initslope) { /* sufficient reduction */
224       if (monitor) {
225         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
226         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, lambda=%18.16e\n", (double)lambda));
227         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
228       }
229     } else {
230       for (count = 0; count < max_it; count++) {
231         if (lambda <= minlambda) {
232           if (monitor) {
233             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
234             PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: unable to find good step length! After %" PetscInt_FMT " tries \n", count));
235             if (!objective) {
236               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));
237             } else {
238               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));
239             }
240             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
241           }
242           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
243           PetscFunctionReturn(PETSC_SUCCESS);
244         }
245         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
246           /* Fit points with cubic */
247           t1 = g - f - lambda * initslope;
248           t2 = gprev - f - lambdaprev * initslope;
249           a  = (t1 / (lambda * lambda) - t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
250           b  = (-lambdaprev * t1 / (lambda * lambda) + lambda * t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
251           d  = b * b - 3 * a * initslope;
252           if (d < 0.0) d = 0.0;
253           if (a == 0.0) lambdatemp = -initslope / (2.0 * b);
254           else lambdatemp = (-b + PetscSqrtReal(d)) / (3.0 * a);
255         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
256           lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
257         } else if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { /* Just backtrack */
258           lambdatemp = .5 * lambda;
259         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "Line search order %" PetscInt_FMT " for type bt", linesearch->order);
260         lambdaprev = lambda;
261         gprev      = g;
262 
263         lambda = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
264         PetscCall(VecWAXPY(W, -lambda, Y, X));
265         if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
266         if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
267           PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while looking for good step length! %" PetscInt_FMT " \n", count));
268           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));
269           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
270           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
271           PetscFunctionReturn(PETSC_SUCCESS);
272         }
273         if (objective) {
274           PetscCall(SNESComputeObjective(snes, W, &g));
275         } else {
276           PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
277           if (linesearch->ops->vinorm) {
278             gnorm = fnorm;
279             PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
280           } else {
281             PetscCall(VecNorm(G, NORM_2, &gnorm));
282           }
283           g = 0.5 * PetscSqr(gnorm);
284         }
285         if (PetscIsInfOrNanReal(g)) {
286           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
287           PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
288           PetscFunctionReturn(PETSC_SUCCESS);
289         }
290         if (g <= f + lambda * alpha * initslope) { /* is reduction enough? */
291           if (monitor) {
292             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
293             if (!objective) {
294               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: %s step, current gnorm %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)gnorm, (double)lambda));
295               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
296             } else {
297               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: %s step, obj %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)g, (double)lambda));
298               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
299             }
300           }
301           break;
302         } else if (monitor) {
303           PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
304           if (!objective) {
305             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));
306             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
307           } else {
308             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));
309             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
310           }
311         }
312       }
313     }
314   }
315 
316   /* postcheck */
317   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
318   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
319   if (changed_y) {
320     if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
321     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
322   }
323   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
324     PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
325     if (linesearch->ops->vinorm) {
326       gnorm = fnorm;
327       PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
328     } else {
329       PetscCall(VecNorm(G, NORM_2, &gnorm));
330     }
331     PetscCall(VecNorm(Y, NORM_2, &ynorm));
332     if (PetscIsInfOrNanReal(gnorm)) {
333       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
334       PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
335       PetscFunctionReturn(PETSC_SUCCESS);
336     }
337   }
338 
339   /* copy the solution over */
340   PetscCall(VecCopy(W, X));
341   PetscCall(VecCopy(G, F));
342   PetscCall(VecNorm(X, NORM_2, &xnorm));
343   PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm));
344   PetscFunctionReturn(PETSC_SUCCESS);
345 }
346 
347 static PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
348 {
349   PetscBool          isascii;
350   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
351 
352   PetscFunctionBegin;
353   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
354   if (isascii) {
355     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
356       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n"));
357     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
358       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n"));
359     }
360     PetscCall(PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha));
361   }
362   PetscFunctionReturn(PETSC_SUCCESS);
363 }
364 
365 static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
366 {
367   PetscFunctionBegin;
368   PetscCall(PetscFree(linesearch->data));
369   PetscFunctionReturn(PETSC_SUCCESS);
370 }
371 
372 static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch, PetscOptionItems PetscOptionsObject)
373 {
374   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
375 
376   PetscFunctionBegin;
377   PetscOptionsHeadBegin(PetscOptionsObject, "SNESLineSearch BT options");
378   PetscCall(PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL));
379   PetscOptionsHeadEnd();
380   PetscFunctionReturn(PETSC_SUCCESS);
381 }
382 
383 /*MC
384    SNESLINESEARCHBT - Backtracking line search {cite}`dennis:83`.
385 
386    This line search finds the minimum of a polynomial fitting either $1/2 ||F(x_k + \lambda Y_k)||_2^2$,
387    or the objective function $G(x_k + \lambda Y_k)$ if it is provided with `SNESSetObjective()`.
388    If this fit does not satisfy the sufficient decrease conditions, the interval shrinks
389    and the fit is reattempted at most `max_it` times or until $\lambda$ is below `minlambda`.
390 
391    Options Database Keys:
392 +  -snes_linesearch_alpha <1e\-4>      - slope descent parameter
393 .  -snes_linesearch_damping <1.0>      - initial `lambda` on entry to the line search
394 .  -snes_linesearch_max_it <40>        - maximum number of shrinking iterations in the line search
395 .  -snes_linesearch_minlambda <1e\-12> - minimum `lambda` (scaling of solution update) allowed
396 -  -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
397 
398    Level: advanced
399 
400    Note:
401    This line search will always produce a step that is less than or equal to, in length, the full step size.
402 
403 .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
404 M*/
405 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
406 {
407   SNESLineSearch_BT *bt;
408 
409   PetscFunctionBegin;
410   linesearch->ops->apply          = SNESLineSearchApply_BT;
411   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
412   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
413   linesearch->ops->reset          = NULL;
414   linesearch->ops->view           = SNESLineSearchView_BT;
415   linesearch->ops->setup          = NULL;
416 
417   PetscCall(PetscNew(&bt));
418 
419   linesearch->data   = (void *)bt;
420   linesearch->max_it = 40;
421   linesearch->order  = SNES_LINESEARCH_ORDER_CUBIC;
422   bt->alpha          = 1e-4;
423   PetscFunctionReturn(PETSC_SUCCESS);
424 }
425