xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision 3eec2f6aa40ee6962d3fbcd0607a6f3de30c0607)
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 /*@
92f4102e2SPeter Brune    SNESLineSearchBTSetAlpha - Sets the descent parameter, alpha, in the BT linesearch 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 
171a4f838cSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
182f4102e2SPeter Brune @*/
192f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha)
202f4102e2SPeter Brune {
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;
262f4102e2SPeter Brune   PetscFunctionReturn(0);
272f4102e2SPeter Brune }
282f4102e2SPeter Brune 
292f4102e2SPeter Brune /*@
302f4102e2SPeter Brune    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant.
312f4102e2SPeter Brune 
322f4102e2SPeter Brune    Input Parameters:
332f4102e2SPeter Brune .  linesearch - linesearch context
348e557f58SPeter Brune 
358e557f58SPeter Brune    Output Parameters:
362f4102e2SPeter Brune .  alpha - The descent parameter
372f4102e2SPeter Brune 
382f4102e2SPeter Brune    Level: intermediate
392f4102e2SPeter Brune 
401a4f838cSPeter Brune .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
412f4102e2SPeter Brune @*/
422f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
432f4102e2SPeter Brune {
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;
492f4102e2SPeter Brune   PetscFunctionReturn(0);
502f4102e2SPeter Brune }
512f4102e2SPeter Brune 
52f1c6b773SPeter Brune static PetscErrorCode  SNESLineSearchApply_BT(SNESLineSearch linesearch)
5371b4ebd2SPeter Brune {
5471b4ebd2SPeter Brune   PetscBool         changed_y,changed_w;
5571b4ebd2SPeter Brune   PetscErrorCode    ierr;
5671b4ebd2SPeter Brune   Vec               X,F,Y,W,G;
5771b4ebd2SPeter Brune   SNES              snes;
588a430afbSPeter Brune   PetscReal         fnorm, xnorm, ynorm, gnorm;
593b288129SPeter Brune   PetscReal         lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol;
60dfcd3828SPeter Brune   PetscReal         t1,t2,a,b,d;
618a430afbSPeter Brune   PetscReal         f;
628a430afbSPeter Brune   PetscReal         g,gprev;
6371b4ebd2SPeter Brune   PetscViewer       monitor;
6471b4ebd2SPeter Brune   PetscInt          max_its,count;
65d7ca06c0SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT*)linesearch->data;
6671b4ebd2SPeter Brune   Mat               jac;
67bd4a8a71SBarry Smith   PetscErrorCode    (*objective)(SNES,Vec,PetscReal*,void*);
6871b4ebd2SPeter Brune 
6971b4ebd2SPeter Brune   PetscFunctionBegin;
70f1c6b773SPeter Brune   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
71f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
72f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
73f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
74dcf2fd19SBarry Smith   ierr = SNESLineSearchGetDefaultMonitor(linesearch, &monitor);CHKERRQ(ierr);
750298fd71SBarry Smith   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,NULL,NULL,NULL,&max_its);CHKERRQ(ierr);
760298fd71SBarry Smith   ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr);
770298fd71SBarry Smith   ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr);
7871b4ebd2SPeter Brune   alpha = bt->alpha;
7971b4ebd2SPeter Brune 
800298fd71SBarry Smith   ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr);
81ce94432eSBarry Smith   if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
82c69d1a72SBarry Smith 
837b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
84422a814eSBarry Smith   ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
8571b4ebd2SPeter Brune 
863b288129SPeter Brune   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
873b288129SPeter Brune   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
883b288129SPeter Brune   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
893b288129SPeter Brune   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
903b288129SPeter Brune 
9171b4ebd2SPeter Brune   if (ynorm == 0.0) {
9271b4ebd2SPeter Brune     if (monitor) {
93ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
9471b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
95ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
9671b4ebd2SPeter Brune     }
9771b4ebd2SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
9871b4ebd2SPeter Brune     ierr = VecCopy(F,G);CHKERRQ(ierr);
9971cee68bSPeter Brune     ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
100e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
10171b4ebd2SPeter Brune     PetscFunctionReturn(0);
10271b4ebd2SPeter Brune   }
10371b4ebd2SPeter Brune   if (ynorm > maxstep) {        /* Step too big, so scale back */
10471b4ebd2SPeter Brune     if (monitor) {
105ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
106c69d1a72SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr);
107ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10871b4ebd2SPeter Brune     }
10971b4ebd2SPeter Brune     ierr  = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
11071b4ebd2SPeter Brune     ynorm = maxstep;
11171b4ebd2SPeter Brune   }
1128a430afbSPeter Brune 
1138a430afbSPeter Brune   /* if the SNES has an objective set, use that instead of the function value */
114bd4a8a71SBarry Smith   if (objective) {
1158a430afbSPeter Brune     ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr);
1168a430afbSPeter Brune   } else {
1178a430afbSPeter Brune     f = fnorm*fnorm;
1188a430afbSPeter Brune   }
1198a430afbSPeter Brune 
1208a430afbSPeter Brune   /* compute the initial slope */
121bd4a8a71SBarry Smith   if (objective) {
1228a430afbSPeter Brune     /* slope comes from the function (assumed to be the gradient of the objective */
12367392de3SBarry Smith     ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr);
1248a430afbSPeter Brune   } else {
1258a430afbSPeter Brune     /* slope comes from the normal equations */
12671b4ebd2SPeter Brune     ierr = MatMult(jac,Y,W);CHKERRQ(ierr);
12767392de3SBarry Smith     ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr);
12871b4ebd2SPeter Brune     if (initslope > 0.0)  initslope = -initslope;
12971b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1308a430afbSPeter Brune   }
13171b4ebd2SPeter Brune 
132df019d78SBarry Smith   while (PETSC_TRUE) {
133c427506bSPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1349bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1359bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1369bd66eb0SPeter Brune     }
137e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
13871b4ebd2SPeter Brune       ierr         = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
13971b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
140e9b602ebSSatish Balay       ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
14171b4ebd2SPeter Brune       PetscFunctionReturn(0);
14271b4ebd2SPeter Brune     }
1438a430afbSPeter Brune 
144bd4a8a71SBarry Smith     if (objective) {
1458a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
1468a430afbSPeter Brune     } else {
147ce5a860aSPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
1489bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
1499bd66eb0SPeter Brune         gnorm = fnorm;
1509bd66eb0SPeter Brune         ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1519bd66eb0SPeter Brune       } else {
15271b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1538a430afbSPeter Brune       }
15479e7e81bSJed Brown       g = PetscSqr(gnorm);
1559bd66eb0SPeter Brune     }
156dcf2fd19SBarry Smith     ierr = SNESLineSearchMonitor(linesearch);CHKERRQ(ierr);
1579bd66eb0SPeter Brune 
158df019d78SBarry Smith     if (!PetscIsInfOrNanReal(g)) break;
159df019d78SBarry Smith     if (monitor) {
160df019d78SBarry Smith       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
161df019d78SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: objective function at lambdas = %g is Inf or Nan, cutting lambda\n",(double)lambda);CHKERRQ(ierr);
162df019d78SBarry Smith       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
163df019d78SBarry Smith     }
164df019d78SBarry Smith     if (lambda <= minlambda) {
1657e49c775SBarry Smith       SNESCheckFunctionNorm(snes,g);
166c98378a5SBarry Smith     }
167df019d78SBarry Smith     lambda = .5*lambda;
168df019d78SBarry Smith   }
169df019d78SBarry Smith 
170bd4a8a71SBarry Smith   if (!objective) {
171c69d1a72SBarry Smith     ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1728a430afbSPeter Brune   }
1738a430afbSPeter Brune   if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */
17471b4ebd2SPeter Brune     if (monitor) {
175ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
176bd4a8a71SBarry Smith       if (!objective) {
177c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1788a430afbSPeter Brune       } else {
179df019d78SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: obj %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr);
1808a430afbSPeter Brune       }
181ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
18271b4ebd2SPeter Brune     }
18371b4ebd2SPeter Brune   } else {
184c21ba15cSPeter Brune     /* Since the full step didn't work and the step is tiny, quit */
185c21ba15cSPeter Brune     if (stol*xnorm > ynorm) {
186ced04f9dSPeter Brune       ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
187efa57af7SBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
188c61ba1e2SPeter Brune       if (monitor) {
189c61ba1e2SPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
190efa57af7SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n",(double)ynorm,(double)stol*xnorm);CHKERRQ(ierr);
191c61ba1e2SPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
192c61ba1e2SPeter Brune       }
193c21ba15cSPeter Brune       PetscFunctionReturn(0);
194c21ba15cSPeter Brune     }
19571b4ebd2SPeter Brune     /* Fit points with quadratic */
1968a430afbSPeter Brune     lambdatemp = -initslope/(g - f - 2.0*lambda*initslope);
19771b4ebd2SPeter Brune     lambdaprev = lambda;
1988a430afbSPeter Brune     gprev      = g;
19971b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
20071b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
20171b4ebd2SPeter Brune     else                         lambda = lambdatemp;
20271b4ebd2SPeter Brune 
20371b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2049bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2059bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2069bd66eb0SPeter Brune     }
207e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
20871b4ebd2SPeter Brune       ierr         = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
20971b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
210e9b602ebSSatish Balay       ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
21171b4ebd2SPeter Brune       PetscFunctionReturn(0);
21271b4ebd2SPeter Brune     }
213bd4a8a71SBarry Smith     if (objective) {
2148a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2158a430afbSPeter Brune     } else {
216ce5a860aSPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
2179bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
2189bd66eb0SPeter Brune         gnorm = fnorm;
2199bd66eb0SPeter Brune         ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2209bd66eb0SPeter Brune       } else {
22171b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2229bd66eb0SPeter Brune       }
223f1ab69e3SPeter Brune       g = PetscSqr(gnorm);
2248a430afbSPeter Brune     }
225c98378a5SBarry Smith     if (PetscIsInfOrNanReal(g)) {
226422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
227c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
228c98378a5SBarry Smith       PetscFunctionReturn(0);
229c98378a5SBarry Smith     }
23071b4ebd2SPeter Brune     if (monitor) {
231ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
232bd4a8a71SBarry Smith       if (!objective) {
233c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr);
2348a430afbSPeter Brune       } else {
2358a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr);
2368a430afbSPeter Brune       }
237ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
23871b4ebd2SPeter Brune     }
2398a430afbSPeter Brune     if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */
24071b4ebd2SPeter Brune       if (monitor) {
241ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24271b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
243ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24471b4ebd2SPeter Brune       }
24571b4ebd2SPeter Brune     } else {
24671b4ebd2SPeter Brune       /* Fit points with cubic */
24771b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
24871b4ebd2SPeter Brune         if (lambda <= minlambda) {
24971b4ebd2SPeter Brune           if (monitor) {
250ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
25171b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
252bd4a8a71SBarry Smith             if (!objective) {
2534a2f8832SBarry Smith               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
254c69d1a72SBarry Smith                                                          (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2558a430afbSPeter Brune             } else {
2564a2f8832SBarry Smith               ierr = 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",
2578a430afbSPeter Brune                                                          (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2588a430afbSPeter Brune             }
259ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26071b4ebd2SPeter Brune           }
261e9b602ebSSatish Balay           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
26271b4ebd2SPeter Brune           PetscFunctionReturn(0);
26371b4ebd2SPeter Brune         }
264b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
2658a430afbSPeter Brune           t1 = .5*(g - f) - lambda*initslope;
2668a430afbSPeter Brune           t2 = .5*(gprev  - f) - lambdaprev*initslope;
26771b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
26871b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
26971b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
27071b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
271f5af7f23SKarl Rupp           if (a == 0.0) lambdatemp = -initslope/(2.0*b);
272f5af7f23SKarl Rupp           else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
273f5af7f23SKarl Rupp 
274b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
2758a430afbSPeter Brune           lambdatemp = -initslope/(g - f - 2.0*initslope);
276ce94432eSBarry Smith         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt");
27771b4ebd2SPeter Brune         lambdaprev = lambda;
2788a430afbSPeter Brune         gprev      = g;
27971b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
28071b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
28171b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
28271b4ebd2SPeter Brune         ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
283b7b2e573SPeter Brune         if (linesearch->ops->viproject) {
284b7b2e573SPeter Brune           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
285b7b2e573SPeter Brune         }
286e71169deSBarry Smith         if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
28771b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
288bd4a8a71SBarry Smith           if (!objective) {
28971b4ebd2SPeter Brune             ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
290c69d1a72SBarry Smith                               (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
2918a430afbSPeter Brune           }
292e9b602ebSSatish Balay           ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
29371b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
29471b4ebd2SPeter Brune           PetscFunctionReturn(0);
29571b4ebd2SPeter Brune         }
296bd4a8a71SBarry Smith         if (objective) {
2978a430afbSPeter Brune           ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2988a430afbSPeter Brune         } else {
299ce5a860aSPeter Brune           ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3009bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
3019bd66eb0SPeter Brune             gnorm = fnorm;
3029bd66eb0SPeter Brune             ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3039bd66eb0SPeter Brune           } else {
30471b4ebd2SPeter Brune             ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3058a430afbSPeter Brune           }
306f1ab69e3SPeter Brune           g = PetscSqr(gnorm);
3079bd66eb0SPeter Brune         }
3084a2f8832SBarry Smith         if (PetscIsInfOrNanReal(g)) {
309422a814eSBarry Smith           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
310c61ba1e2SPeter Brune           ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
311c98378a5SBarry Smith           PetscFunctionReturn(0);
312c98378a5SBarry Smith         }
3138a430afbSPeter Brune         if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */
31471b4ebd2SPeter Brune           if (monitor) {
315ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
316bd4a8a71SBarry Smith             if (!objective) {
317b000cd8dSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
318c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
31971b4ebd2SPeter Brune               } else {
320c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32171b4ebd2SPeter Brune               }
322ea5d4fccSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3238a430afbSPeter Brune             } else {
3248a430afbSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3258a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3268a430afbSPeter Brune               } else {
3278a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3288a430afbSPeter Brune               }
3298a430afbSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3308a430afbSPeter Brune             }
33171b4ebd2SPeter Brune           }
33271b4ebd2SPeter Brune           break;
333f5af7f23SKarl Rupp         } else if (monitor) {
334ea5d4fccSPeter Brune           ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
335bd4a8a71SBarry Smith           if (!objective) {
336b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
337c69d1a72SBarry Smith               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
33871b4ebd2SPeter Brune             } else {
339c69d1a72SBarry Smith               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
34071b4ebd2SPeter Brune             }
341ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3428a430afbSPeter Brune           } else {
3438a430afbSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3448a430afbSPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3458a430afbSPeter Brune             } else {
3468a430afbSPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3478a430afbSPeter Brune             }
3488a430afbSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3498a430afbSPeter Brune           }
35071b4ebd2SPeter Brune         }
35171b4ebd2SPeter Brune       }
35271b4ebd2SPeter Brune     }
35371b4ebd2SPeter Brune   }
35471b4ebd2SPeter Brune 
35571b4ebd2SPeter Brune   /* postcheck */
3560c1ac483SGlenn Hammond   /* update Y to lambda*Y so that W is consistent with  X - lambda*Y */
3570c1ac483SGlenn Hammond   ierr = VecScale(Y,lambda);CHKERRQ(ierr);
3587b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
35971b4ebd2SPeter Brune   if (changed_y) {
3600c1ac483SGlenn Hammond     ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);
3619bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
3629bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
3639bd66eb0SPeter Brune     }
36471b4ebd2SPeter Brune   }
365bd4a8a71SBarry Smith   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
366ce5a860aSPeter Brune     ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3679bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3689bd66eb0SPeter Brune       gnorm = fnorm;
3699bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3709bd66eb0SPeter Brune     } else {
3719bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3729bd66eb0SPeter Brune     }
3739bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
374c98378a5SBarry Smith     if (PetscIsInfOrNanReal(gnorm)) {
375422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
376c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
377c98378a5SBarry Smith       PetscFunctionReturn(0);
378c98378a5SBarry Smith     }
379c427506bSPeter Brune   }
38071b4ebd2SPeter Brune 
38171b4ebd2SPeter Brune   /* copy the solution over */
38271b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
383c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
384c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
385f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
386f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
38771b4ebd2SPeter Brune   PetscFunctionReturn(0);
38871b4ebd2SPeter Brune }
38971b4ebd2SPeter Brune 
3907f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
3917f1410a3SPeter Brune {
3927f1410a3SPeter Brune   PetscErrorCode    ierr;
3937f1410a3SPeter Brune   PetscBool         iascii;
394d7ca06c0SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT*)linesearch->data;
395075cc632SBarry Smith 
3967f1410a3SPeter Brune   PetscFunctionBegin;
397251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
3987f1410a3SPeter Brune   if (iascii) {
399b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
4007f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
401b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
4027f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
4037f1410a3SPeter Brune     }
4047904a332SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha);CHKERRQ(ierr);
4057f1410a3SPeter Brune   }
4067f1410a3SPeter Brune   PetscFunctionReturn(0);
4077f1410a3SPeter Brune }
4087f1410a3SPeter Brune 
409f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
410731c067cSJed Brown {
411731c067cSJed Brown   PetscErrorCode ierr;
41271b4ebd2SPeter Brune 
41371b4ebd2SPeter Brune   PetscFunctionBegin;
414731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
41571b4ebd2SPeter Brune   PetscFunctionReturn(0);
41671b4ebd2SPeter Brune }
41771b4ebd2SPeter Brune 
4184416b707SBarry Smith static PetscErrorCode SNESLineSearchSetFromOptions_BT(PetscOptionItems *PetscOptionsObject,SNESLineSearch linesearch)
41971b4ebd2SPeter Brune {
42071b4ebd2SPeter Brune   PetscErrorCode    ierr;
421eb23ba39SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT*)linesearch->data;
42271b4ebd2SPeter Brune 
4236e111a19SKarl Rupp   PetscFunctionBegin;
424e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESLineSearch BT options");CHKERRQ(ierr);
4250298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr);
42671b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
42771b4ebd2SPeter Brune   PetscFunctionReturn(0);
42871b4ebd2SPeter Brune }
42971b4ebd2SPeter Brune 
430954494b2SJed Brown /*MC
431db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
432954494b2SJed Brown 
433db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
434eb23ba39SBarry Smith    function or the objective function if it is provided with SNESSetObjective(). If this fit does not satisfy the conditions for progress, the interval shrinks
435db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
436954494b2SJed Brown 
437cd7522eaSPeter Brune    Options Database Keys:
438*3eec2f6aSPierre Jolivet +  -snes_linesearch_alpha <1e\-4> - slope descent parameter
439db609ea7SPeter Brune .  -snes_linesearch_damping <1.0> - initial step length
440eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - if the length the initial step is larger than this then the
441e1bc860dSBarry Smith                                        step is scaled back to be of this length at the beginning of the line search
442db609ea7SPeter Brune .  -snes_linesearch_max_it <40> - maximum number of shrinking step
443*3eec2f6aSPierre Jolivet .  -snes_linesearch_minlambda <1e\-12> - minimum step length allowed
444cd7522eaSPeter Brune -  -snes_linesearch_order <cubic,quadratic> - order of the approximation
445cd7522eaSPeter Brune 
446954494b2SJed Brown    Level: advanced
447954494b2SJed Brown 
448db609ea7SPeter Brune    Notes:
449db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
450db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
451e55accfeSBarryFSmith 
452e55accfeSBarryFSmith    This line search will always produce a step that is less than or equal to, in length, the full step size.
453db609ea7SPeter Brune 
454f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
455954494b2SJed Brown M*/
4568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
45771b4ebd2SPeter Brune {
45871b4ebd2SPeter Brune 
459f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
46071b4ebd2SPeter Brune   PetscErrorCode    ierr;
46171b4ebd2SPeter Brune 
46271b4ebd2SPeter Brune   PetscFunctionBegin;
463f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
464f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
465f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
4660298fd71SBarry Smith   linesearch->ops->reset          = NULL;
4677f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
4680298fd71SBarry Smith   linesearch->ops->setup          = NULL;
46971b4ebd2SPeter Brune 
470b00a9115SJed Brown   ierr = PetscNewLog(linesearch,&bt);CHKERRQ(ierr);
471f5af7f23SKarl Rupp 
47271b4ebd2SPeter Brune   linesearch->data    = (void*)bt;
47371b4ebd2SPeter Brune   linesearch->max_its = 40;
474b000cd8dSPeter Brune   linesearch->order   = SNES_LINESEARCH_ORDER_CUBIC;
47571b4ebd2SPeter Brune   bt->alpha           = 1e-4;
47671b4ebd2SPeter Brune   PetscFunctionReturn(0);
47771b4ebd2SPeter Brune }
478