xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision a99ef635ef3fb31080710bc8e0dbbde292199c98)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
371b4ebd2SPeter Brune 
471b4ebd2SPeter Brune typedef struct {
571b4ebd2SPeter Brune   PetscReal alpha; /* sufficient decrease parameter */
6f1c6b773SPeter Brune } SNESLineSearch_BT;
771b4ebd2SPeter Brune 
82f4102e2SPeter Brune /*@
9420bcc1bSBarry Smith   SNESLineSearchBTSetAlpha - Sets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` `SNESLineSearch` variant.
102f4102e2SPeter Brune 
112f4102e2SPeter Brune   Input Parameters:
122f4102e2SPeter Brune + linesearch - linesearch context
132f4102e2SPeter Brune - alpha      - The descent parameter
142f4102e2SPeter Brune 
152f4102e2SPeter Brune   Level: intermediate
162f4102e2SPeter Brune 
17420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetTolerances()`, `SNESLINESEARCHBT`, `SNESLineSearchBTGetAlpha()`
182f4102e2SPeter Brune @*/
19d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha)
20d71ae5a4SJacob Faibussowitsch {
21d7ca06c0SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
226e111a19SKarl Rupp 
232f4102e2SPeter Brune   PetscFunctionBegin;
242f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
252f4102e2SPeter Brune   bt->alpha = alpha;
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272f4102e2SPeter Brune }
282f4102e2SPeter Brune 
292f4102e2SPeter Brune /*@
30420bcc1bSBarry Smith   SNESLineSearchBTGetAlpha - Gets the descent parameter, `alpha`, in the `SNESLINESEARCHBT` variant that was set with `SNESLineSearchBTSetAlpha()`
312f4102e2SPeter Brune 
32f6dfbefdSBarry Smith   Input Parameter:
332f4102e2SPeter Brune . linesearch - linesearch context
348e557f58SPeter Brune 
35f6dfbefdSBarry Smith   Output Parameter:
362f4102e2SPeter Brune . alpha - The descent parameter
372f4102e2SPeter Brune 
382f4102e2SPeter Brune   Level: intermediate
392f4102e2SPeter Brune 
400241b274SPierre Jolivet .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchGetLambda()`, `SNESLineSearchGetTolerances()`, `SNESLINESEARCHBT`, `SNESLineSearchBTSetAlpha()`
412f4102e2SPeter Brune @*/
42d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
43d71ae5a4SJacob Faibussowitsch {
44d7ca06c0SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
456e111a19SKarl Rupp 
462f4102e2SPeter Brune   PetscFunctionBegin;
472f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1);
482f4102e2SPeter Brune   *alpha = bt->alpha;
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
502f4102e2SPeter Brune }
512f4102e2SPeter Brune 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch)
53d71ae5a4SJacob Faibussowitsch {
546b72add0SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
5571b4ebd2SPeter Brune   PetscBool          changed_y, changed_w;
5671b4ebd2SPeter Brune   Vec                X, F, Y, W, G;
5771b4ebd2SPeter Brune   SNES               snes;
588a430afbSPeter Brune   PetscReal          fnorm, xnorm, ynorm, gnorm;
59*a99ef635SJonas Heinzmann   PetscReal          lambda, lambdatemp, lambdaprev, minlambda, initslope, alpha, stol;
60dfcd3828SPeter Brune   PetscReal          t1, t2, a, b, d;
618a430afbSPeter Brune   PetscReal          f;
628a430afbSPeter Brune   PetscReal          g, gprev;
6371b4ebd2SPeter Brune   PetscViewer        monitor;
64*a99ef635SJonas Heinzmann   PetscInt           max_it, count;
6571b4ebd2SPeter Brune   Mat                jac;
666b72add0SBarry Smith   SNESObjectiveFn   *objective;
6732661dccSStefano Zampini   const char *const  ordStr[] = {"Linear", "Quadratic", "Cubic"};
6871b4ebd2SPeter Brune 
6971b4ebd2SPeter Brune   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G));
71a68bbae5SBarry Smith   PetscCall(SNESLineSearchGetNorms(linesearch, NULL, &fnorm, NULL));
729566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
739566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
749566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
75*a99ef635SJonas Heinzmann   PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, NULL, NULL, NULL, NULL, &max_it));
769566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, NULL, NULL, &stol, NULL, NULL));
779566063dSJacob Faibussowitsch   PetscCall(SNESGetObjective(snes, &objective, NULL));
7871b4ebd2SPeter Brune   alpha = bt->alpha;
7971b4ebd2SPeter Brune 
809566063dSJacob Faibussowitsch   PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
8108401ef6SPierre Jolivet   PetscCheck(jac || objective, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
82c69d1a72SBarry Smith 
839566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
849566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
8571b4ebd2SPeter Brune 
869566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(Y, NORM_2, &ynorm));
879566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(X, NORM_2, &xnorm));
889566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(Y, NORM_2, &ynorm));
899566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
903b288129SPeter Brune 
9171b4ebd2SPeter Brune   if (ynorm == 0.0) {
9271b4ebd2SPeter Brune     if (monitor) {
939566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
949566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Initial direction and size is 0\n"));
959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
9671b4ebd2SPeter Brune     }
979566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, W));
989566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, G));
999566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
1009566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
1013ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
10271b4ebd2SPeter Brune   }
1038a430afbSPeter Brune 
1048a430afbSPeter Brune   /* if the SNES has an objective set, use that instead of the function value */
105bd4a8a71SBarry Smith   if (objective) {
1069566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, X, &f));
1078a430afbSPeter Brune   } else {
108e3ce9940SStefano Zampini     f = 0.5 * PetscSqr(fnorm);
1098a430afbSPeter Brune   }
1108a430afbSPeter Brune 
1118a430afbSPeter Brune   /* compute the initial slope */
112bd4a8a71SBarry Smith   if (objective) {
11332661dccSStefano Zampini     /* slope comes from the function (assumed to be the gradient of the objective) */
1149566063dSJacob Faibussowitsch     PetscCall(VecDotRealPart(Y, F, &initslope));
1158a430afbSPeter Brune   } else {
1168a430afbSPeter Brune     /* slope comes from the normal equations */
1179566063dSJacob Faibussowitsch     PetscCall(MatMult(jac, Y, W));
1189566063dSJacob Faibussowitsch     PetscCall(VecDotRealPart(F, W, &initslope));
11971b4ebd2SPeter Brune     if (initslope > 0.0) initslope = -initslope;
12071b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1218a430afbSPeter Brune   }
12271b4ebd2SPeter Brune 
123df019d78SBarry Smith   while (PETSC_TRUE) {
1249566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(W, -lambda, Y, X));
1251baa6e33SBarry Smith     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
126e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
1279566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while checking full step length!\n"));
12871b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1299566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
1303ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
13171b4ebd2SPeter Brune     }
1328a430afbSPeter Brune 
133bd4a8a71SBarry Smith     if (objective) {
1349566063dSJacob Faibussowitsch       PetscCall(SNESComputeObjective(snes, W, &g));
1358a430afbSPeter Brune     } else {
1369566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
1379bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
1389bd66eb0SPeter Brune         gnorm = fnorm;
1399566063dSJacob Faibussowitsch         PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
1409bd66eb0SPeter Brune       } else {
1419566063dSJacob Faibussowitsch         PetscCall(VecNorm(G, NORM_2, &gnorm));
1428a430afbSPeter Brune       }
143e3ce9940SStefano Zampini       g = 0.5 * PetscSqr(gnorm);
1449bd66eb0SPeter Brune     }
1459566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchMonitor(linesearch));
1469bd66eb0SPeter Brune 
147df019d78SBarry Smith     if (!PetscIsInfOrNanReal(g)) break;
148df019d78SBarry Smith     if (monitor) {
1499566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: objective function at lambdas = %g is Inf or Nan, cutting lambda\n", (double)lambda));
1519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
152df019d78SBarry Smith     }
153ad540459SPierre Jolivet     if (lambda <= minlambda) SNESCheckFunctionNorm(snes, g);
15432661dccSStefano Zampini     lambda *= .5;
155df019d78SBarry Smith   }
156df019d78SBarry Smith 
15748a46eb9SPierre Jolivet   if (!objective) PetscCall(PetscInfo(snes, "Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
158e3ce9940SStefano Zampini   if (g <= f + lambda * alpha * initslope) { /* Sufficient reduction or step tolerance convergence */
15971b4ebd2SPeter Brune     if (monitor) {
1609566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
161bd4a8a71SBarry Smith       if (!objective) {
1629566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm));
1638a430afbSPeter Brune       } else {
164e3ce9940SStefano Zampini         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Using full step: old obj %14.12e new obj %14.12e\n", (double)f, (double)g));
1658a430afbSPeter Brune       }
1669566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
16771b4ebd2SPeter Brune     }
16871b4ebd2SPeter Brune   } else {
169c21ba15cSPeter Brune     if (stol * xnorm > ynorm) {
170a68bbae5SBarry Smith       /* Since the full step didn't give sufficient decrease and the step is tiny, exit */
1719566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
1729566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
173c61ba1e2SPeter Brune       if (monitor) {
1749566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
17563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n", (double)ynorm, (double)(stol * xnorm)));
1769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
177c61ba1e2SPeter Brune       }
1783ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
179c21ba15cSPeter Brune     }
18032661dccSStefano Zampini     /* Here to avoid -Wmaybe-uninitiliazed warnings */
18171b4ebd2SPeter Brune     lambdaprev = lambda;
1828a430afbSPeter Brune     gprev      = g;
18332661dccSStefano Zampini     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR) {
18432661dccSStefano Zampini       /* Fit points with quadratic */
185e3ce9940SStefano Zampini       lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
186e3ce9940SStefano Zampini       lambda     = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
18771b4ebd2SPeter Brune 
1889566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(W, -lambda, Y, X));
1891baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
190e71169deSBarry Smith       if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
19163a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while attempting quadratic backtracking! %" PetscInt_FMT " \n", snes->nfuncs));
19271b4ebd2SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1939566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
1943ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
19571b4ebd2SPeter Brune       }
196bd4a8a71SBarry Smith       if (objective) {
1979566063dSJacob Faibussowitsch         PetscCall(SNESComputeObjective(snes, W, &g));
1988a430afbSPeter Brune       } else {
1999566063dSJacob Faibussowitsch         PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
2009bd66eb0SPeter Brune         if (linesearch->ops->vinorm) {
2019bd66eb0SPeter Brune           gnorm = fnorm;
2029566063dSJacob Faibussowitsch           PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
2039bd66eb0SPeter Brune         } else {
2049566063dSJacob Faibussowitsch           PetscCall(VecNorm(G, NORM_2, &gnorm));
2059bd66eb0SPeter Brune         }
206e3ce9940SStefano Zampini         g = 0.5 * PetscSqr(gnorm);
2078a430afbSPeter Brune       }
208c98378a5SBarry Smith       if (PetscIsInfOrNanReal(g)) {
2099566063dSJacob Faibussowitsch         PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
2109566063dSJacob Faibussowitsch         PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
2113ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
212c98378a5SBarry Smith       }
21371b4ebd2SPeter Brune       if (monitor) {
2149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
215bd4a8a71SBarry Smith         if (!objective) {
2169566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: gnorm after quadratic fit %14.12e\n", (double)gnorm));
2178a430afbSPeter Brune         } else {
2189566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: obj after quadratic fit %14.12e\n", (double)g));
2198a430afbSPeter Brune         }
2209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
22171b4ebd2SPeter Brune       }
22232661dccSStefano Zampini     }
223e3ce9940SStefano Zampini     if (linesearch->order != SNES_LINESEARCH_ORDER_LINEAR && g <= f + lambda * alpha * initslope) { /* sufficient reduction */
22471b4ebd2SPeter Brune       if (monitor) {
2259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
2269566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Quadratically determined step, lambda=%18.16e\n", (double)lambda));
2279566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
22871b4ebd2SPeter Brune       }
22971b4ebd2SPeter Brune     } else {
230*a99ef635SJonas Heinzmann       for (count = 0; count < max_it; count++) {
23171b4ebd2SPeter Brune         if (lambda <= minlambda) {
23271b4ebd2SPeter Brune           if (monitor) {
2339566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
23463a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: unable to find good step length! After %" PetscInt_FMT " tries \n", count));
235bd4a8a71SBarry Smith             if (!objective) {
2369371c9d4SSatish Balay               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));
2378a430afbSPeter Brune             } else {
238e3ce9940SStefano Zampini               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));
2398a430afbSPeter Brune             }
2409566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
24171b4ebd2SPeter Brune           }
2429566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
2433ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
24471b4ebd2SPeter Brune         }
245b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
246e3ce9940SStefano Zampini           /* Fit points with cubic */
247e3ce9940SStefano Zampini           t1 = g - f - lambda * initslope;
248e3ce9940SStefano Zampini           t2 = gprev - f - lambdaprev * initslope;
24971b4ebd2SPeter Brune           a  = (t1 / (lambda * lambda) - t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
25071b4ebd2SPeter Brune           b  = (-lambdaprev * t1 / (lambda * lambda) + lambda * t2 / (lambdaprev * lambdaprev)) / (lambda - lambdaprev);
25171b4ebd2SPeter Brune           d  = b * b - 3 * a * initslope;
25271b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
253f5af7f23SKarl Rupp           if (a == 0.0) lambdatemp = -initslope / (2.0 * b);
254f5af7f23SKarl Rupp           else lambdatemp = (-b + PetscSqrtReal(d)) / (3.0 * a);
255b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
256e3ce9940SStefano Zampini           lambdatemp = -initslope * PetscSqr(lambda) / (2.0 * (g - f - lambda * initslope));
25732661dccSStefano Zampini         } else if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) { /* Just backtrack */
25832661dccSStefano Zampini           lambdatemp = .5 * lambda;
25932661dccSStefano Zampini         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "Line search order %" PetscInt_FMT " for type bt", linesearch->order);
26071b4ebd2SPeter Brune         lambdaprev = lambda;
2618a430afbSPeter Brune         gprev      = g;
262e3ce9940SStefano Zampini 
263e3ce9940SStefano Zampini         lambda = PetscClipInterval(lambdatemp, .1 * lambda, .5 * lambda);
2649566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(W, -lambda, Y, X));
2651baa6e33SBarry Smith         if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
266e71169deSBarry Smith         if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
26763a3b9bcSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations, while looking for good step length! %" PetscInt_FMT " \n", count));
26848a46eb9SPierre Jolivet           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));
2699566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
27071b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
2713ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
27271b4ebd2SPeter Brune         }
273bd4a8a71SBarry Smith         if (objective) {
2749566063dSJacob Faibussowitsch           PetscCall(SNESComputeObjective(snes, W, &g));
2758a430afbSPeter Brune         } else {
2769566063dSJacob Faibussowitsch           PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
2779bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
2789bd66eb0SPeter Brune             gnorm = fnorm;
2799566063dSJacob Faibussowitsch             PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
2809bd66eb0SPeter Brune           } else {
2819566063dSJacob Faibussowitsch             PetscCall(VecNorm(G, NORM_2, &gnorm));
2828a430afbSPeter Brune           }
283e3ce9940SStefano Zampini           g = 0.5 * PetscSqr(gnorm);
2849bd66eb0SPeter Brune         }
2854a2f8832SBarry Smith         if (PetscIsInfOrNanReal(g)) {
2869566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
2879566063dSJacob Faibussowitsch           PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
2883ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
289c98378a5SBarry Smith         }
290e3ce9940SStefano Zampini         if (g <= f + lambda * alpha * initslope) { /* is reduction enough? */
29171b4ebd2SPeter Brune           if (monitor) {
2929566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
293bd4a8a71SBarry Smith             if (!objective) {
29432661dccSStefano Zampini               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: %s step, current gnorm %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)gnorm, (double)lambda));
2959566063dSJacob Faibussowitsch               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
2968a430afbSPeter Brune             } else {
29732661dccSStefano Zampini               PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: %s step, obj %14.12e lambda=%18.16e\n", ordStr[linesearch->order - 1], (double)g, (double)lambda));
2989566063dSJacob Faibussowitsch               PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
2998a430afbSPeter Brune             }
30071b4ebd2SPeter Brune           }
30171b4ebd2SPeter Brune           break;
302f5af7f23SKarl Rupp         } else if (monitor) {
3039566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
304bd4a8a71SBarry Smith           if (!objective) {
30532661dccSStefano Zampini             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));
3069566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
3078a430afbSPeter Brune           } else {
30832661dccSStefano Zampini             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));
3099566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
3108a430afbSPeter Brune           }
31171b4ebd2SPeter Brune         }
31271b4ebd2SPeter Brune       }
31371b4ebd2SPeter Brune     }
31471b4ebd2SPeter Brune   }
31571b4ebd2SPeter Brune 
31671b4ebd2SPeter Brune   /* postcheck */
3176b095a85SStefano Zampini   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
3189566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
31971b4ebd2SPeter Brune   if (changed_y) {
3206b095a85SStefano Zampini     if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
3211baa6e33SBarry Smith     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
32271b4ebd2SPeter Brune   }
323bd4a8a71SBarry Smith   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
3249566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
3259bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3269bd66eb0SPeter Brune       gnorm = fnorm;
3279566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->vinorm)(snes, G, W, &gnorm));
3289bd66eb0SPeter Brune     } else {
3299566063dSJacob Faibussowitsch       PetscCall(VecNorm(G, NORM_2, &gnorm));
3309bd66eb0SPeter Brune     }
3319566063dSJacob Faibussowitsch     PetscCall(VecNorm(Y, NORM_2, &ynorm));
332c98378a5SBarry Smith     if (PetscIsInfOrNanReal(gnorm)) {
3339566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
3349566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Aborted due to Nan or Inf in function evaluation\n"));
3353ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
336c98378a5SBarry Smith     }
337c427506bSPeter Brune   }
33871b4ebd2SPeter Brune 
33971b4ebd2SPeter Brune   /* copy the solution over */
3409566063dSJacob Faibussowitsch   PetscCall(VecCopy(W, X));
3419566063dSJacob Faibussowitsch   PetscCall(VecCopy(G, F));
3429566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm));
3439566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm));
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34571b4ebd2SPeter Brune }
34671b4ebd2SPeter Brune 
34766976f2fSJacob Faibussowitsch static PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
348d71ae5a4SJacob Faibussowitsch {
3497f1410a3SPeter Brune   PetscBool          iascii;
350d7ca06c0SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
351075cc632SBarry Smith 
3527f1410a3SPeter Brune   PetscFunctionBegin;
3539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3547f1410a3SPeter Brune   if (iascii) {
355b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n"));
357b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
3589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n"));
3597f1410a3SPeter Brune     }
3609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha));
3617f1410a3SPeter Brune   }
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3637f1410a3SPeter Brune }
3647f1410a3SPeter Brune 
365d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
366d71ae5a4SJacob Faibussowitsch {
36771b4ebd2SPeter Brune   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCall(PetscFree(linesearch->data));
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37071b4ebd2SPeter Brune }
37171b4ebd2SPeter Brune 
372ce78bad3SBarry Smith static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch, PetscOptionItems PetscOptionsObject)
373d71ae5a4SJacob Faibussowitsch {
374eb23ba39SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT *)linesearch->data;
37571b4ebd2SPeter Brune 
3766e111a19SKarl Rupp   PetscFunctionBegin;
377d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNESLineSearch BT options");
3789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL));
379d0609cedSBarry Smith   PetscOptionsHeadEnd();
3803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38171b4ebd2SPeter Brune }
38271b4ebd2SPeter Brune 
383954494b2SJed Brown /*MC
3841d27aa22SBarry Smith    SNESLINESEARCHBT - Backtracking line search {cite}`dennis:83`.
385954494b2SJed Brown 
386*a99ef635SJonas Heinzmann    This line search finds the minimum of a polynomial fitting either $1/2 ||F(x_k + \lambda Y_k)||_2^2$,
387*a99ef635SJonas Heinzmann    or the objective function $G(x_k + \lambda Y_k)$ if it is provided with `SNESSetObjective()`.
388*a99ef635SJonas Heinzmann    If this fit does not satisfy the sufficient decrease conditions, the interval shrinks
38989dfbbd5SBarry Smith    and the fit is reattempted at most `max_it` times or until $\lambda$ is below `minlambda`.
390954494b2SJed Brown 
391cd7522eaSPeter Brune    Options Database Keys:
3923eec2f6aSPierre Jolivet +  -snes_linesearch_alpha <1e\-4>      - slope descent parameter
393*a99ef635SJonas Heinzmann .  -snes_linesearch_damping <1.0>      - initial `lambda` on entry to the line search
394*a99ef635SJonas Heinzmann .  -snes_linesearch_max_it <40>        - maximum number of shrinking iterations in the line search
395*a99ef635SJonas Heinzmann .  -snes_linesearch_minlambda <1e\-12> - minimum `lambda` (scaling of solution update) allowed
396*a99ef635SJonas Heinzmann -  -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
397cd7522eaSPeter Brune 
398954494b2SJed Brown    Level: advanced
399954494b2SJed Brown 
400f6dfbefdSBarry Smith    Note:
401e55accfeSBarryFSmith    This line search will always produce a step that is less than or equal to, in length, the full step size.
402db609ea7SPeter Brune 
403420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
404954494b2SJed Brown M*/
405d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
406d71ae5a4SJacob Faibussowitsch {
407f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
40871b4ebd2SPeter Brune 
40971b4ebd2SPeter Brune   PetscFunctionBegin;
410f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
411f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
412f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
4130298fd71SBarry Smith   linesearch->ops->reset          = NULL;
4147f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
4150298fd71SBarry Smith   linesearch->ops->setup          = NULL;
41671b4ebd2SPeter Brune 
4174dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&bt));
418f5af7f23SKarl Rupp 
41971b4ebd2SPeter Brune   linesearch->data   = (void *)bt;
420*a99ef635SJonas Heinzmann   linesearch->max_it = 40;
421b000cd8dSPeter Brune   linesearch->order  = SNES_LINESEARCH_ORDER_CUBIC;
42271b4ebd2SPeter Brune   bt->alpha          = 1e-4;
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42471b4ebd2SPeter Brune }
425