xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision dfcd382804734d69f5ee5e440b71b7097e4d653f)
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;
66*dfcd3828SPeter Brune   PetscReal      lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha;
67*dfcd3828SPeter Brune   PetscReal      t1,t2,a,b,d;
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);
85*dfcd3828SPeter Brune   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&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 */
957b1df9c1SPeter Brune   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&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      = MatMult(jac,Y,W);CHKERRQ(ierr);
12071b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX)
12171b4ebd2SPeter Brune   ierr      = VecDot(F,W,&cinitslope);CHKERRQ(ierr);
12271b4ebd2SPeter Brune   initslope = PetscRealPart(cinitslope);
12371b4ebd2SPeter Brune #else
12471b4ebd2SPeter Brune   ierr      = VecDot(F,W,&initslope);CHKERRQ(ierr);
12571b4ebd2SPeter Brune #endif
12671b4ebd2SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
12771b4ebd2SPeter Brune   if (initslope == 0.0) initslope = -1.0;
12871b4ebd2SPeter Brune 
129c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1309bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1319bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1329bd66eb0SPeter Brune   }
13371b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
13471b4ebd2SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
13571b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
136f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
13771b4ebd2SPeter Brune     PetscFunctionReturn(0);
13871b4ebd2SPeter Brune   }
13971b4ebd2SPeter Brune   ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
140c427506bSPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
141c427506bSPeter Brune   if (domainerror) {
142f1c6b773SPeter Brune     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
14371b4ebd2SPeter Brune     PetscFunctionReturn(0);
14471b4ebd2SPeter Brune   }
1459bd66eb0SPeter Brune   if (linesearch->ops->vinorm) {
1469bd66eb0SPeter Brune     gnorm = fnorm;
1479bd66eb0SPeter Brune     ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1489bd66eb0SPeter Brune   } else {
14971b4ebd2SPeter Brune     ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1509bd66eb0SPeter Brune   }
1519bd66eb0SPeter Brune 
15271b4ebd2SPeter Brune   if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
15371b4ebd2SPeter Brune   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr);
15471b4ebd2SPeter Brune   if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */
15571b4ebd2SPeter Brune     if (monitor) {
156ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
15771b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr);
158ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
15971b4ebd2SPeter Brune     }
16071b4ebd2SPeter Brune   } else {
16171b4ebd2SPeter Brune     /* Fit points with quadratic */
16271b4ebd2SPeter Brune     lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope);
16371b4ebd2SPeter Brune     lambdaprev = lambda;
16471b4ebd2SPeter Brune     gnormprev  = gnorm;
16571b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
16671b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
16771b4ebd2SPeter Brune     else                         lambda = lambdatemp;
16871b4ebd2SPeter Brune 
16971b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1709bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1719bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1729bd66eb0SPeter Brune     }
17371b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
17471b4ebd2SPeter Brune       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
17571b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
176f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
17771b4ebd2SPeter Brune       PetscFunctionReturn(0);
17871b4ebd2SPeter Brune     }
17971b4ebd2SPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
1806188f407SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1816188f407SPeter Brune     if (domainerror) {
18271b4ebd2SPeter Brune       PetscFunctionReturn(0);
18371b4ebd2SPeter Brune     }
1849bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1859bd66eb0SPeter Brune       gnorm = fnorm;
1869bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1879bd66eb0SPeter Brune     } else {
18871b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1899bd66eb0SPeter Brune     }
19071b4ebd2SPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
19171b4ebd2SPeter Brune     if (monitor) {
192ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19371b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr);
194ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19571b4ebd2SPeter Brune     }
19671b4ebd2SPeter Brune     if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
19771b4ebd2SPeter Brune       if (monitor) {
198ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
19971b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
200ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
20171b4ebd2SPeter Brune       }
20271b4ebd2SPeter Brune     } else {
20371b4ebd2SPeter Brune       /* Fit points with cubic */
20471b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
20571b4ebd2SPeter Brune         if (lambda <= minlambda) {
20671b4ebd2SPeter Brune           if (monitor) {
207ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
20871b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
20971b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,
21071b4ebd2SPeter Brune                                           "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
21171b4ebd2SPeter Brune                                           fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr);
212ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
21371b4ebd2SPeter Brune           }
214f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
21571b4ebd2SPeter Brune           PetscFunctionReturn(0);
21671b4ebd2SPeter Brune         }
217b000cd8dSPeter Brune         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
21871b4ebd2SPeter Brune           t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope;
21971b4ebd2SPeter Brune           t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
22071b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
22171b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
22271b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
22371b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
22471b4ebd2SPeter Brune           if (a == 0.0) {
22571b4ebd2SPeter Brune             lambdatemp = -initslope/(2.0*b);
22671b4ebd2SPeter Brune           } else {
22771b4ebd2SPeter Brune             lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
22871b4ebd2SPeter Brune           }
229b000cd8dSPeter Brune         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
23071b4ebd2SPeter Brune           lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope);
23159405d5eSPeter Brune         } else {
23259405d5eSPeter Brune           SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt");
23371b4ebd2SPeter Brune         }
23471b4ebd2SPeter Brune           lambdaprev = lambda;
23571b4ebd2SPeter Brune           gnormprev  = gnorm;
23671b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
23771b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
23871b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
23971b4ebd2SPeter Brune         ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
24071b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
24171b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
24271b4ebd2SPeter Brune           ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
24371b4ebd2SPeter Brune                             fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr);
244f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
24571b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
24671b4ebd2SPeter Brune           PetscFunctionReturn(0);
24771b4ebd2SPeter Brune         }
24871b4ebd2SPeter Brune         ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
249c427506bSPeter Brune         ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
250c427506bSPeter Brune         if (domainerror) {
251f1c6b773SPeter Brune           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
25271b4ebd2SPeter Brune           PetscFunctionReturn(0);
25371b4ebd2SPeter Brune         }
2549bd66eb0SPeter Brune         if (linesearch->ops->vinorm) {
2559bd66eb0SPeter Brune           gnorm = fnorm;
2569bd66eb0SPeter Brune           ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2579bd66eb0SPeter Brune         } else {
25871b4ebd2SPeter Brune           ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2599bd66eb0SPeter Brune         }
26071b4ebd2SPeter Brune         if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
26171b4ebd2SPeter Brune         if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
26271b4ebd2SPeter Brune           if (monitor) {
263ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
264b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
26571b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
26671b4ebd2SPeter Brune             } else {
26771b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
26871b4ebd2SPeter Brune             }
269ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
27071b4ebd2SPeter Brune           }
27171b4ebd2SPeter Brune           break;
27271b4ebd2SPeter Brune         } else {
27371b4ebd2SPeter Brune           if (monitor) {
274ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
275b000cd8dSPeter Brune             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
27671b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
27771b4ebd2SPeter Brune             } else {
27871b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
27971b4ebd2SPeter Brune             }
280ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
28171b4ebd2SPeter Brune           }
28271b4ebd2SPeter Brune         }
28371b4ebd2SPeter Brune       }
28471b4ebd2SPeter Brune     }
28571b4ebd2SPeter Brune   }
28671b4ebd2SPeter Brune 
28771b4ebd2SPeter Brune   /* postcheck */
2887b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
28971b4ebd2SPeter Brune   if (changed_y) {
29071b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2919bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2929bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2939bd66eb0SPeter Brune     }
29471b4ebd2SPeter Brune   }
295c427506bSPeter Brune   if (changed_y || changed_w) { /* recompute the function if the step has changed */
296c427506bSPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
29771b4ebd2SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
29871b4ebd2SPeter Brune     if (domainerror) {
299f1c6b773SPeter Brune       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
30071b4ebd2SPeter Brune       PetscFunctionReturn(0);
30171b4ebd2SPeter Brune     }
3029bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
3039bd66eb0SPeter Brune       gnorm = fnorm;
3049bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
3059bd66eb0SPeter Brune     } else {
3069bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
3079bd66eb0SPeter Brune     }
3089bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
309c427506bSPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3109bd66eb0SPeter Brune 
311c427506bSPeter Brune   }
31271b4ebd2SPeter Brune 
31371b4ebd2SPeter Brune   /* copy the solution over */
31471b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
315c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
316c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
317f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
318f1c6b773SPeter Brune   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
31971b4ebd2SPeter Brune   PetscFunctionReturn(0);
32071b4ebd2SPeter Brune }
32171b4ebd2SPeter Brune 
32271b4ebd2SPeter Brune #undef __FUNCT__
3237f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT"
3247f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
3257f1410a3SPeter Brune {
3267f1410a3SPeter Brune   PetscErrorCode    ierr;
3277f1410a3SPeter Brune   PetscBool         iascii;
3287f1410a3SPeter Brune   SNESLineSearch_BT *bt;
3297f1410a3SPeter Brune   PetscFunctionBegin;
330251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
3317f1410a3SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
3327f1410a3SPeter Brune   if (iascii) {
333b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
3347f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
335b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
3367f1410a3SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
3377f1410a3SPeter Brune     }
338007b6e36SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", bt->alpha);CHKERRQ(ierr);
3397f1410a3SPeter Brune   }
3407f1410a3SPeter Brune   PetscFunctionReturn(0);
3417f1410a3SPeter Brune }
3427f1410a3SPeter Brune 
3437f1410a3SPeter Brune 
3447f1410a3SPeter Brune #undef __FUNCT__
345f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT"
346f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
347731c067cSJed Brown {
348731c067cSJed Brown   PetscErrorCode ierr;
34971b4ebd2SPeter Brune 
35071b4ebd2SPeter Brune   PetscFunctionBegin;
351731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
35271b4ebd2SPeter Brune   PetscFunctionReturn(0);
35371b4ebd2SPeter Brune }
35471b4ebd2SPeter Brune 
35571b4ebd2SPeter Brune 
35671b4ebd2SPeter Brune #undef __FUNCT__
357f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
358f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch)
35971b4ebd2SPeter Brune {
36071b4ebd2SPeter Brune 
36171b4ebd2SPeter Brune   PetscErrorCode       ierr;
362f1c6b773SPeter Brune   SNESLineSearch_BT    *bt;
36371b4ebd2SPeter Brune   PetscFunctionBegin;
36471b4ebd2SPeter Brune 
365f1c6b773SPeter Brune   bt = (SNESLineSearch_BT*)linesearch->data;
36671b4ebd2SPeter Brune 
367f1c6b773SPeter Brune   ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr);
3687a35526eSPeter Brune   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr);
36971b4ebd2SPeter Brune 
37071b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
37171b4ebd2SPeter Brune   PetscFunctionReturn(0);
37271b4ebd2SPeter Brune }
37371b4ebd2SPeter Brune 
37471b4ebd2SPeter Brune 
37571b4ebd2SPeter Brune #undef __FUNCT__
376f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT"
377954494b2SJed Brown /*MC
378db609ea7SPeter Brune    SNESLINESEARCHBT - Backtracking line search.
379954494b2SJed Brown 
380db609ea7SPeter Brune    This line search finds the minimum of a polynomial fitting of the L2 norm of the
381db609ea7SPeter Brune    function. If this fit does not satisfy the conditions for progress, the interval shrinks
382db609ea7SPeter Brune    and the fit is reattempted at most max_it times or until lambda is below minlambda.
383954494b2SJed Brown 
384cd7522eaSPeter Brune    Options Database Keys:
385cd7522eaSPeter Brune +  -snes_linesearch_alpha<1e-4> - slope descent parameter
386db609ea7SPeter Brune .  -snes_linesearch_damping<1.0> - initial step length
387db609ea7SPeter Brune .  -snes_linesearch_max_it<40> - maximum number of shrinking step
388db609ea7SPeter Brune .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
389cd7522eaSPeter Brune -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
390cd7522eaSPeter Brune 
391954494b2SJed Brown    Level: advanced
392954494b2SJed Brown 
393db609ea7SPeter Brune    Notes:
394db609ea7SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
395db609ea7SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
396db609ea7SPeter Brune 
397f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping
398954494b2SJed Brown 
399f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
400954494b2SJed Brown M*/
401f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
40271b4ebd2SPeter Brune {
40371b4ebd2SPeter Brune 
404f1c6b773SPeter Brune   SNESLineSearch_BT  *bt;
40571b4ebd2SPeter Brune   PetscErrorCode ierr;
40671b4ebd2SPeter Brune 
40771b4ebd2SPeter Brune   PetscFunctionBegin;
408f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_BT;
409f1c6b773SPeter Brune   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
410f1c6b773SPeter Brune   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
41171b4ebd2SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
4127f1410a3SPeter Brune   linesearch->ops->view           = SNESLineSearchView_BT;
41371b4ebd2SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
41471b4ebd2SPeter Brune 
415f1c6b773SPeter Brune   ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr);
41671b4ebd2SPeter Brune   linesearch->data = (void *)bt;
41771b4ebd2SPeter Brune   linesearch->max_its = 40;
418b000cd8dSPeter Brune   linesearch->order = SNES_LINESEARCH_ORDER_CUBIC;
41971b4ebd2SPeter Brune   bt->alpha = 1e-4;
42071b4ebd2SPeter Brune 
42171b4ebd2SPeter Brune   PetscFunctionReturn(0);
42271b4ebd2SPeter Brune }
423