Actual source code: linesearchbt.c

  1: #include <petsc/private/linesearchimpl.h>
  2: #include <petsc/private/snesimpl.h>

  4: typedef struct {
  5:   PetscReal alpha; /* sufficient decrease parameter */
  6: } SNESLineSearch_BT;

  8: /*@
  9:   SNESLineSearchBTSetAlpha - Sets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` `SNESLineSearch` variant.

 11:   Input Parameters:
 12: + linesearch - linesearch context
 13: - alpha      - The descent parameter

 15:   Level: intermediate

 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;

 24:   PetscFunctionBegin;
 26:   PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, SNESLINESEARCHBT, &isbt));
 27:   if (isbt) bt->alpha = alpha;
 28:   PetscFunctionReturn(PETSC_SUCCESS);
 29: }

 31: /*@
 32:   SNESLineSearchBTGetAlpha - Gets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` variant that was set with `SNESLineSearchBTSetAlpha()`

 34:   Input Parameter:
 35: . linesearch - linesearch context

 37:   Output Parameter:
 38: . alpha - The descent parameter

 40:   Level: intermediate

 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;

 49:   PetscFunctionBegin;
 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: }

 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"};

 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;

 85:   PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
 86:   PetscCheck(jac || objective, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");

 88:   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
 89:   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));

 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));

 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:   }

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:     SNESLineSearchCheckObjectiveDomainError(snes, f);
113:   } else {
114:     f = 0.5 * PetscSqr(fnorm);
115:   }

117:   /* compute the initial slope */
118:   if (objective) {
119:     /* slope comes from the function (assumed to be the gradient of the objective) */
120:     PetscCall(VecDotRealPart(Y, F, &initslope));
121:   } else {
122:     /* slope comes from the normal equations */
123:     PetscCall(MatMult(jac, Y, W));
124:     PetscCall(VecDotRealPart(F, W, &initslope));
125:     if (initslope > 0.0) initslope = -initslope;
126:     if (initslope == 0.0) initslope = -1.0;
127:   }

129:   /* repeatedly cut lambda until the norm or objective function is not infinity or NaN or lambda is too small */
130:   while (PETSC_TRUE) {
131:     PetscCall(VecWAXPY(W, -lambda, Y, X));
132:     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
133:     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
134:       PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while checking full step length!\n"));
135:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
136:       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
137:       PetscFunctionReturn(PETSC_SUCCESS);
138:     }

140:     if (objective) PetscCall(SNESComputeObjective(snes, W, &g));
141:     else {
142:       PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
143:       if (linesearch->ops->vinorm) {
144:         gnorm = fnorm;
145:         PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
146:       } else {
147:         PetscCall(VecNorm(G, NORM_2, &gnorm));
148:       }
149:       g = 0.5 * PetscSqr(gnorm);
150:     }
151:     PetscCall(SNESLineSearchMonitor(linesearch));

153:     if (!PetscIsInfOrNanReal(g)) break;
154:     if (monitor) {
155:       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
156:       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: objective function at lambdas = %g is infinity or NaN, cutting lambda\n", (double)lambda));
157:       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
158:     }
159:     if (lambda <= minlambda) SNESCheckFunctionDomainError(snes, g);
160:     lambda *= .5;
161:   }

