xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision e9b602ebd060280b77de97e0e037db14ab93a4a7)
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);
83f1c6b773SPeter Brune   ierr = SNESLineSearchGetMonitor(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 
8971b4ebd2SPeter Brune   alpha = bt->alpha;
9071b4ebd2SPeter Brune 
910298fd71SBarry Smith   ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr);
928a430afbSPeter Brune 
93ce94432eSBarry Smith   if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
94c69d1a72SBarry Smith 
9571b4ebd2SPeter Brune   /* precheck */
967b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
97422a814eSBarry Smith   ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
9871b4ebd2SPeter Brune 
993b288129SPeter Brune   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
1003b288129SPeter Brune   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
1013b288129SPeter Brune   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
1023b288129SPeter Brune   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
1033b288129SPeter Brune 
10471b4ebd2SPeter Brune   if (ynorm == 0.0) {
10571b4ebd2SPeter Brune     if (monitor) {
106ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10771b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
108ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10971b4ebd2SPeter Brune     }
11071b4ebd2SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
11171b4ebd2SPeter Brune     ierr = VecCopy(F,G);CHKERRQ(ierr);
11271cee68bSPeter Brune     ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
113*e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
11471b4ebd2SPeter Brune     PetscFunctionReturn(0);
11571b4ebd2SPeter Brune   }
11671b4ebd2SPeter Brune   if (ynorm > maxstep) {        /* Step too big, so scale back */
11771b4ebd2SPeter Brune     if (monitor) {
118ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
119c69d1a72SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr);
120ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
12171b4ebd2SPeter Brune     }
12271b4ebd2SPeter Brune     ierr  = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
12371b4ebd2SPeter Brune     ynorm = maxstep;
12471b4ebd2SPeter Brune   }
1258a430afbSPeter Brune 
1268a430afbSPeter Brune   /* if the SNES has an objective set, use that instead of the function value */
127bd4a8a71SBarry Smith   if (objective) {
1288a430afbSPeter Brune     ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr);
1298a430afbSPeter Brune   } else {
1308a430afbSPeter Brune     f = fnorm*fnorm;
1318a430afbSPeter Brune   }
1328a430afbSPeter Brune 
1338a430afbSPeter Brune   /* compute the initial slope */
134bd4a8a71SBarry Smith   if (objective) {
1358a430afbSPeter Brune     /* slope comes from the function (assumed to be the gradient of the objective */
13667392de3SBarry Smith     ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr);
1378a430afbSPeter Brune   } else {
1388a430afbSPeter Brune     /* slope comes from the normal equations */
13971b4ebd2SPeter Brune     ierr = MatMult(jac,Y,W);CHKERRQ(ierr);
14067392de3SBarry Smith     ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr);
14171b4ebd2SPeter Brune     if (initslope > 0.0)  initslope = -initslope;
14271b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1438a430afbSPeter Brune   }
14471b4ebd2SPeter Brune 
145c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1469bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1479bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1489bd66eb0SPeter Brune   }
14971b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
15071b4ebd2SPeter Brune     ierr         = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
15171b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
152*e9b602ebSSatish Balay     ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
15371b4ebd2SPeter Brune     PetscFunctionReturn(0);
15471b4ebd2SPeter Brune   }
1558a430afbSPeter Brune 
156bd4a8a71SBarry Smith   if (objective) {
1578a430afbSPeter Brune     ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
1588a430afbSPeter Brune   } else {
159ce5a860aSPeter Brune     ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
1609bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1619bd66eb0SPeter Brune       gnorm = fnorm;
1629bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1639bd66eb0SPeter Brune     } else {
16471b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1658a430afbSPeter Brune     }
16679e7e81bSJed Brown     g = PetscSqr(gnorm);
1679bd66eb0SPeter Brune   }
1689bd66eb0SPeter Brune 
169c98378a5SBarry Smith   if (PetscIsInfOrNanReal(g)) {
170422a814eSBarry Smith     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
171c61ba1e2SPeter Brune     ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
172c98378a5SBarry Smith     PetscFunctionReturn(0);
173c98378a5SBarry Smith   }
174bd4a8a71SBarry Smith   if (!objective) {
175c69d1a72SBarry Smith     ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1768a430afbSPeter Brune   }
1778a430afbSPeter Brune   if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */
17871b4ebd2SPeter Brune     if (monitor) {
179ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
180bd4a8a71SBarry Smith       if (!objective) {
181c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1828a430afbSPeter Brune       } else {
1838a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr);
1848a430afbSPeter Brune       }
185ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
18671b4ebd2SPeter Brune     }
18771b4ebd2SPeter Brune   } else {
188c21ba15cSPeter Brune     /* Since the full step didn't work and the step is tiny, quit */
189c21ba15cSPeter Brune     if (stol*xnorm > ynorm) {
190ced04f9dSPeter Brune       ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
191*e9b602ebSSatish Balay       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
192c61ba1e2SPeter Brune       if (monitor) {
193c61ba1e2SPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19408ed2907SPeter 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);
195c61ba1e2SPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
196c61ba1e2SPeter Brune       }
197c21ba15cSPeter Brune       PetscFunctionReturn(0);
198c21ba15cSPeter Brune     }
19971b4ebd2SPeter Brune     /* Fit points with quadratic */
2008a430afbSPeter Brune     lambdatemp = -initslope/(g - f - 2.0*lambda*initslope);
20171b4ebd2SPeter Brune     lambdaprev = lambda;
2028a430afbSPeter Brune     gprev      = g;
20371b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
20471b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
20571b4ebd2SPeter Brune     else                         lambda = lambdatemp;
20671b4ebd2SPeter Brune 
20771b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2089bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2099bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2109bd66eb0SPeter Brune     }
21171b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
21271b4ebd2SPeter Brune       ierr         = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
21371b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
214*e9b602ebSSatish Balay       ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
21571b4ebd2SPeter Brune       PetscFunctionReturn(0);
21671b4ebd2SPeter Brune     }
217bd4a8a71SBarry Smith     if (objective) {
2188a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2198a430afbSPeter Brune     } else {
220ce5a860aSPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
2219bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
2229bd66eb0SPeter Brune         gnorm = fnorm;
2239bd66eb0SPeter Brune         ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2249bd66eb0SPeter Brune       } else {
22571b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2269bd66eb0SPeter Brune       }
227f1ab69e3SPeter Brune       g = PetscSqr(gnorm);
2288a430afbSPeter Brune     }
229c98378a5SBarry Smith     if (PetscIsInfOrNanReal(g)) {
230422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
231c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
232c98378a5SBarry Smith       PetscFunctionReturn(0);
233c98378a5SBarry Smith     }
23471b4ebd2SPeter Brune     if (monitor) {
235ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
236bd4a8a71SBarry Smith       if (!objective) {
237c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr);
2388a430afbSPeter Brune       } else {
2398a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr);
2408a430afbSPeter Brune       }
241ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24271b4ebd2SPeter Brune     }
2438a430afbSPeter Brune     if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */
24471b4ebd2SPeter Brune       if (monitor) {
245ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24671b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
247ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24871b4ebd2SPeter Brune       }
24971b4ebd2SPeter Brune     } else {
25071b4ebd2SPeter Brune       /* Fit points with cubic */
25171b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
25271b4ebd2SPeter Brune         if (lambda <= minlambda) {
25371b4ebd2SPeter Brune           if (monitor) {
254ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
25571b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
256bd4a8a71SBarry Smith             if (!objective) {
25771b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,
25871b4ebd2SPeter Brune                                             "    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 {
2618a430afbSPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,
2628a430afbSPeter Brune                                             "    Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
2638a430afbSPeter Brune                                             (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2648a430afbSPeter Brune             }
265ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26671b4ebd2SPeter Brune           }
267*e9b602ebSSatish Balay           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
26871b4ebd2SPeter Brune           PetscFunctionReturn(0);
26971b4ebd2SPeter Brune         }
270b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
2718a430afbSPeter Brune           t1 = .5*(g - f) - lambda*initslope;
2728a430afbSPeter Brune           t2 = .5*(gprev  - f) - lambdaprev*initslope;
27371b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27471b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27571b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
27671b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
277f5af7f23SKarl Rupp           if (a == 0.0) lambdatemp = -initslope/(2.0*b);
278f5af7f23SKarl Rupp           else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
279f5af7f23SKarl Rupp 
280b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
2818a430afbSPeter Brune           lambdatemp = -initslope/(g - f - 2.0*initslope);
282ce94432eSBarry Smith         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt");
28371b4ebd2SPeter Brune         lambdaprev = lambda;
2848a430afbSPeter Brune         gprev      = g;
28571b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
28671b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
28771b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
28871b4ebd2SPeter Brune         ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
289b7b2e573SPeter Brune         if (linesearch->ops->viproject) {
290b7b2e573SPeter Brune           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
291b7b2e573SPeter Brune         }
29271b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
29371b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
294bd4a8a71SBarry Smith           if (!objective) {
29571b4ebd2SPeter Brune             ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
296c69d1a72SBarry Smith                               (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
2978a430afbSPeter Brune           }
298*e9b602ebSSatish Balay           ierr         = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION);CHKERRQ(ierr);
29971b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
30071b4ebd2SPeter Brune           PetscFunctionReturn(0);
30171b4ebd2SPeter Brune         }
302bd4a8a71SBarry Smith         if (objective) {
3038a430afbSPeter Brune           ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
3048a430afbSPeter Brune         } else {
305ce5a860aSPeter Brune           ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3069bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
3079bd66eb0SPeter Brune             gnorm = fnorm;
3089bd66eb0SPeter Brune             ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3099bd66eb0SPeter Brune           } else {
31071b4ebd2SPeter Brune             ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3118a430afbSPeter Brune           }
312f1ab69e3SPeter Brune           g = PetscSqr(gnorm);
3139bd66eb0SPeter Brune         }
314c98378a5SBarry Smith         if (PetscIsInfOrNanReal(gnorm)) {
315422a814eSBarry Smith           ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
316c61ba1e2SPeter Brune           ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
317c98378a5SBarry Smith           PetscFunctionReturn(0);
318c98378a5SBarry Smith         }
3198a430afbSPeter Brune         if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */
32071b4ebd2SPeter Brune           if (monitor) {
321ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
322bd4a8a71SBarry Smith             if (!objective) {
323b000cd8dSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
324c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32571b4ebd2SPeter Brune               } else {
326c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32771b4ebd2SPeter Brune               }
328ea5d4fccSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3298a430afbSPeter Brune             } else {
3308a430afbSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3318a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3328a430afbSPeter Brune               } else {
3338a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3348a430afbSPeter Brune               }
3358a430afbSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3368a430afbSPeter Brune             }
33771b4ebd2SPeter Brune           }
33871b4ebd2SPeter Brune           break;
339f5af7f23SKarl Rupp         } else if (monitor) {
340ea5d4fccSPeter Brune           ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
341bd4a8a71SBarry Smith           if (!objective) {
342b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
343c69d1a72SBarry 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);
34471b4ebd2SPeter Brune             } else {
345c69d1a72SBarry 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);
34671b4ebd2SPeter Brune             }
347ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3488a430afbSPeter Brune           } else {
3498a430afbSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3508a430afbSPeter 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);
3518a430afbSPeter Brune             } else {
3528a430afbSPeter 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);
3538a430afbSPeter Brune             }
3548a430afbSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3558a430afbSPeter Brune           }
35671b4ebd2SPeter Brune         }
35771b4ebd2SPeter Brune       }
35871b4ebd2SPeter Brune     }
35971b4ebd2SPeter Brune   }
36071b4ebd2SPeter Brune 
36171b4ebd2SPeter Brune   /* postcheck */
3627b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
36371b4ebd2SPeter Brune   if (changed_y) {
36471b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
3659bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
3669bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
3679bd66eb0SPeter Brune     }
36871b4ebd2SPeter Brune   }
369bd4a8a71SBarry Smith   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
370ce5a860aSPeter Brune     ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr);
3719bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3729bd66eb0SPeter Brune       gnorm = fnorm;
3739bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3749bd66eb0SPeter Brune     } else {
3759bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3769bd66eb0SPeter Brune     }
3779bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
378c98378a5SBarry Smith     if (PetscIsInfOrNanReal(gnorm)) {
379422a814eSBarry Smith       ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_FAILED_NANORINF);CHKERRQ(ierr);
380c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
381c98378a5SBarry Smith       PetscFunctionReturn(0);
382c98378a5SBarry Smith     }
383c427506bSPeter Brune   }
38471b4ebd2SPeter Brune 
38571b4ebd2SPeter Brune   /* copy the solution over */
38671b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
387c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
388c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
389f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
390f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
39171b4ebd2SPeter Brune   PetscFunctionReturn(0);
39271b4ebd2SPeter Brune }
39371b4ebd2SPeter Brune 
39471b4ebd2SPeter Brune #undef __FUNCT__
3957f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT"
3967f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
3977f1410a3SPeter Brune {
3987f1410a3SPeter Brune   PetscErrorCode    ierr;
3997f1410a3SPeter Brune   PetscBool         iascii;
4007f1410a3SPeter Brune   SNESLineSearch_BT *bt;
401075cc632SBarry Smith 
4027f1410a3SPeter Brune   PetscFunctionBegin;
403251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4047f1410a3SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
4057f1410a3SPeter Brune   if (iascii) {
406b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
4077f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
408b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
4097f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
4107f1410a3SPeter Brune     }
4117904a332SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", (double)bt->alpha);CHKERRQ(ierr);
4127f1410a3SPeter Brune   }
4137f1410a3SPeter Brune   PetscFunctionReturn(0);
4147f1410a3SPeter Brune }
4157f1410a3SPeter Brune 
4167f1410a3SPeter Brune 
4177f1410a3SPeter Brune #undef __FUNCT__
418f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT"
419f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
420731c067cSJed Brown {
421731c067cSJed Brown   PetscErrorCode ierr;
42271b4ebd2SPeter Brune 
42371b4ebd2SPeter Brune   PetscFunctionBegin;
424731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
42571b4ebd2SPeter Brune   PetscFunctionReturn(0);
42671b4ebd2SPeter Brune }
42771b4ebd2SPeter Brune 
42871b4ebd2SPeter Brune 
42971b4ebd2SPeter Brune #undef __FUNCT__
430f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
4318c34d3f5SBarry Smith static PetscErrorCode SNESLineSearchSetFromOptions_BT(PetscOptions *PetscOptionsObject,SNESLineSearch linesearch)
43271b4ebd2SPeter Brune {
43371b4ebd2SPeter Brune 
43471b4ebd2SPeter Brune   PetscErrorCode    ierr;
435eb23ba39SBarry Smith   SNESLineSearch_BT *bt = (SNESLineSearch_BT*)linesearch->data;
43671b4ebd2SPeter Brune 
4376e111a19SKarl Rupp   PetscFunctionBegin;
438e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESLineSearch BT options");CHKERRQ(ierr);
4390298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr);
44071b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
44171b4ebd2SPeter Brune   PetscFunctionReturn(0);
44271b4ebd2SPeter Brune }
44371b4ebd2SPeter Brune 
44471b4ebd2SPeter Brune 
44571b4ebd2SPeter Brune #undef __FUNCT__
446f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT"
447954494b2SJed Brown /*MC
448db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
449954494b2SJed Brown 
450db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
451eb23ba39SBarry 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
452db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
453954494b2SJed Brown 
454cd7522eaSPeter Brune    Options Database Keys:
455cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
456db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
457eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - if the length the initial step is larger than this then the
458e1bc860dSBarry Smith                                        step is scaled back to be of this length at the beginning of the line search
459db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
460db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
461cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
462cd7522eaSPeter Brune 
463954494b2SJed Brown    Level: advanced
464954494b2SJed Brown 
465db609ea7SPeter Brune    Notes:
466db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
467db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
468db609ea7SPeter Brune 
469f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
470954494b2SJed Brown 
471f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
472954494b2SJed Brown M*/
4738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
47471b4ebd2SPeter Brune {
47571b4ebd2SPeter Brune 
476f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
47771b4ebd2SPeter Brune   PetscErrorCode    ierr;
47871b4ebd2SPeter Brune 
47971b4ebd2SPeter Brune   PetscFunctionBegin;
480f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
481f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
482f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
4830298fd71SBarry Smith   linesearch->ops->reset          = NULL;
4847f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
4850298fd71SBarry Smith   linesearch->ops->setup          = NULL;
48671b4ebd2SPeter Brune 
487b00a9115SJed Brown   ierr = PetscNewLog(linesearch,&bt);CHKERRQ(ierr);
488f5af7f23SKarl Rupp 
48971b4ebd2SPeter Brune   linesearch->data    = (void*)bt;
49071b4ebd2SPeter Brune   linesearch->max_its = 40;
491b000cd8dSPeter Brune   linesearch->order   = SNES_LINESEARCH_ORDER_CUBIC;
49271b4ebd2SPeter Brune   bt->alpha           = 1e-4;
49371b4ebd2SPeter Brune   PetscFunctionReturn(0);
49471b4ebd2SPeter Brune }
495