xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision c61ba1e2242b7baea15736f193d31bad7aa27007)
13d443117SPeter Brune #include <petsc-private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
2b45d2f2cSJed Brown #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   PetscBool         domainerror;
7371b4ebd2SPeter Brune   PetscViewer       monitor;
7471b4ebd2SPeter Brune   PetscInt          max_its,count;
75f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
7671b4ebd2SPeter Brune   Mat               jac;
77bd4a8a71SBarry Smith   PetscErrorCode    (*objective)(SNES,Vec,PetscReal*,void*);
7871b4ebd2SPeter Brune 
7971b4ebd2SPeter Brune   PetscFunctionBegin;
80f1c6b773SPeter Brune   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
81f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
82f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
83f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
84f1c6b773SPeter Brune   ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr);
850298fd71SBarry Smith   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,NULL,NULL,NULL,&max_its);CHKERRQ(ierr);
860298fd71SBarry Smith   ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr);
870298fd71SBarry Smith   ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr);
88f1c6b773SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
8971b4ebd2SPeter Brune 
9071b4ebd2SPeter Brune   alpha = bt->alpha;
9171b4ebd2SPeter Brune 
920298fd71SBarry Smith   ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr);
938a430afbSPeter Brune 
94ce94432eSBarry Smith   if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
95c69d1a72SBarry Smith 
9671b4ebd2SPeter Brune   /* precheck */
977b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
98f1c6b773SPeter Brune   ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
9971b4ebd2SPeter Brune 
1003b288129SPeter Brune   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
1013b288129SPeter Brune   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
1023b288129SPeter Brune   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
1033b288129SPeter Brune   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
1043b288129SPeter Brune 
10571b4ebd2SPeter Brune   if (ynorm == 0.0) {
10671b4ebd2SPeter Brune     if (monitor) {
107ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10871b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
109ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
11071b4ebd2SPeter Brune     }
11171b4ebd2SPeter Brune     ierr = VecCopy(X,W);CHKERRQ(ierr);
11271b4ebd2SPeter Brune     ierr = VecCopy(F,G);CHKERRQ(ierr);
11371cee68bSPeter Brune     ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
114f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
11571b4ebd2SPeter Brune     PetscFunctionReturn(0);
11671b4ebd2SPeter Brune   }
11771b4ebd2SPeter Brune   if (ynorm > maxstep) {        /* Step too big, so scale back */
11871b4ebd2SPeter Brune     if (monitor) {
119ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
120c69d1a72SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr);
121ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
12271b4ebd2SPeter Brune     }
12371b4ebd2SPeter Brune     ierr  = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
12471b4ebd2SPeter Brune     ynorm = maxstep;
12571b4ebd2SPeter Brune   }
1268a430afbSPeter Brune 
1278a430afbSPeter Brune   /* if the SNES has an objective set, use that instead of the function value */
128bd4a8a71SBarry Smith   if (objective) {
1298a430afbSPeter Brune     ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr);
1308a430afbSPeter Brune   } else {
1318a430afbSPeter Brune     f = fnorm*fnorm;
1328a430afbSPeter Brune   }
1338a430afbSPeter Brune 
1348a430afbSPeter Brune   /* compute the initial slope */
135bd4a8a71SBarry Smith   if (objective) {
1368a430afbSPeter Brune     /* slope comes from the function (assumed to be the gradient of the objective */
13767392de3SBarry Smith     ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr);
1388a430afbSPeter Brune   } else {
1398a430afbSPeter Brune     /* slope comes from the normal equations */
14071b4ebd2SPeter Brune     ierr = MatMult(jac,Y,W);CHKERRQ(ierr);
14167392de3SBarry Smith     ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr);
14271b4ebd2SPeter Brune     if (initslope > 0.0)  initslope = -initslope;
14371b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1448a430afbSPeter Brune   }
14571b4ebd2SPeter Brune 
146c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1479bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1489bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1499bd66eb0SPeter Brune   }
15071b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
15171b4ebd2SPeter Brune     ierr         = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
15271b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
153f1c6b773SPeter Brune     ierr         = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
15471b4ebd2SPeter Brune     PetscFunctionReturn(0);
15571b4ebd2SPeter Brune   }
1568a430afbSPeter Brune 
157bd4a8a71SBarry Smith   if (objective) {
1588a430afbSPeter Brune     ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
1598a430afbSPeter Brune   } else {
16071b4ebd2SPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
161c427506bSPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
162c427506bSPeter Brune     if (domainerror) {
163f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
16471b4ebd2SPeter Brune       PetscFunctionReturn(0);
16571b4ebd2SPeter Brune     }
1669bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1679bd66eb0SPeter Brune       gnorm = fnorm;
1689bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1699bd66eb0SPeter Brune     } else {
17071b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1718a430afbSPeter Brune     }
17279e7e81bSJed Brown     g = PetscSqr(gnorm);
1739bd66eb0SPeter Brune   }
1749bd66eb0SPeter Brune 
175c98378a5SBarry Smith   if (PetscIsInfOrNanReal(g)) {
176c98378a5SBarry Smith     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
177*c61ba1e2SPeter Brune     ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
178c98378a5SBarry Smith     PetscFunctionReturn(0);
179c98378a5SBarry Smith   }
180bd4a8a71SBarry Smith   if (!objective) {
181*c61ba1e2SPeter Brune     if (monitor) {
182*c61ba1e2SPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
183*c61ba1e2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
184*c61ba1e2SPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
185*c61ba1e2SPeter Brune     }
1868a430afbSPeter Brune   }
1878a430afbSPeter Brune   if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */
18871b4ebd2SPeter Brune     if (monitor) {
189ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
190bd4a8a71SBarry Smith       if (!objective) {
191c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1928a430afbSPeter Brune       } else {
1938a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr);
1948a430afbSPeter Brune       }
195ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19671b4ebd2SPeter Brune     }
19771b4ebd2SPeter Brune   } else {
198c21ba15cSPeter Brune     /* Since the full step didn't work and the step is tiny, quit */
199c21ba15cSPeter Brune     if (stol*xnorm > ynorm) {
200ced04f9dSPeter Brune       ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
201c21ba15cSPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
202*c61ba1e2SPeter Brune       if (monitor) {
203*c61ba1e2SPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
204*c61ba1e2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr);
205*c61ba1e2SPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
206*c61ba1e2SPeter Brune       }
207c21ba15cSPeter Brune       PetscFunctionReturn(0);
208c21ba15cSPeter Brune     }
20971b4ebd2SPeter Brune     /* Fit points with quadratic */
2108a430afbSPeter Brune     lambdatemp = -initslope/(g - f - 2.0*lambda*initslope);
21171b4ebd2SPeter Brune     lambdaprev = lambda;
2128a430afbSPeter Brune     gprev      = g;
21371b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
21471b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
21571b4ebd2SPeter Brune     else                         lambda = lambdatemp;
21671b4ebd2SPeter Brune 
21771b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2189bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2199bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2209bd66eb0SPeter Brune     }
22171b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
22271b4ebd2SPeter Brune       ierr         = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
22371b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
224f1c6b773SPeter Brune       ierr         = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
22571b4ebd2SPeter Brune       PetscFunctionReturn(0);
22671b4ebd2SPeter Brune     }
227bd4a8a71SBarry Smith     if (objective) {
2288a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2298a430afbSPeter Brune     } else {
23071b4ebd2SPeter Brune       ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
2316188f407SPeter Brune       ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
232f5af7f23SKarl Rupp       if (domainerror) PetscFunctionReturn(0);
233f5af7f23SKarl Rupp 
2349bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
2359bd66eb0SPeter Brune         gnorm = fnorm;
2369bd66eb0SPeter Brune         ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2379bd66eb0SPeter Brune       } else {
23871b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2399bd66eb0SPeter Brune       }
240f1ab69e3SPeter Brune       g = PetscSqr(gnorm);
2418a430afbSPeter Brune     }
242c98378a5SBarry Smith     if (PetscIsInfOrNanReal(g)) {
243c98378a5SBarry Smith       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
244*c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
245c98378a5SBarry Smith       PetscFunctionReturn(0);
246c98378a5SBarry Smith     }
24771b4ebd2SPeter Brune     if (monitor) {
248ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
249bd4a8a71SBarry Smith       if (!objective) {
250c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr);
2518a430afbSPeter Brune       } else {
2528a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr);
2538a430afbSPeter Brune       }
254ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
25571b4ebd2SPeter Brune     }
2568a430afbSPeter Brune     if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */
25771b4ebd2SPeter Brune       if (monitor) {
258ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
25971b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
260ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26171b4ebd2SPeter Brune       }
26271b4ebd2SPeter Brune     } else {
26371b4ebd2SPeter Brune       /* Fit points with cubic */
26471b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
26571b4ebd2SPeter Brune         if (lambda <= minlambda) {
26671b4ebd2SPeter Brune           if (monitor) {
267ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26871b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
269bd4a8a71SBarry Smith             if (!objective) {
27071b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,
27171b4ebd2SPeter Brune                                             "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
272c69d1a72SBarry Smith                                             (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2738a430afbSPeter Brune             } else {
2748a430afbSPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,
2758a430afbSPeter Brune                                             "    Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
2768a430afbSPeter Brune                                             (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2778a430afbSPeter Brune             }
278ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
27971b4ebd2SPeter Brune           }
280f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
28171b4ebd2SPeter Brune           PetscFunctionReturn(0);
28271b4ebd2SPeter Brune         }
283b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
2848a430afbSPeter Brune           t1 = .5*(g - f) - lambda*initslope;
2858a430afbSPeter Brune           t2 = .5*(gprev  - f) - lambdaprev*initslope;
28671b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
28771b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
28871b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
28971b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
290f5af7f23SKarl Rupp           if (a == 0.0) lambdatemp = -initslope/(2.0*b);
291f5af7f23SKarl Rupp           else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
292f5af7f23SKarl Rupp 
293b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
2948a430afbSPeter Brune           lambdatemp = -initslope/(g - f - 2.0*initslope);
295ce94432eSBarry Smith         } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt");
29671b4ebd2SPeter Brune         lambdaprev = lambda;
2978a430afbSPeter Brune         gprev      = g;
29871b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
29971b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
30071b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
30171b4ebd2SPeter Brune         ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
302b7b2e573SPeter Brune         if (linesearch->ops->viproject) {
303b7b2e573SPeter Brune           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
304b7b2e573SPeter Brune         }
30571b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
30671b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
307bd4a8a71SBarry Smith           if (!objective) {
30871b4ebd2SPeter Brune             ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
309c69d1a72SBarry Smith                               (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
3108a430afbSPeter Brune           }
311f1c6b773SPeter Brune           ierr         = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
31271b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
31371b4ebd2SPeter Brune           PetscFunctionReturn(0);
31471b4ebd2SPeter Brune         }
315bd4a8a71SBarry Smith         if (objective) {
3168a430afbSPeter Brune           ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
3178a430afbSPeter Brune         } else {
31871b4ebd2SPeter Brune           ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
319c427506bSPeter Brune           ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
320c427506bSPeter Brune           if (domainerror) {
321f1c6b773SPeter Brune             ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
32271b4ebd2SPeter Brune             PetscFunctionReturn(0);
32371b4ebd2SPeter Brune           }
3249bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
3259bd66eb0SPeter Brune             gnorm = fnorm;
3269bd66eb0SPeter Brune             ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3279bd66eb0SPeter Brune           } else {
32871b4ebd2SPeter Brune             ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3298a430afbSPeter Brune           }
330f1ab69e3SPeter Brune           g = PetscSqr(gnorm);
3319bd66eb0SPeter Brune         }
332c98378a5SBarry Smith         if (PetscIsInfOrNanReal(gnorm)) {
333c98378a5SBarry Smith           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
334*c61ba1e2SPeter Brune           ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
335c98378a5SBarry Smith           PetscFunctionReturn(0);
336c98378a5SBarry Smith         }
3378a430afbSPeter Brune         if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */
33871b4ebd2SPeter Brune           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: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
34371b4ebd2SPeter Brune               } else {
344c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.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: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3508a430afbSPeter Brune               } else {
3518a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, obj %14.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           break;
357f5af7f23SKarl Rupp         } else if (monitor) {
358ea5d4fccSPeter Brune           ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
359bd4a8a71SBarry Smith           if (!objective) {
360b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
361c69d1a72SBarry 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);
36271b4ebd2SPeter Brune             } else {
363c69d1a72SBarry 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);
36471b4ebd2SPeter Brune             }
365ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3668a430afbSPeter Brune           } else {
3678a430afbSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3688a430afbSPeter 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);
3698a430afbSPeter Brune             } else {
3708a430afbSPeter 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);
3718a430afbSPeter Brune             }
3728a430afbSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3738a430afbSPeter Brune           }
37471b4ebd2SPeter Brune         }
37571b4ebd2SPeter Brune       }
37671b4ebd2SPeter Brune     }
37771b4ebd2SPeter Brune   }
37871b4ebd2SPeter Brune 
37971b4ebd2SPeter Brune   /* postcheck */
3807b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
38171b4ebd2SPeter Brune   if (changed_y) {
38271b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
3839bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
3849bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
3859bd66eb0SPeter Brune     }
38671b4ebd2SPeter Brune   }
387bd4a8a71SBarry Smith   if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */
388c427506bSPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
38971b4ebd2SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
39071b4ebd2SPeter Brune     if (domainerror) {
391f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
39271b4ebd2SPeter Brune       PetscFunctionReturn(0);
39371b4ebd2SPeter Brune     }
3949bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3959bd66eb0SPeter Brune       gnorm = fnorm;
3969bd66eb0SPeter Brune       ierr  = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3979bd66eb0SPeter Brune     } else {
3989bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3999bd66eb0SPeter Brune     }
4009bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
401c98378a5SBarry Smith     if (PetscIsInfOrNanReal(gnorm)) {
402c98378a5SBarry Smith       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
403*c61ba1e2SPeter Brune       ierr = PetscInfo(snes,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr);
404c98378a5SBarry Smith       PetscFunctionReturn(0);
405c98378a5SBarry Smith     }
406c427506bSPeter Brune   }
40771b4ebd2SPeter Brune 
40871b4ebd2SPeter Brune   /* copy the solution over */
40971b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
410c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
411c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
412f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
413f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
41471b4ebd2SPeter Brune   PetscFunctionReturn(0);
41571b4ebd2SPeter Brune }
41671b4ebd2SPeter Brune 
41771b4ebd2SPeter Brune #undef __FUNCT__
4187f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT"
4197f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
4207f1410a3SPeter Brune {
4217f1410a3SPeter Brune   PetscErrorCode    ierr;
4227f1410a3SPeter Brune   PetscBool         iascii;
4237f1410a3SPeter Brune   SNESLineSearch_BT *bt;
424075cc632SBarry Smith 
4257f1410a3SPeter Brune   PetscFunctionBegin;
426251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4277f1410a3SPeter Brune   bt   = (SNESLineSearch_BT*)linesearch->data;
4287f1410a3SPeter Brune   if (iascii) {
429b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
4307f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
431b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
4327f1410a3SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
4337f1410a3SPeter Brune     }
434007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", bt->alpha);CHKERRQ(ierr);
4357f1410a3SPeter Brune   }
4367f1410a3SPeter Brune   PetscFunctionReturn(0);
4377f1410a3SPeter Brune }
4387f1410a3SPeter Brune 
4397f1410a3SPeter Brune 
4407f1410a3SPeter Brune #undef __FUNCT__
441f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT"
442f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
443731c067cSJed Brown {
444731c067cSJed Brown   PetscErrorCode ierr;
44571b4ebd2SPeter Brune 
44671b4ebd2SPeter Brune   PetscFunctionBegin;
447731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
44871b4ebd2SPeter Brune   PetscFunctionReturn(0);
44971b4ebd2SPeter Brune }
45071b4ebd2SPeter Brune 
45171b4ebd2SPeter Brune 
45271b4ebd2SPeter Brune #undef __FUNCT__
453f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
454f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch)
45571b4ebd2SPeter Brune {
45671b4ebd2SPeter Brune 
45771b4ebd2SPeter Brune   PetscErrorCode    ierr;
458f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
45971b4ebd2SPeter Brune 
4606e111a19SKarl Rupp   PetscFunctionBegin;
461f1c6b773SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
46271b4ebd2SPeter Brune 
463f1c6b773SPeter Brune   ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr);
4640298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr);
46571b4ebd2SPeter Brune 
46671b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
46771b4ebd2SPeter Brune   PetscFunctionReturn(0);
46871b4ebd2SPeter Brune }
46971b4ebd2SPeter Brune 
47071b4ebd2SPeter Brune 
47171b4ebd2SPeter Brune #undef __FUNCT__
472f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT"
473954494b2SJed Brown /*MC
474db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
475954494b2SJed Brown 
476db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
477db609ea7SPeter Brune    function. If this fit does not satisfy the conditions for progress, the interval shrinks
478db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
479954494b2SJed Brown 
480cd7522eaSPeter Brune    Options Database Keys:
481cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
482db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
483db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
484db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
485cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
486cd7522eaSPeter Brune 
487954494b2SJed Brown    Level: advanced
488954494b2SJed Brown 
489db609ea7SPeter Brune    Notes:
490db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
491db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
492db609ea7SPeter Brune 
493f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
494954494b2SJed Brown 
495f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
496954494b2SJed Brown M*/
4978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
49871b4ebd2SPeter Brune {
49971b4ebd2SPeter Brune 
500f1c6b773SPeter Brune   SNESLineSearch_BT *bt;
50171b4ebd2SPeter Brune   PetscErrorCode    ierr;
50271b4ebd2SPeter Brune 
50371b4ebd2SPeter Brune   PetscFunctionBegin;
504f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
505f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
506f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
5070298fd71SBarry Smith   linesearch->ops->reset          = NULL;
5087f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
5090298fd71SBarry Smith   linesearch->ops->setup          = NULL;
51071b4ebd2SPeter Brune 
511f1c6b773SPeter Brune   ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr);
512f5af7f23SKarl Rupp 
51371b4ebd2SPeter Brune   linesearch->data    = (void*)bt;
51471b4ebd2SPeter Brune   linesearch->max_its = 40;
515b000cd8dSPeter Brune   linesearch->order   = SNES_LINESEARCH_ORDER_CUBIC;
51671b4ebd2SPeter Brune   bt->alpha           = 1e-4;
51771b4ebd2SPeter Brune   PetscFunctionReturn(0);
51871b4ebd2SPeter Brune }
519