163:   if (!objective) PetscCall(PetscInfo(snes, "Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
164:   if (g <= f + lambda * alpha * initslope) { /* Sufficient reduction or step tolerance convergence */
165:     if (monitor) {
166:       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
167:       if (!objective) {
168:         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
169:       } else {
170:         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: old obj %14.12e new obj %14.12e\n", (double)f, (double)g));
171:       }
172:       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
173:     }
174:   } else {
175:     if (stol * xnorm > ynorm) {
176:       /* Since the full step didn't give sufficient decrease and the step is tiny, exit */
177:       PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
178:       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
179:       if (monitor) {
180:         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
181:         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n", (double)ynorm, (double)(stol * xnorm)));
182:         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
183:       }
184:       PetscFunctionReturn(PETSC_SUCCESS);
185:     }
186:     /* Here to avoid -Wmaybe-uninitiliazed warnings */
187:     lambdaprev = lambda;
188:     gprev      = g;
189:     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR) {
190:       /* Fit points with quadratic */
191:       lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
192:       lambda     = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);

194:       PetscCall(VecWAXPY(W, -lambda, Y, X));
195:       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
196:       if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
197:         PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while attempting quadratic backtracking! %" PetscInt_FMT " \n", snes->nfuncs));
198:         snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
199:         PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
200:         PetscFunctionReturn(PETSC_SUCCESS);
201:       }
202:       if (objective) {
203:         PetscCall(SNESComputeObjective(snes, W, &g));
204:         SNESLineSearchCheckObjectiveDomainError(snes, g);
205:       } else {
206:         PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
207:         if (linesearch->ops->vinorm) {
208:           gnorm = fnorm;
209:           PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
210:         } else {
211:           PetscCall(VecNorm(G, NORM_2, &gnorm));
212:         }
213:         SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm);
214:         g = 0.5 * PetscSqr(gnorm);
215:       }
216:       if (PetscIsInfOrNanReal(g)) {
217:         PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
218:         PetscCall(PetscInfo(snes, "Aborted due to infinity or NaN in function evaluation\n"));
219:         PetscFunctionReturn(PETSC_SUCCESS);
220:       }
221:       if (monitor) {
222:         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
223:         if (!objective) {
224:           PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: gnorm after quadratic fit %14.12e\n", (double)gnorm));
225:         } else {
226:           PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: obj after quadratic fit %14.12e\n", (double)g));
227:         }
228:         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
229:       }
230:     }
231:     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR && g <= f + lambda * alpha * initslope) { /* sufficient reduction */
232:       if (monitor) {
233:         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
234:         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, lambda=%18.16e\n", (double)lambda));
235:         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
236:       }
237:     } else {
238:       for (count = 0; count < max_it; count++) {
239:         if (lambda <= minlambda) {
240:           if (monitor) {
241:             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
242:             PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: unable to find good step length! After %" PetscInt_FMT " tries \n", count));
243:             if (!objective) {
244:               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));
245:             } else {
246:               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));
247:             }
248:             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
249:           }
250:           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
251:           PetscFunctionReturn(PETSC_SUCCESS);
252:         }
253:         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
254:           /* Fit points with cubic */
255:           t1 = g - f - lambda * initslope;
256:           t2 = gprev - f - lambdaprev * initslope;
257:           a  = (t1 / (lambda * lambda) - t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
258:           b  = (-lambdaprev * t1 / (lambda * lambda) + lambda * t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
259:           d  = b * b - 3 * a * initslope;
260:           if (d < 0.0) d = 0.0;
261:           if (a == 0.0) lambdatemp = -initslope / (2.0 * b);
262:           else lambdatemp = (-b + PetscSqrtReal(d)) / (3.0 * a);
263:         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
264:           lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
265:         } else if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { /* Just backtrack */
266:           lambdatemp = .5 * lambda;
267:         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "Line search order %" PetscInt_FMT " for type bt", linesearch->order);
268:         lambdaprev = lambda;
269:         gprev      = g;

271:         lambda = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
272:         PetscCall(VecWAXPY(W, -lambda, Y, X));
273:         if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
274:         if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
275:           PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while looking for good step length! %" PetscInt_FMT " \n", count));
276:           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));
277:           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
278:           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
279:           PetscFunctionReturn(PETSC_SUCCESS);
280:         }
281:         if (objective) {
282:           PetscCall(SNESComputeObjective(snes, W, &g));
283:           SNESLineSearchCheckObjectiveDomainError(snes, 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:           SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm);
293:           g = 0.5 * PetscSqr(gnorm);
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:   }

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:     SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm);
337:     PetscCall(VecNorm(Y, NORM_2, &ynorm));
338:   }

340:   /* copy the solution over */
341:   PetscCall(VecCopy(W, X));
342:   PetscCall(VecCopy(G, F));
343:   PetscCall(VecNorm(X, NORM_2, &xnorm));
344:   PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm));
345:   PetscFunctionReturn(PETSC_SUCCESS);
346: }

348: static PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
349: {
350:   PetscBool          isascii;
351:   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;

353:   PetscFunctionBegin;
354:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
355:   if (isascii) {
356:     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
357:       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n"));
358:     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
359:       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n"));
360:     }
361:     PetscCall(PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha));
362:   }
363:   PetscFunctionReturn(PETSC_SUCCESS);
364: }

366: static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
367: {
368:   PetscFunctionBegin;
369:   PetscCall(PetscFree(linesearch->data));
370:   PetscFunctionReturn(PETSC_SUCCESS);
371: }

373: static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch, PetscOptionItems PetscOptionsObject)
374: {
375:   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;

377:   PetscFunctionBegin;
378:   PetscOptionsHeadBegin(PetscOptionsObject, "SNESLineSearch BT options");
379:   PetscCall(PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL));
380:   PetscOptionsHeadEnd();
381:   PetscFunctionReturn(PETSC_SUCCESS);
382: }

384: /*MC
385:    SNESLINESEARCHBT - Backtracking line search {cite}`dennis:83`.

387:    This line search finds the minimum of a polynomial fitting either $1/2 ||F(x_k + \lambda Y_k)||_2^2$,
388:    or the objective function $G(x_k + \lambda Y_k)$ if it is provided with `SNESSetObjective()`.
389:    If this fit does not satisfy the sufficient decrease conditions, the interval shrinks
390:    and the fit is reattempted at most `max_it` times or until $\lambda$ is below `minlambda`.

392:    Options Database Keys:
393: +  -snes_linesearch_alpha 1e\-4      - slope descent parameter
394: .  -snes_linesearch_damping 1.0      - initial `lambda` on entry to the line search
395: .  -snes_linesearch_max_it 40        - maximum number of shrinking iterations in the line search
396: .  -snes_linesearch_minlambda 1e\-12 - minimum `lambda` (scaling of solution update) allowed
397: -  -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

399:    Level: advanced

401:    Note:
402:    This line search will always produce a step that is less than or equal to, in length, the full step size.

404: .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
405: M*/
406: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
407: {
408:   SNESLineSearch_BT *bt;

410:   PetscFunctionBegin;
411:   linesearch->ops->apply          = SNESLineSearchApply_BT;
412:   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
413:   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
414:   linesearch->ops->reset          = NULL;
415:   linesearch->ops->view           = SNESLineSearchView_BT;
416:   linesearch->ops->setup          = NULL;

418:   PetscCall(PetscNew(&bt));

420:   linesearch->data   = (void *)bt;
421:   linesearch->max_it = 40;
422:   linesearch->order  = SNES_LINESEARCH_ORDER_CUBIC;
423:   bt->alpha          = 1e-4;
424:   PetscFunctionReturn(PETSC_SUCCESS);
425: }