xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision 251f4c6745e69799cb0300fd80598c58f96eb452)
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;
6571b4ebd2SPeter Brune   PetscReal      fnorm, xnorm, ynorm, gnorm, gnormprev;
6671b4ebd2SPeter Brune   PetscReal      lambda, lambdatemp, lambdaprev, minlambda, maxstep, rellength, initslope, alpha;
6771b4ebd2SPeter Brune   PetscReal      t1, t2, a, b, d, steptol;
6871b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX)
6971b4ebd2SPeter Brune   PetscScalar    cinitslope;
7071b4ebd2SPeter Brune #endif
7171b4ebd2SPeter Brune   PetscBool      domainerror;
7271b4ebd2SPeter Brune   PetscViewer    monitor;
7371b4ebd2SPeter Brune   PetscInt       max_its, count;
74f1c6b773SPeter Brune   SNESLineSearch_BT  *bt;
7571b4ebd2SPeter Brune   Mat            jac;
7671b4ebd2SPeter Brune 
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);
85f1c6b773SPeter Brune   ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_its);
86f1c6b773SPeter Brune   bt = (SNESLineSearch_BT *)linesearch->data;
8771b4ebd2SPeter Brune 
8871b4ebd2SPeter Brune   alpha = bt->alpha;
8971b4ebd2SPeter Brune 
9071b4ebd2SPeter Brune   ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
9171b4ebd2SPeter Brune   if (!jac) {
92f1c6b773SPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
9371b4ebd2SPeter Brune   }
9471b4ebd2SPeter Brune   /* precheck */
95f1c6b773SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
96f1c6b773SPeter Brune   ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
9771b4ebd2SPeter Brune 
9871b4ebd2SPeter Brune   ierr = VecNorm(Y, NORM_2, &ynorm);CHKERRQ(ierr);
9971b4ebd2SPeter Brune   if (ynorm == 0.0) {
10071b4ebd2SPeter Brune     if (monitor) {
101ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10271b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
103ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
10471b4ebd2SPeter Brune     }
10571b4ebd2SPeter Brune     ierr   = VecCopy(X,W);CHKERRQ(ierr);
10671b4ebd2SPeter Brune     ierr   = VecCopy(F,G);CHKERRQ(ierr);
107f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
10871b4ebd2SPeter Brune     PetscFunctionReturn(0);
10971b4ebd2SPeter Brune   }
11071b4ebd2SPeter Brune   if (ynorm > maxstep) {	/* Step too big, so scale back */
11171b4ebd2SPeter Brune     if (monitor) {
112ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
11371b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);CHKERRQ(ierr);
114ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
11571b4ebd2SPeter Brune     }
11671b4ebd2SPeter Brune     ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
11771b4ebd2SPeter Brune     ynorm = maxstep;
11871b4ebd2SPeter Brune   }
11971b4ebd2SPeter Brune   ierr      = VecMaxPointwiseDivide(Y,X,&rellength);CHKERRQ(ierr);
12071b4ebd2SPeter Brune   minlambda = steptol/rellength;
12171b4ebd2SPeter Brune   ierr      = MatMult(jac,Y,W);CHKERRQ(ierr);
12271b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX)
12371b4ebd2SPeter Brune   ierr      = VecDot(F,W,&cinitslope);CHKERRQ(ierr);
12471b4ebd2SPeter Brune   initslope = PetscRealPart(cinitslope);
12571b4ebd2SPeter Brune #else
12671b4ebd2SPeter Brune   ierr      = VecDot(F,W,&initslope);CHKERRQ(ierr);
12771b4ebd2SPeter Brune #endif
12871b4ebd2SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
12971b4ebd2SPeter Brune   if (initslope == 0.0) initslope = -1.0;
13071b4ebd2SPeter Brune 
131c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1329bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1339bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1349bd66eb0SPeter Brune   }
13571b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
13671b4ebd2SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
13771b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
138f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
13971b4ebd2SPeter Brune     PetscFunctionReturn(0);
14071b4ebd2SPeter Brune   }
14171b4ebd2SPeter Brune   ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
142c427506bSPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
143c427506bSPeter Brune   if (domainerror) {
144f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
14571b4ebd2SPeter Brune     PetscFunctionReturn(0);
14671b4ebd2SPeter Brune   }
1479bd66eb0SPeter Brune   if (linesearch->ops->vinorm) {
1489bd66eb0SPeter Brune     gnorm = fnorm;
1499bd66eb0SPeter Brune     ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1509bd66eb0SPeter Brune   } else {
15171b4ebd2SPeter Brune     ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1529bd66eb0SPeter Brune   }
1539bd66eb0SPeter Brune 
15471b4ebd2SPeter Brune   if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
15571b4ebd2SPeter Brune   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr);
15671b4ebd2SPeter Brune   if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */
15771b4ebd2SPeter Brune     if (monitor) {
158ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
15971b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr);
160ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
16171b4ebd2SPeter Brune     }
16271b4ebd2SPeter Brune   } else {
16371b4ebd2SPeter Brune     /* Fit points with quadratic */
16471b4ebd2SPeter Brune     lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope);
16571b4ebd2SPeter Brune     lambdaprev = lambda;
16671b4ebd2SPeter Brune     gnormprev  = gnorm;
16771b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
16871b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
16971b4ebd2SPeter Brune     else                         lambda = lambdatemp;
17071b4ebd2SPeter Brune 
17171b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1729bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1739bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1749bd66eb0SPeter Brune     }
17571b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
17671b4ebd2SPeter Brune       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
17771b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
178f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
17971b4ebd2SPeter Brune       PetscFunctionReturn(0);
18071b4ebd2SPeter Brune     }
18171b4ebd2SPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
1826188f407SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1836188f407SPeter Brune     if (domainerror) {
18471b4ebd2SPeter Brune       PetscFunctionReturn(0);
18571b4ebd2SPeter Brune     }
1869bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1879bd66eb0SPeter Brune       gnorm = fnorm;
1889bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1899bd66eb0SPeter Brune     } else {
19071b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1919bd66eb0SPeter Brune     }
19271b4ebd2SPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
19371b4ebd2SPeter Brune     if (monitor) {
194ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19571b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr);
196ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19771b4ebd2SPeter Brune     }
19871b4ebd2SPeter Brune     if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
19971b4ebd2SPeter Brune       if (monitor) {
200ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
20171b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
202ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
20371b4ebd2SPeter Brune       }
20471b4ebd2SPeter Brune     } else {
20571b4ebd2SPeter Brune       /* Fit points with cubic */
20671b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
20771b4ebd2SPeter Brune         if (lambda <= minlambda) {
20871b4ebd2SPeter Brune           if (monitor) {
209ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
21071b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
21171b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,
21271b4ebd2SPeter Brune                                           "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
21371b4ebd2SPeter Brune                                           fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr);
214ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
21571b4ebd2SPeter Brune           }
216f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
21771b4ebd2SPeter Brune           PetscFunctionReturn(0);
21871b4ebd2SPeter Brune         }
219b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
22071b4ebd2SPeter Brune           t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope;
22171b4ebd2SPeter Brune           t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
22271b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
22371b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
22471b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
22571b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
22671b4ebd2SPeter Brune           if (a == 0.0) {
22771b4ebd2SPeter Brune             lambdatemp = -initslope/(2.0*b);
22871b4ebd2SPeter Brune           } else {
22971b4ebd2SPeter Brune             lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
23071b4ebd2SPeter Brune           }
231b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
23271b4ebd2SPeter Brune           lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope);
23359405d5eSPeter Brune         } else {
23459405d5eSPeter Brune           SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt");
23571b4ebd2SPeter Brune         }
23671b4ebd2SPeter Brune           lambdaprev = lambda;
23771b4ebd2SPeter Brune           gnormprev  = gnorm;
23871b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
23971b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
24071b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
24171b4ebd2SPeter Brune         ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
24271b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
24371b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
24471b4ebd2SPeter Brune           ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
24571b4ebd2SPeter Brune                             fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr);
246f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
24771b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
24871b4ebd2SPeter Brune           PetscFunctionReturn(0);
24971b4ebd2SPeter Brune         }
25071b4ebd2SPeter Brune         ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
251c427506bSPeter Brune         ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
252c427506bSPeter Brune         if (domainerror) {
253f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
25471b4ebd2SPeter Brune           PetscFunctionReturn(0);
25571b4ebd2SPeter Brune         }
2569bd66eb0SPeter Brune         if (linesearch->ops->vinorm) {
2579bd66eb0SPeter Brune           gnorm = fnorm;
2589bd66eb0SPeter Brune           ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2599bd66eb0SPeter Brune         } else {
26071b4ebd2SPeter Brune           ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2619bd66eb0SPeter Brune         }
26271b4ebd2SPeter Brune         if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
26371b4ebd2SPeter Brune         if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
26471b4ebd2SPeter Brune           if (monitor) {
265ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
266b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
26771b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
26871b4ebd2SPeter Brune             } else {
26971b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
27071b4ebd2SPeter Brune             }
271ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
27271b4ebd2SPeter Brune           }
27371b4ebd2SPeter Brune           break;
27471b4ebd2SPeter Brune         } else {
27571b4ebd2SPeter Brune           if (monitor) {
276ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
277b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
27871b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
27971b4ebd2SPeter Brune             } else {
28071b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
28171b4ebd2SPeter Brune             }
282ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
28371b4ebd2SPeter Brune           }
28471b4ebd2SPeter Brune         }
28571b4ebd2SPeter Brune       }
28671b4ebd2SPeter Brune     }
28771b4ebd2SPeter Brune   }
28871b4ebd2SPeter Brune 
28971b4ebd2SPeter Brune   /* postcheck */
290f1c6b773SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
29171b4ebd2SPeter Brune   if (changed_y) {
29271b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2939bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2949bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2959bd66eb0SPeter Brune     }
29671b4ebd2SPeter Brune   }
297c427506bSPeter Brune   if (changed_y || changed_w) { /* recompute the function if the step has changed */
298c427506bSPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
29971b4ebd2SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
30071b4ebd2SPeter Brune     if (domainerror) {
301f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
30271b4ebd2SPeter Brune       PetscFunctionReturn(0);
30371b4ebd2SPeter Brune     }
3049bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3059bd66eb0SPeter Brune       gnorm = fnorm;
3069bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3079bd66eb0SPeter Brune     } else {
3089bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3099bd66eb0SPeter Brune     }
3109bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
311c427506bSPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3129bd66eb0SPeter Brune 
313c427506bSPeter Brune   }
31471b4ebd2SPeter Brune 
31571b4ebd2SPeter Brune   /* copy the solution over */
31671b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
317c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
318c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
319f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
320f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
32171b4ebd2SPeter Brune   PetscFunctionReturn(0);
32271b4ebd2SPeter Brune }
32371b4ebd2SPeter Brune 
32471b4ebd2SPeter Brune #undef __FUNCT__
3257f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT"
3267f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
3277f1410a3SPeter Brune {
3287f1410a3SPeter Brune   PetscErrorCode    ierr;
3297f1410a3SPeter Brune   PetscBool         iascii;
3307f1410a3SPeter Brune   SNESLineSearch_BT *bt;
3317f1410a3SPeter Brune   PetscFunctionBegin;
332*251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
3337f1410a3SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
3347f1410a3SPeter Brune   if (iascii) {
335b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3367f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
337b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
3387f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
3397f1410a3SPeter Brune     }
3407f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%G\n", bt->alpha);CHKERRQ(ierr);
3417f1410a3SPeter Brune   }
3427f1410a3SPeter Brune   PetscFunctionReturn(0);
3437f1410a3SPeter Brune }
3447f1410a3SPeter Brune 
3457f1410a3SPeter Brune 
3467f1410a3SPeter Brune #undef __FUNCT__
347f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT"
348f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
349731c067cSJed Brown {
350731c067cSJed Brown   PetscErrorCode ierr;
35171b4ebd2SPeter Brune 
35271b4ebd2SPeter Brune   PetscFunctionBegin;
353731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
35471b4ebd2SPeter Brune   PetscFunctionReturn(0);
35571b4ebd2SPeter Brune }
35671b4ebd2SPeter Brune 
35771b4ebd2SPeter Brune 
35871b4ebd2SPeter Brune #undef __FUNCT__
359f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
360f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch)
36171b4ebd2SPeter Brune {
36271b4ebd2SPeter Brune 
36371b4ebd2SPeter Brune   PetscErrorCode       ierr;
364f1c6b773SPeter Brune   SNESLineSearch_BT    *bt;
36571b4ebd2SPeter Brune   PetscFunctionBegin;
36671b4ebd2SPeter Brune 
367f1c6b773SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
36871b4ebd2SPeter Brune 
369f1c6b773SPeter Brune   ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr);
3707a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr);
37171b4ebd2SPeter Brune 
37271b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
37371b4ebd2SPeter Brune   PetscFunctionReturn(0);
37471b4ebd2SPeter Brune }
37571b4ebd2SPeter Brune 
37671b4ebd2SPeter Brune 
37771b4ebd2SPeter Brune #undef __FUNCT__
378f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT"
379954494b2SJed Brown /*MC
380db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
381954494b2SJed Brown 
382db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
383db609ea7SPeter Brune    function. If this fit does not satisfy the conditions for progress, the interval shrinks
384db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
385954494b2SJed Brown 
386cd7522eaSPeter Brune    Options Database Keys:
387cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
388db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
389db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
390db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
391cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
392cd7522eaSPeter Brune 
393954494b2SJed Brown    Level: advanced
394954494b2SJed Brown 
395db609ea7SPeter Brune    Notes:
396db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
397db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
398db609ea7SPeter Brune 
399f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
400954494b2SJed Brown 
401f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
402954494b2SJed Brown M*/
403f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
40471b4ebd2SPeter Brune {
40571b4ebd2SPeter Brune 
406f1c6b773SPeter Brune   SNESLineSearch_BT  *bt;
40771b4ebd2SPeter Brune   PetscErrorCode ierr;
40871b4ebd2SPeter Brune 
40971b4ebd2SPeter Brune   PetscFunctionBegin;
410f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
411f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
412f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
41371b4ebd2SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
4147f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
41571b4ebd2SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
41671b4ebd2SPeter Brune 
417f1c6b773SPeter Brune   ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr);
41871b4ebd2SPeter Brune   linesearch->data = (void *)bt;
41971b4ebd2SPeter Brune   linesearch->max_its = 40;
420b000cd8dSPeter Brune   linesearch->order = SNES_LINESEARCH_ORDER_CUBIC;
42171b4ebd2SPeter Brune   bt->alpha = 1e-4;
42271b4ebd2SPeter Brune 
42371b4ebd2SPeter Brune   PetscFunctionReturn(0);
42471b4ebd2SPeter Brune }
425