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