xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision f1ab69e3b3a0b45fb61c0e722fed61bbd3dde8fd)
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;
242f4102e2SPeter Brune   PetscFunctionBegin;
252f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
262f4102e2SPeter Brune   bt = (SNESLineSearch_BT *)linesearch->data;
272f4102e2SPeter Brune   bt->alpha = alpha;
282f4102e2SPeter Brune   PetscFunctionReturn(0);
292f4102e2SPeter Brune }
302f4102e2SPeter Brune 
312f4102e2SPeter Brune 
322f4102e2SPeter Brune #undef __FUNCT__
332f4102e2SPeter Brune #define __FUNCT__ "SNESLineSearchBTGetAlpha"
342f4102e2SPeter Brune /*@
352f4102e2SPeter Brune    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant.
362f4102e2SPeter Brune 
372f4102e2SPeter Brune    Input Parameters:
382f4102e2SPeter Brune .  linesearch - linesearch context
398e557f58SPeter Brune 
408e557f58SPeter Brune    Output Parameters:
412f4102e2SPeter Brune .  alpha - The descent parameter
422f4102e2SPeter Brune 
432f4102e2SPeter Brune    Level: intermediate
442f4102e2SPeter Brune 
451a4f838cSPeter Brune .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
462f4102e2SPeter Brune @*/
472f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
482f4102e2SPeter Brune {
492f4102e2SPeter Brune   SNESLineSearch_BT  *bt;
502f4102e2SPeter Brune   PetscFunctionBegin;
512f4102e2SPeter Brune   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
522f4102e2SPeter Brune   bt = (SNESLineSearch_BT *)linesearch->data;
532f4102e2SPeter Brune   *alpha = bt->alpha;
542f4102e2SPeter Brune   PetscFunctionReturn(0);
552f4102e2SPeter Brune }
562f4102e2SPeter Brune 
572f4102e2SPeter Brune #undef __FUNCT__
58f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_BT"
59f1c6b773SPeter Brune static PetscErrorCode  SNESLineSearchApply_BT(SNESLineSearch linesearch)
6071b4ebd2SPeter Brune {
6171b4ebd2SPeter Brune   PetscBool      changed_y,changed_w;
6271b4ebd2SPeter Brune   PetscErrorCode ierr;
6371b4ebd2SPeter Brune   Vec            X,F,Y,W,G;
6471b4ebd2SPeter Brune   SNES           snes;
658a430afbSPeter Brune   PetscReal      fnorm, xnorm, ynorm, gnorm;
663b288129SPeter Brune   PetscReal      lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol;
67dfcd3828SPeter Brune   PetscReal      t1,t2,a,b,d;
688a430afbSPeter Brune   PetscReal      f;
698a430afbSPeter Brune   PetscReal      g,gprev;
7071b4ebd2SPeter Brune   PetscScalar    cinitslope;
7171b4ebd2SPeter Brune   PetscBool      domainerror;
7271b4ebd2SPeter Brune   PetscViewer    monitor;
7371b4ebd2SPeter Brune   PetscInt       max_its,count;
74f1c6b773SPeter Brune   SNESLineSearch_BT  *bt;
7571b4ebd2SPeter Brune   Mat            jac;
768a430afbSPeter Brune   SNESObjective  obj;
7771b4ebd2SPeter Brune 
7871b4ebd2SPeter Brune   PetscFunctionBegin;
7971b4ebd2SPeter Brune 
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);
85dfcd3828SPeter Brune   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,PETSC_NULL,PETSC_NULL,PETSC_NULL,&max_its);
863b288129SPeter Brune   ierr = SNESGetTolerances(snes,PETSC_NULL,PETSC_NULL,&stol,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
878a430afbSPeter Brune   ierr = SNESGetObjective(snes,&obj,PETSC_NULL);CHKERRQ(ierr);
88f1c6b773SPeter Brune   bt = (SNESLineSearch_BT *)linesearch->data;
8971b4ebd2SPeter Brune 
9071b4ebd2SPeter Brune   alpha = bt->alpha;
9171b4ebd2SPeter Brune 
9271b4ebd2SPeter Brune   ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
938a430afbSPeter Brune 
948a430afbSPeter Brune   if (!jac && !obj) SETERRQ(((PetscObject)linesearch)->comm, 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);
113f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);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 */
1278a430afbSPeter Brune   if (obj) {
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 */
1348a430afbSPeter Brune   if (obj) {
1358a430afbSPeter Brune     /* slope comes from the function (assumed to be the gradient of the objective */
1368a430afbSPeter Brune     ierr = VecDot(Y,F,&cinitslope);CHKERRQ(ierr);
1378a430afbSPeter Brune     initslope = PetscRealPart(cinitslope);
1388a430afbSPeter Brune   } else {
1398a430afbSPeter Brune     /* slope comes from the normal equations */
14071b4ebd2SPeter Brune     ierr      = MatMult(jac,Y,W);CHKERRQ(ierr);
14171b4ebd2SPeter Brune     ierr      = VecDot(F,W,&cinitslope);CHKERRQ(ierr);
14271b4ebd2SPeter Brune     initslope = PetscRealPart(cinitslope);
14371b4ebd2SPeter Brune     if (initslope > 0.0)  initslope = -initslope;
14471b4ebd2SPeter Brune     if (initslope == 0.0) initslope = -1.0;
1458a430afbSPeter Brune   }
14671b4ebd2SPeter Brune 
147c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1489bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1499bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1509bd66eb0SPeter Brune   }
15171b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
15271b4ebd2SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
15371b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
154f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
15571b4ebd2SPeter Brune     PetscFunctionReturn(0);
15671b4ebd2SPeter Brune   }
1578a430afbSPeter Brune 
1588a430afbSPeter Brune   if (obj) {
1598a430afbSPeter Brune     ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
1608a430afbSPeter Brune   } else {
16171b4ebd2SPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
162c427506bSPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
163c427506bSPeter Brune     if (domainerror) {
164f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
16571b4ebd2SPeter Brune       PetscFunctionReturn(0);
16671b4ebd2SPeter Brune     }
1679bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1689bd66eb0SPeter Brune       gnorm = fnorm;
1699bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1709bd66eb0SPeter Brune     } else {
17171b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1728a430afbSPeter Brune     }
173*f1ab69e3SPeter Brune     g = PetscSqr(gnorm);
1749bd66eb0SPeter Brune   }
1759bd66eb0SPeter Brune 
1768a430afbSPeter Brune   if (PetscIsInfOrNanReal(g)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1778a430afbSPeter Brune   if (!obj) {
178c69d1a72SBarry Smith     ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1798a430afbSPeter Brune   }
1808a430afbSPeter Brune   if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */
18171b4ebd2SPeter Brune     if (monitor) {
182ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
1838a430afbSPeter Brune       if (!obj) {
184c69d1a72SBarry Smith         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
1858a430afbSPeter Brune       } else {
1868a430afbSPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr);
1878a430afbSPeter Brune       }
188ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
18971b4ebd2SPeter Brune     }
19071b4ebd2SPeter Brune   } else {
191c21ba15cSPeter Brune     /* Since the full step didn't work and the step is tiny, quit */
192c21ba15cSPeter Brune     if (stol*xnorm > ynorm) {
193c21ba15cSPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
194c21ba15cSPeter Brune       ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr);
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;
212f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
21371b4ebd2SPeter Brune       PetscFunctionReturn(0);
21471b4ebd2SPeter Brune     }
2158a430afbSPeter Brune     if (obj) {
2168a430afbSPeter Brune       ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
2178a430afbSPeter Brune     } else {
21871b4ebd2SPeter Brune       ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
2196188f407SPeter Brune       ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
2206188f407SPeter Brune       if (domainerror) {
22171b4ebd2SPeter Brune         PetscFunctionReturn(0);
22271b4ebd2SPeter Brune       }
2239bd66eb0SPeter Brune       if (linesearch->ops->vinorm) {
2249bd66eb0SPeter Brune         gnorm = fnorm;
2259bd66eb0SPeter Brune         ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2269bd66eb0SPeter Brune       } else {
22771b4ebd2SPeter Brune         ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2289bd66eb0SPeter Brune       }
229*f1ab69e3SPeter Brune       g = PetscSqr(gnorm);
2308a430afbSPeter Brune     }
2318a430afbSPeter Brune     if (PetscIsInfOrNanReal(g)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
23271b4ebd2SPeter Brune     if (monitor) {
233ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
2348a430afbSPeter Brune       if (!obj) {
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);
2548a430afbSPeter Brune             if (!obj) {
25571b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,
25671b4ebd2SPeter Brune                                           "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
257c69d1a72SBarry Smith                                           (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2588a430afbSPeter Brune             } else {
2598a430afbSPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,
2608a430afbSPeter Brune                                             "    Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
2618a430afbSPeter Brune                                             (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
2628a430afbSPeter Brune             }
263ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
26471b4ebd2SPeter Brune           }
265f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
26671b4ebd2SPeter Brune           PetscFunctionReturn(0);
26771b4ebd2SPeter Brune         }
268b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
2698a430afbSPeter Brune           t1 = .5*(g - f) - lambda*initslope;
2708a430afbSPeter Brune           t2 = .5*(gprev  - f) - lambdaprev*initslope;
27171b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27271b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
27371b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
27471b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
27571b4ebd2SPeter Brune           if (a == 0.0) {
27671b4ebd2SPeter Brune             lambdatemp = -initslope/(2.0*b);
27771b4ebd2SPeter Brune           } else {
27871b4ebd2SPeter Brune             lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
27971b4ebd2SPeter Brune           }
280b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
2818a430afbSPeter Brune           lambdatemp = -initslope/(g - f - 2.0*initslope);
28259405d5eSPeter Brune         } else {
28359405d5eSPeter Brune           SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt");
28471b4ebd2SPeter Brune         }
28571b4ebd2SPeter Brune           lambdaprev = lambda;
2868a430afbSPeter Brune           gprev  = g;
28771b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
28871b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
28971b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
29071b4ebd2SPeter Brune         ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
291b7b2e573SPeter Brune         if (linesearch->ops->viproject) {
292b7b2e573SPeter Brune           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
293b7b2e573SPeter Brune         }
29471b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
29571b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
2968a430afbSPeter Brune           if (!obj) {
29771b4ebd2SPeter Brune             ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
298c69d1a72SBarry Smith                               (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
2998a430afbSPeter Brune           }
300f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
30171b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
30271b4ebd2SPeter Brune           PetscFunctionReturn(0);
30371b4ebd2SPeter Brune         }
3048a430afbSPeter Brune         if (obj) {
3058a430afbSPeter Brune           ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr);
3068a430afbSPeter Brune         } else {
30771b4ebd2SPeter Brune           ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
308c427506bSPeter Brune           ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
309c427506bSPeter Brune           if (domainerror) {
310f1c6b773SPeter Brune             ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
31171b4ebd2SPeter Brune             PetscFunctionReturn(0);
31271b4ebd2SPeter Brune           }
3139bd66eb0SPeter Brune           if (linesearch->ops->vinorm) {
3149bd66eb0SPeter Brune             gnorm = fnorm;
3159bd66eb0SPeter Brune             ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3169bd66eb0SPeter Brune           } else {
31771b4ebd2SPeter Brune             ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3188a430afbSPeter Brune           }
319*f1ab69e3SPeter Brune           g = PetscSqr(gnorm);
3209bd66eb0SPeter Brune         }
32171b4ebd2SPeter Brune         if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3228a430afbSPeter Brune         if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */
32371b4ebd2SPeter Brune           if (monitor) {
324ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3258a430afbSPeter Brune             if (!obj) {
326b000cd8dSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
327c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
32871b4ebd2SPeter Brune               } else {
329c69d1a72SBarry Smith                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
33071b4ebd2SPeter Brune               }
331ea5d4fccSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3328a430afbSPeter Brune             } else {
3338a430afbSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3348a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3358a430afbSPeter Brune               } else {
3368a430afbSPeter Brune                 ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr);
3378a430afbSPeter Brune               }
3388a430afbSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3398a430afbSPeter Brune             }
34071b4ebd2SPeter Brune           }
34171b4ebd2SPeter Brune           break;
34271b4ebd2SPeter Brune         } else {
34371b4ebd2SPeter Brune           if (monitor) {
344ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3458a430afbSPeter Brune             if (!obj) {
346b000cd8dSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
347c69d1a72SBarry 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);
34871b4ebd2SPeter Brune               } else {
349c69d1a72SBarry 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);
35071b4ebd2SPeter Brune               }
351ea5d4fccSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3528a430afbSPeter Brune             } else {
3538a430afbSPeter Brune               if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3548a430afbSPeter 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);
3558a430afbSPeter Brune               } else {
3568a430afbSPeter 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);
3578a430afbSPeter Brune               }
3588a430afbSPeter Brune               ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
3598a430afbSPeter Brune             }
36071b4ebd2SPeter Brune           }
36171b4ebd2SPeter Brune         }
36271b4ebd2SPeter Brune       }
36371b4ebd2SPeter Brune     }
36471b4ebd2SPeter Brune   }
36571b4ebd2SPeter Brune 
36671b4ebd2SPeter Brune   /* postcheck */
3677b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
36871b4ebd2SPeter Brune   if (changed_y) {
36971b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
3709bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
3719bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
3729bd66eb0SPeter Brune     }
37371b4ebd2SPeter Brune   }
3748a430afbSPeter Brune   if (changed_y || changed_w || obj) { /* recompute the function norm if the step has changed or the objective isn't the norm */
375c427506bSPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
37671b4ebd2SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
37771b4ebd2SPeter Brune     if (domainerror) {
378f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
37971b4ebd2SPeter Brune       PetscFunctionReturn(0);
38071b4ebd2SPeter Brune     }
3819bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3829bd66eb0SPeter Brune       gnorm = fnorm;
3839bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3849bd66eb0SPeter Brune     } else {
3859bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3869bd66eb0SPeter Brune     }
3879bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
388c427506bSPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3899bd66eb0SPeter Brune 
390c427506bSPeter Brune   }
39171b4ebd2SPeter Brune 
39271b4ebd2SPeter Brune   /* copy the solution over */
39371b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
394c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
395c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
396f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
397f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
39871b4ebd2SPeter Brune   PetscFunctionReturn(0);
39971b4ebd2SPeter Brune }
40071b4ebd2SPeter Brune 
40171b4ebd2SPeter Brune #undef __FUNCT__
4027f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT"
4037f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
4047f1410a3SPeter Brune {
4057f1410a3SPeter Brune   PetscErrorCode    ierr;
4067f1410a3SPeter Brune   PetscBool         iascii;
4077f1410a3SPeter Brune   SNESLineSearch_BT *bt;
4087f1410a3SPeter Brune   PetscFunctionBegin;
409251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4107f1410a3SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
4117f1410a3SPeter Brune   if (iascii) {
412b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
4137f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
414b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
4157f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
4167f1410a3SPeter Brune     }
417007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", bt->alpha);CHKERRQ(ierr);
4187f1410a3SPeter Brune   }
4197f1410a3SPeter Brune   PetscFunctionReturn(0);
4207f1410a3SPeter Brune }
4217f1410a3SPeter Brune 
4227f1410a3SPeter Brune 
4237f1410a3SPeter Brune #undef __FUNCT__
424f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT"
425f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
426731c067cSJed Brown {
427731c067cSJed Brown   PetscErrorCode ierr;
42871b4ebd2SPeter Brune 
42971b4ebd2SPeter Brune   PetscFunctionBegin;
430731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
43171b4ebd2SPeter Brune   PetscFunctionReturn(0);
43271b4ebd2SPeter Brune }
43371b4ebd2SPeter Brune 
43471b4ebd2SPeter Brune 
43571b4ebd2SPeter Brune #undef __FUNCT__
436f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
437f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch)
43871b4ebd2SPeter Brune {
43971b4ebd2SPeter Brune 
44071b4ebd2SPeter Brune   PetscErrorCode       ierr;
441f1c6b773SPeter Brune   SNESLineSearch_BT    *bt;
44271b4ebd2SPeter Brune   PetscFunctionBegin;
44371b4ebd2SPeter Brune 
444f1c6b773SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
44571b4ebd2SPeter Brune 
446f1c6b773SPeter Brune   ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr);
4477a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr);
44871b4ebd2SPeter Brune 
44971b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
45071b4ebd2SPeter Brune   PetscFunctionReturn(0);
45171b4ebd2SPeter Brune }
45271b4ebd2SPeter Brune 
45371b4ebd2SPeter Brune 
45471b4ebd2SPeter Brune #undef __FUNCT__
455f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT"
456954494b2SJed Brown /*MC
457db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
458954494b2SJed Brown 
459db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
460db609ea7SPeter Brune    function. If this fit does not satisfy the conditions for progress, the interval shrinks
461db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
462954494b2SJed Brown 
463cd7522eaSPeter Brune    Options Database Keys:
464cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
465db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
466db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
467db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
468cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
469cd7522eaSPeter Brune 
470954494b2SJed Brown    Level: advanced
471954494b2SJed Brown 
472db609ea7SPeter Brune    Notes:
473db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
474db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
475db609ea7SPeter Brune 
476f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
477954494b2SJed Brown 
478f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
479954494b2SJed Brown M*/
480f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
48171b4ebd2SPeter Brune {
48271b4ebd2SPeter Brune 
483f1c6b773SPeter Brune   SNESLineSearch_BT  *bt;
48471b4ebd2SPeter Brune   PetscErrorCode ierr;
48571b4ebd2SPeter Brune 
48671b4ebd2SPeter Brune   PetscFunctionBegin;
487f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
488f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
489f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
49071b4ebd2SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
4917f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
49271b4ebd2SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
49371b4ebd2SPeter Brune 
494f1c6b773SPeter Brune   ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr);
49571b4ebd2SPeter Brune   linesearch->data = (void *)bt;
49671b4ebd2SPeter Brune   linesearch->max_its = 40;
497b000cd8dSPeter Brune   linesearch->order = SNES_LINESEARCH_ORDER_CUBIC;
49871b4ebd2SPeter Brune   bt->alpha = 1e-4;
49971b4ebd2SPeter Brune 
50071b4ebd2SPeter Brune   PetscFunctionReturn(0);
50171b4ebd2SPeter Brune }
502