xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision 4a2f8832e687075714ce232dc09a9bd05f920d30)
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 
871b4ebd2SPeter Brune #undef __FUNCT__
92f4102e2SPeter Brune #define __FUNCT__ "SNESLineSearchBTSetAlpha"
102f4102e2SPeter Brune /*@
112f4102e2SPeter Brune    SNESLineSearchBTSetAlpha - Sets the descent parameter, alpha, in the BT linesearch variant.
122f4102e2SPeter Brune 
132f4102e2SPeter Brune    Input Parameters:
142f4102e2SPeter Brune +  linesearch - linesearch context
152f4102e2SPeter Brune -  alpha - The descent parameter
162f4102e2SPeter Brune 
172f4102e2SPeter Brune    Level: intermediate
182f4102e2SPeter Brune 
191a4f838cSPeter Brune .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
202f4102e2SPeter Brune @*/
212f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha)
222f4102e2SPeter Brune {
232f4102e2SPeter Brune   SNESLineSearch_BT *bt;
246e111a19SKarl Rupp 
252f4102e2SPeter Brune   PetscFunctionBegin;
262f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
272f4102e2SPeter Brune   bt        = (SNESLineSearch_BT*)linesearch->data;
282f4102e2SPeter Brune   bt->alpha = alpha;
292f4102e2SPeter Brune   PetscFunctionReturn(0);
302f4102e2SPeter Brune }
312f4102e2SPeter Brune 
322f4102e2SPeter Brune 
332f4102e2SPeter Brune #undef __FUNCT__
342f4102e2SPeter Brune #define __FUNCT__ "SNESLineSearchBTGetAlpha"
352f4102e2SPeter Brune /*@
362f4102e2SPeter Brune    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant.
372f4102e2SPeter Brune 
382f4102e2SPeter Brune    Input Parameters:
392f4102e2SPeter Brune .  linesearch - linesearch context
408e557f58SPeter Brune 
418e557f58SPeter Brune    Output Parameters:
422f4102e2SPeter Brune .  alpha - The descent parameter
432f4102e2SPeter Brune 
442f4102e2SPeter Brune    Level: intermediate
452f4102e2SPeter Brune 
461a4f838cSPeter Brune .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
472f4102e2SPeter Brune @*/
482f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
492f4102e2SPeter Brune {
502f4102e2SPeter Brune   SNESLineSearch_BT *bt;
516e111a19SKarl Rupp 
522f4102e2SPeter Brune   PetscFunctionBegin;
532f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
542f4102e2SPeter Brune   bt     = (SNESLineSearch_BT*)linesearch->data;
552f4102e2SPeter Brune   *alpha = bt->alpha;
562f4102e2SPeter Brune   PetscFunctionReturn(0);
572f4102e2SPeter Brune }
582f4102e2SPeter Brune 
592f4102e2SPeter Brune #undef __FUNCT__
60f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_BT"
61f1c6b773SPeter Brune static PetscErrorCode  SNESLineSearchApply_BT(SNESLineSearch linesearch)
6271b4ebd2SPeter Brune {
6371b4ebd2SPeter Brune   PetscBool         changed_y,changed_w;
6471b4ebd2SPeter Brune   PetscErrorCode    ierr;
6571b4ebd2SPeter Brune   Vec               X,F,Y,W,G;
6671b4ebd2SPeter Brune   SNES              snes;
678a430afbSPeter Brune   PetscReal         fnorm, xnorm, ynorm, gnorm;
683b288129SPeter Brune   PetscReal         lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol;
69dfcd3828SPeter Brune   PetscReal         t1,t2,a,b,d;
708a430afbSPeter Brune   PetscReal         f;
718a430afbSPeter Brune   PetscReal         g,gprev;
7271b4ebd2SPeter Brune   PetscViewer       monitor;
7371b4ebd2SPeter Brune   PetscInt          max_its,count;
74f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
7571b4ebd2SPeter Brune   Mat               jac;
76bd4a8a71SBarry Smith   PetscErrorCode    (*objective)(SNES,Vec,PetscReal*,void*);
7771b4ebd2SPeter Brune 
7871b4ebd2SPeter Brune   PetscFunctionBegin;
79f1c6b773SPeter Brune   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
80f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
81f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
82f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
83dcf2fd19SBarry Smith   ierr = SNESLineSearchGetDefaultMonitor(linesearch, &monitor);CHKERRQ(ierr);
840298fd71SBarry Smith   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,NULL,NULL,NULL,&max_its);CHKERRQ(ierr);
850298fd71SBarry Smith   ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr);
860298fd71SBarry Smith   ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr);
87f1c6b773SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
8871b4ebd2SPeter Brune   alpha = bt->alpha;
8971b4ebd2SPeter Brune 
900298fd71SBarry Smith   ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr);
91ce94432eSBarry Smith   if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
92c69d1a72SBarry Smith 
937b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
94422a814eSBarry Smith   ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
9571b4ebd2SPeter Brune 
963b288129SPeter Brune   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
973b288129SPeter Brune   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
983b288129SPeter Brune   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
993b288129SPeter Brune   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
1003b288129SPeter Brune 
10171b4ebd2SPeter Brune   if (ynorm == 0.0) {
10271b4ebd2SPeter Brune     if (monitor) {
103ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10471b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
105ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10671b4ebd2SPeter Brune     }
10771b4ebd2SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
10871b4ebd2SPeter Brune     ierr = VecCopy(F,G);CHKERRQ(ierr);
10971cee68bSPeter Brune     ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
110e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
11171b4ebd2SPeter Brune     PetscFunctionReturn(0);
11271b4ebd2SPeter Brune   }
11371b4ebd2SPeter Brune   if (ynorm > maxstep) {        /* Step too big, so scale back */
11471b4ebd2SPeter Brune     if (monitor) {
115ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
116c69d1a72SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr);
117ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
11871b4ebd2SPeter Brune     }
11971b4ebd2SPeter Brune     ierr  = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
12071b4ebd2SPeter Brune     ynorm = maxstep;
12171b4ebd2SPeter Brune   }
1228a430afbSPeter Brune 
1238a430afbSPeter Brune   /* if the SNES has an objective set, use that instead of the function value */
124bd4a8a71SBarry Smith   if (objective) {
1258a430afbSPeter Brune     ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr);
1268a430afbSPeter Brune   } else {
1278a430afbSPeter Brune     f = fnorm*fnorm;
1288a430afbSPeter Brune   }
1298a430afbSPeter Brune 
1308a430afbSPeter Brune   /* compute the initial slope */
131bd4a8a71SBarry Smith   if (objective) {
1328a430afbSPeter Brune     /* slope comes from the function (assumed to be the gradient of the objective */
13367392de3SBarry Smith     ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr);
1348a430afbSPeter Brune   } else {
1358a430afbSPeter Brune     /* slope comes from the normal equations */
13671b4ebd2SPeter Brune     ierr = MatMult(jac,Y,W);CHKERRQ(ierr);
13767392de3SBarry Smith     ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr);
13871b4ebd2SPeter Brune     if (initslope > 0.0)  initslope = -initslope;
13971b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1408a430afbSPeter Brune   }
14171b4ebd2SPeter Brune 
142c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1439bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1449bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1459bd66eb0SPeter Brune   }
14671b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
14771b4ebd2SPeter Brune     ierr         = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
14871b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
149e9b602ebSSatish Balay     ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
15071b4ebd2SPeter Brune     PetscFunctionReturn(0);
15171b4ebd2SPeter Brune   }
1528a430afbSPeter Brune 
153bd4a8a71SBarry Smith   if (objective) {
1548a430afbSPeter Brune     ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
1558a430afbSPeter Brune   } else {
156ce5a860aSPeter Brune     ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
1579bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1589bd66eb0SPeter Brune       gnorm = fnorm;
1599bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1609bd66eb0SPeter Brune     } else {
16171b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1628a430afbSPeter Brune     }
16379e7e81bSJed Brown     g = PetscSqr(gnorm);
1649bd66eb0SPeter Brune   }
165dcf2fd19SBarry Smith   ierr = SNESLineSearchMonitor(linesearch);CHKERRQ(ierr);
1669bd66eb0SPeter Brune 
167c98378a5SBarry Smith   if (PetscIsInfOrNanReal(g)) {
168422a814eSBarry Smith     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
169c61ba1e2SPeter Brune     ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
170c98378a5SBarry Smith     PetscFunctionReturn(0);
171c98378a5SBarry Smith   }
172bd4a8a71SBarry Smith   if (!objective) {
173c69d1a72SBarry Smith     ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1748a430afbSPeter Brune   }
1758a430afbSPeter Brune   if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */
17671b4ebd2SPeter Brune     if (monitor) {
177ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
178bd4a8a71SBarry Smith       if (!objective) {
179c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1808a430afbSPeter Brune       } else {
1818a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr);
1828a430afbSPeter Brune       }
183ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
18471b4ebd2SPeter Brune     }
18571b4ebd2SPeter Brune   } else {
186c21ba15cSPeter Brune     /* Since the full step didn't work and the step is tiny, quit */
187c21ba15cSPeter Brune     if (stol*xnorm > ynorm) {
188ced04f9dSPeter Brune       ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
189e9b602ebSSatish Balay       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
190c61ba1e2SPeter Brune       if (monitor) {
191c61ba1e2SPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19208ed2907SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",(double)ynorm,(double)stol*xnorm);CHKERRQ(ierr);
193c61ba1e2SPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
194c61ba1e2SPeter Brune       }
195c21ba15cSPeter Brune       PetscFunctionReturn(0);
196c21ba15cSPeter Brune     }
19771b4ebd2SPeter Brune     /* Fit points with quadratic */
1988a430afbSPeter Brune     lambdatemp = -initslope/(g - f - 2.0*lambda*initslope);
19971b4ebd2SPeter Brune     lambdaprev = lambda;
2008a430afbSPeter Brune     gprev      = g;
20171b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
20271b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
20371b4ebd2SPeter Brune     else                         lambda = lambdatemp;
20471b4ebd2SPeter Brune 
20571b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2069bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2079bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2089bd66eb0SPeter Brune     }
20971b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
21071b4ebd2SPeter Brune       ierr         = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
21171b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
212e9b602ebSSatish Balay       ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
21371b4ebd2SPeter Brune       PetscFunctionReturn(0);
21471b4ebd2SPeter Brune     }
215bd4a8a71SBarry Smith     if (objective) {
2168a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2178a430afbSPeter Brune     } else {
218ce5a860aSPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
2199bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
2209bd66eb0SPeter Brune         gnorm = fnorm;
2219bd66eb0SPeter Brune         ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2229bd66eb0SPeter Brune       } else {
22371b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2249bd66eb0SPeter Brune       }
225f1ab69e3SPeter Brune       g = PetscSqr(gnorm);
2268a430afbSPeter Brune     }
227c98378a5SBarry Smith     if (PetscIsInfOrNanReal(g)) {
228422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
229c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
230c98378a5SBarry Smith       PetscFunctionReturn(0);
231c98378a5SBarry Smith     }
23271b4ebd2SPeter Brune     if (monitor) {
233ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
234bd4a8a71SBarry Smith       if (!objective) {
235c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr);
2368a430afbSPeter Brune       } else {
2378a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr);
2388a430afbSPeter Brune       }
239ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24071b4ebd2SPeter Brune     }
2418a430afbSPeter Brune     if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */
24271b4ebd2SPeter Brune       if (monitor) {
243ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24471b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
245ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24671b4ebd2SPeter Brune       }
24771b4ebd2SPeter Brune     } else {
24871b4ebd2SPeter Brune       /* Fit points with cubic */
24971b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
25071b4ebd2SPeter Brune         if (lambda <= minlambda) {
25171b4ebd2SPeter Brune           if (monitor) {
252ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
25371b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
254bd4a8a71SBarry Smith             if (!objective) {
255*4a2f8832SBarry 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",
256c69d1a72SBarry Smith                                                          (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2578a430afbSPeter Brune             } else {
258*4a2f8832SBarry 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",
2598a430afbSPeter Brune                                                          (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2608a430afbSPeter Brune             }
261ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26271b4ebd2SPeter Brune           }
263e9b602ebSSatish Balay           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
26471b4ebd2SPeter Brune           PetscFunctionReturn(0);
26571b4ebd2SPeter Brune         }
266b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
2678a430afbSPeter Brune           t1 = .5*(g - f) - lambda*initslope;
2688a430afbSPeter Brune           t2 = .5*(gprev  - f) - lambdaprev*initslope;
26971b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27071b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27171b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
27271b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
273f5af7f23SKarl Rupp           if (a == 0.0) lambdatemp = -initslope/(2.0*b);
274f5af7f23SKarl Rupp           else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
275f5af7f23SKarl Rupp 
276b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
2778a430afbSPeter Brune           lambdatemp = -initslope/(g - f - 2.0*initslope);
278ce94432eSBarry Smith         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt");
27971b4ebd2SPeter Brune         lambdaprev = lambda;
2808a430afbSPeter Brune         gprev      = g;
28171b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
28271b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
28371b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
28471b4ebd2SPeter Brune         ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
285b7b2e573SPeter Brune         if (linesearch->ops->viproject) {
286b7b2e573SPeter Brune           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
287b7b2e573SPeter Brune         }
28871b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
28971b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
290bd4a8a71SBarry Smith           if (!objective) {
29171b4ebd2SPeter Brune             ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
292c69d1a72SBarry Smith                               (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
2938a430afbSPeter Brune           }
294e9b602ebSSatish Balay           ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
29571b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
29671b4ebd2SPeter Brune           PetscFunctionReturn(0);
29771b4ebd2SPeter Brune         }
298bd4a8a71SBarry Smith         if (objective) {
2998a430afbSPeter Brune           ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
3008a430afbSPeter Brune         } else {
301ce5a860aSPeter Brune           ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3029bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
3039bd66eb0SPeter Brune             gnorm = fnorm;
3049bd66eb0SPeter Brune             ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3059bd66eb0SPeter Brune           } else {
30671b4ebd2SPeter Brune             ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3078a430afbSPeter Brune           }
308f1ab69e3SPeter Brune           g = PetscSqr(gnorm);
3099bd66eb0SPeter Brune         }
310*4a2f8832SBarry Smith         if (PetscIsInfOrNanReal(g)) {
311422a814eSBarry Smith           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
312c61ba1e2SPeter Brune           ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
313c98378a5SBarry Smith           PetscFunctionReturn(0);
314c98378a5SBarry Smith         }
3158a430afbSPeter Brune         if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */
31671b4ebd2SPeter Brune           if (monitor) {
317ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
318bd4a8a71SBarry Smith             if (!objective) {
319b000cd8dSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
320c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32171b4ebd2SPeter Brune               } else {
322c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32371b4ebd2SPeter Brune               }
324ea5d4fccSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3258a430afbSPeter Brune             } else {
3268a430afbSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3278a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3288a430afbSPeter Brune               } else {
3298a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3308a430afbSPeter Brune               }
3318a430afbSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3328a430afbSPeter Brune             }
33371b4ebd2SPeter Brune           }
33471b4ebd2SPeter Brune           break;
335f5af7f23SKarl Rupp         } else if (monitor) {
336ea5d4fccSPeter Brune           ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
337bd4a8a71SBarry Smith           if (!objective) {
338b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
339c69d1a72SBarry 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);
34071b4ebd2SPeter Brune             } else {
341c69d1a72SBarry 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);
34271b4ebd2SPeter Brune             }
343ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3448a430afbSPeter Brune           } else {
3458a430afbSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3468a430afbSPeter 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);
3478a430afbSPeter Brune             } else {
3488a430afbSPeter 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);
3498a430afbSPeter Brune             }
3508a430afbSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3518a430afbSPeter Brune           }
35271b4ebd2SPeter Brune         }
35371b4ebd2SPeter Brune       }
35471b4ebd2SPeter Brune     }
35571b4ebd2SPeter Brune   }
35671b4ebd2SPeter Brune 
35771b4ebd2SPeter Brune   /* postcheck */
3587b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
35971b4ebd2SPeter Brune   if (changed_y) {
36071b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,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 
39071b4ebd2SPeter Brune #undef __FUNCT__
3917f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT"
3927f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
3937f1410a3SPeter Brune {
3947f1410a3SPeter Brune   PetscErrorCode    ierr;
3957f1410a3SPeter Brune   PetscBool         iascii;
3967f1410a3SPeter Brune   SNESLineSearch_BT *bt;
397075cc632SBarry Smith 
3987f1410a3SPeter Brune   PetscFunctionBegin;
399251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4007f1410a3SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
4017f1410a3SPeter Brune   if (iascii) {
402b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
4037f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
404b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
4057f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
4067f1410a3SPeter Brune     }
4077904a332SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha);CHKERRQ(ierr);
4087f1410a3SPeter Brune   }
4097f1410a3SPeter Brune   PetscFunctionReturn(0);
4107f1410a3SPeter Brune }
4117f1410a3SPeter Brune 
4127f1410a3SPeter Brune 
4137f1410a3SPeter Brune #undef __FUNCT__
414f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT"
415f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
416731c067cSJed Brown {
417731c067cSJed Brown   PetscErrorCode ierr;
41871b4ebd2SPeter Brune 
41971b4ebd2SPeter Brune   PetscFunctionBegin;
420731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
42171b4ebd2SPeter Brune   PetscFunctionReturn(0);
42271b4ebd2SPeter Brune }
42371b4ebd2SPeter Brune 
42471b4ebd2SPeter Brune 
42571b4ebd2SPeter Brune #undef __FUNCT__
426f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
4274416b707SBarry Smith static PetscErrorCode SNESLineSearchSetFromOptions_BT(PetscOptionItems *PetscOptionsObject,SNESLineSearch linesearch)
42871b4ebd2SPeter Brune {
42971b4ebd2SPeter Brune 
43071b4ebd2SPeter Brune   PetscErrorCode    ierr;
431eb23ba39SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT*)linesearch->data;
43271b4ebd2SPeter Brune 
4336e111a19SKarl Rupp   PetscFunctionBegin;
434e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESLineSearch BT options");CHKERRQ(ierr);
4350298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr);
43671b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
43771b4ebd2SPeter Brune   PetscFunctionReturn(0);
43871b4ebd2SPeter Brune }
43971b4ebd2SPeter Brune 
44071b4ebd2SPeter Brune 
44171b4ebd2SPeter Brune #undef __FUNCT__
442f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT"
443954494b2SJed Brown /*MC
444db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
445954494b2SJed Brown 
446db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
447eb23ba39SBarry 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
448db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
449954494b2SJed Brown 
450cd7522eaSPeter Brune    Options Database Keys:
451cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
452db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
453eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - if the length the initial step is larger than this then the
454e1bc860dSBarry Smith                                        step is scaled back to be of this length at the beginning of the line search
455db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
456db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
457cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
458cd7522eaSPeter Brune 
459954494b2SJed Brown    Level: advanced
460954494b2SJed Brown 
461db609ea7SPeter Brune    Notes:
462db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
463db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
464db609ea7SPeter Brune 
465f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
466954494b2SJed Brown 
467f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
468954494b2SJed Brown M*/
4698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
47071b4ebd2SPeter Brune {
47171b4ebd2SPeter Brune 
472f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
47371b4ebd2SPeter Brune   PetscErrorCode    ierr;
47471b4ebd2SPeter Brune 
47571b4ebd2SPeter Brune   PetscFunctionBegin;
476f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
477f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
478f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
4790298fd71SBarry Smith   linesearch->ops->reset          = NULL;
4807f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
4810298fd71SBarry Smith   linesearch->ops->setup          = NULL;
48271b4ebd2SPeter Brune 
483b00a9115SJed Brown   ierr = PetscNewLog(linesearch,&bt);CHKERRQ(ierr);
484f5af7f23SKarl Rupp 
48571b4ebd2SPeter Brune   linesearch->data    = (void*)bt;
48671b4ebd2SPeter Brune   linesearch->max_its = 40;
487b000cd8dSPeter Brune   linesearch->order   = SNES_LINESEARCH_ORDER_CUBIC;
48871b4ebd2SPeter Brune   bt->alpha           = 1e-4;
48971b4ebd2SPeter Brune   PetscFunctionReturn(0);
49071b4ebd2SPeter Brune }
491