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