xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision efa57af7a86939c2389f3dd5c5e4e558a43c2690)
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 {
212f4102e2SPeter Brune   SNESLineSearch_BT *bt;
226e111a19SKarl Rupp 
232f4102e2SPeter Brune   PetscFunctionBegin;
242f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
252f4102e2SPeter Brune   bt        = (SNESLineSearch_BT*)linesearch->data;
262f4102e2SPeter Brune   bt->alpha = alpha;
272f4102e2SPeter Brune   PetscFunctionReturn(0);
282f4102e2SPeter Brune }
292f4102e2SPeter Brune 
302f4102e2SPeter Brune 
312f4102e2SPeter Brune /*@
322f4102e2SPeter Brune    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant.
332f4102e2SPeter Brune 
342f4102e2SPeter Brune    Input Parameters:
352f4102e2SPeter Brune .  linesearch - linesearch context
368e557f58SPeter Brune 
378e557f58SPeter Brune    Output Parameters:
382f4102e2SPeter Brune .  alpha - The descent parameter
392f4102e2SPeter Brune 
402f4102e2SPeter Brune    Level: intermediate
412f4102e2SPeter Brune 
421a4f838cSPeter Brune .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
432f4102e2SPeter Brune @*/
442f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
452f4102e2SPeter Brune {
462f4102e2SPeter Brune   SNESLineSearch_BT *bt;
476e111a19SKarl Rupp 
482f4102e2SPeter Brune   PetscFunctionBegin;
492f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
502f4102e2SPeter Brune   bt     = (SNESLineSearch_BT*)linesearch->data;
512f4102e2SPeter Brune   *alpha = bt->alpha;
522f4102e2SPeter Brune   PetscFunctionReturn(0);
532f4102e2SPeter Brune }
542f4102e2SPeter Brune 
55f1c6b773SPeter Brune static PetscErrorCode  SNESLineSearchApply_BT(SNESLineSearch linesearch)
5671b4ebd2SPeter Brune {
5771b4ebd2SPeter Brune   PetscBool         changed_y,changed_w;
5871b4ebd2SPeter Brune   PetscErrorCode    ierr;
5971b4ebd2SPeter Brune   Vec               X,F,Y,W,G;
6071b4ebd2SPeter Brune   SNES              snes;
618a430afbSPeter Brune   PetscReal         fnorm, xnorm, ynorm, gnorm;
623b288129SPeter Brune   PetscReal         lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol;
63dfcd3828SPeter Brune   PetscReal         t1,t2,a,b,d;
648a430afbSPeter Brune   PetscReal         f;
658a430afbSPeter Brune   PetscReal         g,gprev;
6671b4ebd2SPeter Brune   PetscViewer       monitor;
6771b4ebd2SPeter Brune   PetscInt          max_its,count;
68f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
6971b4ebd2SPeter Brune   Mat               jac;
70bd4a8a71SBarry Smith   PetscErrorCode    (*objective)(SNES,Vec,PetscReal*,void*);
7171b4ebd2SPeter Brune 
7271b4ebd2SPeter Brune   PetscFunctionBegin;
73f1c6b773SPeter Brune   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
74f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
75f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
76f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
77dcf2fd19SBarry Smith   ierr = SNESLineSearchGetDefaultMonitor(linesearch, &monitor);CHKERRQ(ierr);
780298fd71SBarry Smith   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,NULL,NULL,NULL,&max_its);CHKERRQ(ierr);
790298fd71SBarry Smith   ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr);
800298fd71SBarry Smith   ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr);
81f1c6b773SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
8271b4ebd2SPeter Brune   alpha = bt->alpha;
8371b4ebd2SPeter Brune 
840298fd71SBarry Smith   ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr);
85ce94432eSBarry Smith   if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
86c69d1a72SBarry Smith 
877b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
88422a814eSBarry Smith   ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
8971b4ebd2SPeter Brune 
903b288129SPeter Brune   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
913b288129SPeter Brune   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
923b288129SPeter Brune   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
933b288129SPeter Brune   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
943b288129SPeter Brune 
9571b4ebd2SPeter Brune   if (ynorm == 0.0) {
9671b4ebd2SPeter Brune     if (monitor) {
97ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
9871b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
99ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10071b4ebd2SPeter Brune     }
10171b4ebd2SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
10271b4ebd2SPeter Brune     ierr = VecCopy(F,G);CHKERRQ(ierr);
10371cee68bSPeter Brune     ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
104e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
10571b4ebd2SPeter Brune     PetscFunctionReturn(0);
10671b4ebd2SPeter Brune   }
10771b4ebd2SPeter Brune   if (ynorm > maxstep) {        /* Step too big, so scale back */
10871b4ebd2SPeter Brune     if (monitor) {
109ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
110c69d1a72SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr);
111ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
11271b4ebd2SPeter Brune     }
11371b4ebd2SPeter Brune     ierr  = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
11471b4ebd2SPeter Brune     ynorm = maxstep;
11571b4ebd2SPeter Brune   }
1168a430afbSPeter Brune 
1178a430afbSPeter Brune   /* if the SNES has an objective set, use that instead of the function value */
118bd4a8a71SBarry Smith   if (objective) {
1198a430afbSPeter Brune     ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr);
1208a430afbSPeter Brune   } else {
1218a430afbSPeter Brune     f = fnorm*fnorm;
1228a430afbSPeter Brune   }
1238a430afbSPeter Brune 
1248a430afbSPeter Brune   /* compute the initial slope */
125bd4a8a71SBarry Smith   if (objective) {
1268a430afbSPeter Brune     /* slope comes from the function (assumed to be the gradient of the objective */
12767392de3SBarry Smith     ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr);
1288a430afbSPeter Brune   } else {
1298a430afbSPeter Brune     /* slope comes from the normal equations */
13071b4ebd2SPeter Brune     ierr = MatMult(jac,Y,W);CHKERRQ(ierr);
13167392de3SBarry Smith     ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr);
13271b4ebd2SPeter Brune     if (initslope > 0.0)  initslope = -initslope;
13371b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1348a430afbSPeter Brune   }
13571b4ebd2SPeter Brune 
136df019d78SBarry Smith   while (PETSC_TRUE) {
137c427506bSPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1389bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1399bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1409bd66eb0SPeter Brune     }
14171b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
14271b4ebd2SPeter Brune       ierr         = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
14371b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
144e9b602ebSSatish Balay       ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
14571b4ebd2SPeter Brune       PetscFunctionReturn(0);
14671b4ebd2SPeter Brune     }
1478a430afbSPeter Brune 
148bd4a8a71SBarry Smith     if (objective) {
1498a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
1508a430afbSPeter Brune     } else {
151ce5a860aSPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
1529bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
1539bd66eb0SPeter Brune         gnorm = fnorm;
1549bd66eb0SPeter Brune         ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1559bd66eb0SPeter Brune       } else {
15671b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1578a430afbSPeter Brune       }
15879e7e81bSJed Brown       g = PetscSqr(gnorm);
1599bd66eb0SPeter Brune     }
160dcf2fd19SBarry Smith     ierr = SNESLineSearchMonitor(linesearch);CHKERRQ(ierr);
1619bd66eb0SPeter Brune 
162df019d78SBarry Smith     if (!PetscIsInfOrNanReal(g)) break;
163df019d78SBarry Smith     if (monitor) {
164df019d78SBarry Smith       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
165df019d78SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: objective function at lambdas = %g is Inf or Nan, cutting lambda\n",(double)lambda);CHKERRQ(ierr);
166df019d78SBarry Smith       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
167df019d78SBarry Smith     }
168df019d78SBarry Smith     if (lambda <= minlambda) {
169df019d78SBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
170c98378a5SBarry Smith       PetscFunctionReturn(0);
171c98378a5SBarry Smith     }
172df019d78SBarry Smith     lambda = .5*lambda;
173df019d78SBarry Smith   }
174df019d78SBarry Smith 
175bd4a8a71SBarry Smith   if (!objective) {
176c69d1a72SBarry Smith     ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1778a430afbSPeter Brune   }
1788a430afbSPeter Brune   if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */
17971b4ebd2SPeter Brune     if (monitor) {
180ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
181bd4a8a71SBarry Smith       if (!objective) {
182c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1838a430afbSPeter Brune       } else {
184df019d78SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: obj %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr);
1858a430afbSPeter Brune       }
186ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
18771b4ebd2SPeter Brune     }
18871b4ebd2SPeter Brune   } else {
189c21ba15cSPeter Brune     /* Since the full step didn't work and the step is tiny, quit */
190c21ba15cSPeter Brune     if (stol*xnorm > ynorm) {
191ced04f9dSPeter Brune       ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
192*efa57af7SBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
193c61ba1e2SPeter Brune       if (monitor) {
194c61ba1e2SPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
195*efa57af7SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Ended due to ynorm < stol*xnorm (%14.12e < %14.12e).\n",(double)ynorm,(double)stol*xnorm);CHKERRQ(ierr);
196c61ba1e2SPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
197c61ba1e2SPeter Brune       }
198c21ba15cSPeter Brune       PetscFunctionReturn(0);
199c21ba15cSPeter Brune     }
20071b4ebd2SPeter Brune     /* Fit points with quadratic */
2018a430afbSPeter Brune     lambdatemp = -initslope/(g - f - 2.0*lambda*initslope);
20271b4ebd2SPeter Brune     lambdaprev = lambda;
2038a430afbSPeter Brune     gprev      = g;
20471b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
20571b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
20671b4ebd2SPeter Brune     else                         lambda = lambdatemp;
20771b4ebd2SPeter Brune 
20871b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2099bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2109bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2119bd66eb0SPeter Brune     }
21271b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
21371b4ebd2SPeter Brune       ierr         = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
21471b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
215e9b602ebSSatish Balay       ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
21671b4ebd2SPeter Brune       PetscFunctionReturn(0);
21771b4ebd2SPeter Brune     }
218bd4a8a71SBarry Smith     if (objective) {
2198a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2208a430afbSPeter Brune     } else {
221ce5a860aSPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
2229bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
2239bd66eb0SPeter Brune         gnorm = fnorm;
2249bd66eb0SPeter Brune         ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2259bd66eb0SPeter Brune       } else {
22671b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2279bd66eb0SPeter Brune       }
228f1ab69e3SPeter Brune       g = PetscSqr(gnorm);
2298a430afbSPeter Brune     }
230c98378a5SBarry Smith     if (PetscIsInfOrNanReal(g)) {
231422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
232c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
233c98378a5SBarry Smith       PetscFunctionReturn(0);
234c98378a5SBarry Smith     }
23571b4ebd2SPeter Brune     if (monitor) {
236ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
237bd4a8a71SBarry Smith       if (!objective) {
238c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr);
2398a430afbSPeter Brune       } else {
2408a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr);
2418a430afbSPeter Brune       }
242ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24371b4ebd2SPeter Brune     }
2448a430afbSPeter Brune     if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */
24571b4ebd2SPeter Brune       if (monitor) {
246ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24771b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
248ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24971b4ebd2SPeter Brune       }
25071b4ebd2SPeter Brune     } else {
25171b4ebd2SPeter Brune       /* Fit points with cubic */
25271b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
25371b4ebd2SPeter Brune         if (lambda <= minlambda) {
25471b4ebd2SPeter Brune           if (monitor) {
255ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
25671b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
257bd4a8a71SBarry Smith             if (!objective) {
2584a2f8832SBarry 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",
259c69d1a72SBarry Smith                                                          (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2608a430afbSPeter Brune             } else {
2614a2f8832SBarry 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",
2628a430afbSPeter Brune                                                          (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2638a430afbSPeter Brune             }
264ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26571b4ebd2SPeter Brune           }
266e9b602ebSSatish Balay           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
26771b4ebd2SPeter Brune           PetscFunctionReturn(0);
26871b4ebd2SPeter Brune         }
269b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
2708a430afbSPeter Brune           t1 = .5*(g - f) - lambda*initslope;
2718a430afbSPeter Brune           t2 = .5*(gprev  - f) - lambdaprev*initslope;
27271b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27371b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27471b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
27571b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
276f5af7f23SKarl Rupp           if (a == 0.0) lambdatemp = -initslope/(2.0*b);
277f5af7f23SKarl Rupp           else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
278f5af7f23SKarl Rupp 
279b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
2808a430afbSPeter Brune           lambdatemp = -initslope/(g - f - 2.0*initslope);
281ce94432eSBarry Smith         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt");
28271b4ebd2SPeter Brune         lambdaprev = lambda;
2838a430afbSPeter Brune         gprev      = g;
28471b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
28571b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
28671b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
28771b4ebd2SPeter Brune         ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
288b7b2e573SPeter Brune         if (linesearch->ops->viproject) {
289b7b2e573SPeter Brune           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
290b7b2e573SPeter Brune         }
29171b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
29271b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
293bd4a8a71SBarry Smith           if (!objective) {
29471b4ebd2SPeter Brune             ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
295c69d1a72SBarry Smith                               (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
2968a430afbSPeter Brune           }
297e9b602ebSSatish Balay           ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
29871b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
29971b4ebd2SPeter Brune           PetscFunctionReturn(0);
30071b4ebd2SPeter Brune         }
301bd4a8a71SBarry Smith         if (objective) {
3028a430afbSPeter Brune           ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
3038a430afbSPeter Brune         } else {
304ce5a860aSPeter Brune           ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3059bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
3069bd66eb0SPeter Brune             gnorm = fnorm;
3079bd66eb0SPeter Brune             ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3089bd66eb0SPeter Brune           } else {
30971b4ebd2SPeter Brune             ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3108a430afbSPeter Brune           }
311f1ab69e3SPeter Brune           g = PetscSqr(gnorm);
3129bd66eb0SPeter Brune         }
3134a2f8832SBarry Smith         if (PetscIsInfOrNanReal(g)) {
314422a814eSBarry Smith           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
315c61ba1e2SPeter Brune           ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
316c98378a5SBarry Smith           PetscFunctionReturn(0);
317c98378a5SBarry Smith         }
3188a430afbSPeter Brune         if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */
31971b4ebd2SPeter Brune           if (monitor) {
320ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
321bd4a8a71SBarry Smith             if (!objective) {
322b000cd8dSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
323c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32471b4ebd2SPeter Brune               } else {
325c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32671b4ebd2SPeter Brune               }
327ea5d4fccSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3288a430afbSPeter Brune             } else {
3298a430afbSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3308a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3318a430afbSPeter Brune               } else {
3328a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3338a430afbSPeter Brune               }
3348a430afbSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3358a430afbSPeter Brune             }
33671b4ebd2SPeter Brune           }
33771b4ebd2SPeter Brune           break;
338f5af7f23SKarl Rupp         } else if (monitor) {
339ea5d4fccSPeter Brune           ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
340bd4a8a71SBarry Smith           if (!objective) {
341b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
342c69d1a72SBarry 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);
34371b4ebd2SPeter Brune             } else {
344c69d1a72SBarry 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);
34571b4ebd2SPeter Brune             }
346ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3478a430afbSPeter Brune           } else {
3488a430afbSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3498a430afbSPeter 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);
3508a430afbSPeter Brune             } else {
3518a430afbSPeter 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);
3528a430afbSPeter Brune             }
3538a430afbSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3548a430afbSPeter Brune           }
35571b4ebd2SPeter Brune         }
35671b4ebd2SPeter Brune       }
35771b4ebd2SPeter Brune     }
35871b4ebd2SPeter Brune   }
35971b4ebd2SPeter Brune 
36071b4ebd2SPeter Brune   /* postcheck */
3617b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
36271b4ebd2SPeter Brune   if (changed_y) {
36371b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
3649bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
3659bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
3669bd66eb0SPeter Brune     }
36771b4ebd2SPeter Brune   }
368bd4a8a71SBarry Smith   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
369ce5a860aSPeter Brune     ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3709bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3719bd66eb0SPeter Brune       gnorm = fnorm;
3729bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3739bd66eb0SPeter Brune     } else {
3749bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3759bd66eb0SPeter Brune     }
3769bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
377c98378a5SBarry Smith     if (PetscIsInfOrNanReal(gnorm)) {
378422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
379c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
380c98378a5SBarry Smith       PetscFunctionReturn(0);
381c98378a5SBarry Smith     }
382c427506bSPeter Brune   }
38371b4ebd2SPeter Brune 
38471b4ebd2SPeter Brune   /* copy the solution over */
38571b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
386c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
387c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
388f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
389f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
39071b4ebd2SPeter Brune   PetscFunctionReturn(0);
39171b4ebd2SPeter Brune }
39271b4ebd2SPeter Brune 
3937f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
3947f1410a3SPeter Brune {
3957f1410a3SPeter Brune   PetscErrorCode    ierr;
3967f1410a3SPeter Brune   PetscBool         iascii;
3977f1410a3SPeter Brune   SNESLineSearch_BT *bt;
398075cc632SBarry Smith 
3997f1410a3SPeter Brune   PetscFunctionBegin;
400251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4017f1410a3SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
4027f1410a3SPeter Brune   if (iascii) {
403b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
4047f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
405b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
4067f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
4077f1410a3SPeter Brune     }
4087904a332SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha);CHKERRQ(ierr);
4097f1410a3SPeter Brune   }
4107f1410a3SPeter Brune   PetscFunctionReturn(0);
4117f1410a3SPeter Brune }
4127f1410a3SPeter Brune 
4137f1410a3SPeter Brune 
414f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
415731c067cSJed Brown {
416731c067cSJed Brown   PetscErrorCode ierr;
41771b4ebd2SPeter Brune 
41871b4ebd2SPeter Brune   PetscFunctionBegin;
419731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
42071b4ebd2SPeter Brune   PetscFunctionReturn(0);
42171b4ebd2SPeter Brune }
42271b4ebd2SPeter Brune 
42371b4ebd2SPeter Brune 
4244416b707SBarry Smith static PetscErrorCode SNESLineSearchSetFromOptions_BT(PetscOptionItems *PetscOptionsObject,SNESLineSearch linesearch)
42571b4ebd2SPeter Brune {
42671b4ebd2SPeter Brune 
42771b4ebd2SPeter Brune   PetscErrorCode    ierr;
428eb23ba39SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT*)linesearch->data;
42971b4ebd2SPeter Brune 
4306e111a19SKarl Rupp   PetscFunctionBegin;
431e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESLineSearch BT options");CHKERRQ(ierr);
4320298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr);
43371b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
43471b4ebd2SPeter Brune   PetscFunctionReturn(0);
43571b4ebd2SPeter Brune }
43671b4ebd2SPeter Brune 
43771b4ebd2SPeter Brune 
438954494b2SJed Brown /*MC
439db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
440954494b2SJed Brown 
441db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
442eb23ba39SBarry 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
443db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
444954494b2SJed Brown 
445cd7522eaSPeter Brune    Options Database Keys:
446cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
447db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
448eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - if the length the initial step is larger than this then the
449e1bc860dSBarry Smith                                        step is scaled back to be of this length at the beginning of the line search
450db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
451db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
452cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
453cd7522eaSPeter Brune 
454954494b2SJed Brown    Level: advanced
455954494b2SJed Brown 
456db609ea7SPeter Brune    Notes:
457db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
458db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
459db609ea7SPeter Brune 
460f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
461954494b2SJed Brown 
462f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
463954494b2SJed Brown M*/
4648cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
46571b4ebd2SPeter Brune {
46671b4ebd2SPeter Brune 
467f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
46871b4ebd2SPeter Brune   PetscErrorCode    ierr;
46971b4ebd2SPeter Brune 
47071b4ebd2SPeter Brune   PetscFunctionBegin;
471f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
472f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
473f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
4740298fd71SBarry Smith   linesearch->ops->reset          = NULL;
4757f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
4760298fd71SBarry Smith   linesearch->ops->setup          = NULL;
47771b4ebd2SPeter Brune 
478b00a9115SJed Brown   ierr = PetscNewLog(linesearch,&bt);CHKERRQ(ierr);
479f5af7f23SKarl Rupp 
48071b4ebd2SPeter Brune   linesearch->data    = (void*)bt;
48171b4ebd2SPeter Brune   linesearch->max_its = 40;
482b000cd8dSPeter Brune   linesearch->order   = SNES_LINESEARCH_ORDER_CUBIC;
48371b4ebd2SPeter Brune   bt->alpha           = 1e-4;
48471b4ebd2SPeter Brune   PetscFunctionReturn(0);
48571b4ebd2SPeter Brune }
486