xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision ea5d4fccf296dd2bbd0f9c3a3343651cb1066da7)
171b4ebd2SPeter Brune #include <private/linesearchimpl.h>
271b4ebd2SPeter Brune #include <private/snesimpl.h>
371b4ebd2SPeter Brune 
46188f407SPeter Brune typedef enum {PETSCLINESEARCH_BT_QUADRATIC, PETSCLINESEARCH_BT_CUBIC} PetscLineSearchBTOrder;
571b4ebd2SPeter Brune 
671b4ebd2SPeter Brune typedef struct {
771b4ebd2SPeter Brune   PetscReal        alpha; /* sufficient decrease parameter */
86188f407SPeter Brune   PetscLineSearchBTOrder order;
96188f407SPeter Brune } PetscLineSearch_BT;
1071b4ebd2SPeter Brune 
1171b4ebd2SPeter Brune /*MC
126188f407SPeter Brune    PetscLineSearchBT - Backtracking line searches.
1371b4ebd2SPeter Brune 
1471b4ebd2SPeter Brune    These linesearches try a polynomial fit for the L2 norm of the error
1571b4ebd2SPeter Brune    using the gradient.  Failing that, they step back and try again.
1671b4ebd2SPeter Brune 
1771b4ebd2SPeter Brune    Level: advanced
1871b4ebd2SPeter Brune 
196188f407SPeter Brune .keywords: SNES, PetscLineSearch, damping
2071b4ebd2SPeter Brune 
216188f407SPeter Brune .seealso: PetscLineSearchCreate(), PetscLineSearchSetType()
2271b4ebd2SPeter Brune M*/
2371b4ebd2SPeter Brune 
2471b4ebd2SPeter Brune #undef __FUNCT__
256188f407SPeter Brune #define __FUNCT__ "PetscLineSearchApply_BT"
2671b4ebd2SPeter Brune 
276188f407SPeter Brune PetscErrorCode  PetscLineSearchApply_BT(PetscLineSearch linesearch)
2871b4ebd2SPeter Brune {
2971b4ebd2SPeter Brune   PetscBool      changed_y, changed_w;
3071b4ebd2SPeter Brune   PetscErrorCode ierr;
3171b4ebd2SPeter Brune   Vec            X, F, Y, W, G;
3271b4ebd2SPeter Brune   SNES           snes;
3371b4ebd2SPeter Brune   PetscReal      fnorm, xnorm, ynorm, gnorm, gnormprev;
3471b4ebd2SPeter Brune   PetscReal      lambda, lambdatemp, lambdaprev, minlambda, maxstep, rellength, initslope, alpha;
3571b4ebd2SPeter Brune   PetscReal      t1, t2, a, b, d, steptol;
3671b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX)
3771b4ebd2SPeter Brune   PetscScalar    cinitslope;
3871b4ebd2SPeter Brune #endif
3971b4ebd2SPeter Brune   PetscBool      domainerror;
4071b4ebd2SPeter Brune   PetscViewer    monitor;
4171b4ebd2SPeter Brune   PetscInt       max_its, count;
426188f407SPeter Brune   PetscLineSearch_BT  *bt;
4371b4ebd2SPeter Brune   Mat            jac;
4471b4ebd2SPeter Brune 
4571b4ebd2SPeter Brune 
4671b4ebd2SPeter Brune   PetscFunctionBegin;
4771b4ebd2SPeter Brune 
486188f407SPeter Brune   ierr = PetscLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
496188f407SPeter Brune   ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
506188f407SPeter Brune   ierr = PetscLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
516188f407SPeter Brune   ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
526188f407SPeter Brune   ierr = PetscLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr);
536188f407SPeter Brune   ierr = PetscLineSearchGetTolerances(linesearch, &steptol, &maxstep, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_its);
546188f407SPeter Brune   bt = (PetscLineSearch_BT *)linesearch->data;
5571b4ebd2SPeter Brune 
5671b4ebd2SPeter Brune   alpha = bt->alpha;
5771b4ebd2SPeter Brune 
5871b4ebd2SPeter Brune   ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
5971b4ebd2SPeter Brune   if (!jac) {
606188f407SPeter Brune     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "PetscLineSearchBT requires a Jacobian matrix");
6171b4ebd2SPeter Brune   }
6271b4ebd2SPeter Brune   /* precheck */
636188f407SPeter Brune   ierr = PetscLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr);
646188f407SPeter Brune   ierr = PetscLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
6571b4ebd2SPeter Brune 
6671b4ebd2SPeter Brune   ierr = VecNorm(Y, NORM_2, &ynorm);CHKERRQ(ierr);
6771b4ebd2SPeter Brune   if (ynorm == 0.0) {
6871b4ebd2SPeter Brune     if (monitor) {
69*ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
7071b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
71*ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
7271b4ebd2SPeter Brune     }
7371b4ebd2SPeter Brune     ierr   = VecCopy(X,W);CHKERRQ(ierr);
7471b4ebd2SPeter Brune     ierr   = VecCopy(F,G);CHKERRQ(ierr);
756188f407SPeter Brune     ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
7671b4ebd2SPeter Brune     PetscFunctionReturn(0);
7771b4ebd2SPeter Brune   }
7871b4ebd2SPeter Brune   if (ynorm > maxstep) {	/* Step too big, so scale back */
7971b4ebd2SPeter Brune     if (monitor) {
80*ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
8171b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);CHKERRQ(ierr);
82*ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
8371b4ebd2SPeter Brune     }
8471b4ebd2SPeter Brune     ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
8571b4ebd2SPeter Brune     ynorm = maxstep;
8671b4ebd2SPeter Brune   }
8771b4ebd2SPeter Brune   ierr      = VecMaxPointwiseDivide(Y,X,&rellength);CHKERRQ(ierr);
8871b4ebd2SPeter Brune   minlambda = steptol/rellength;
8971b4ebd2SPeter Brune   ierr      = MatMult(jac,Y,W);CHKERRQ(ierr);
9071b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX)
9171b4ebd2SPeter Brune   ierr      = VecDot(F,W,&cinitslope);CHKERRQ(ierr);
9271b4ebd2SPeter Brune   initslope = PetscRealPart(cinitslope);
9371b4ebd2SPeter Brune #else
9471b4ebd2SPeter Brune   ierr      = VecDot(F,W,&initslope);CHKERRQ(ierr);
9571b4ebd2SPeter Brune #endif
9671b4ebd2SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
9771b4ebd2SPeter Brune   if (initslope == 0.0) initslope = -1.0;
9871b4ebd2SPeter Brune 
99c427506bSPeter Brune   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1009bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1019bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1029bd66eb0SPeter Brune   }
10371b4ebd2SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
10471b4ebd2SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
10571b4ebd2SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1066188f407SPeter Brune     ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
10771b4ebd2SPeter Brune     PetscFunctionReturn(0);
10871b4ebd2SPeter Brune   }
10971b4ebd2SPeter Brune   ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
110c427506bSPeter Brune   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
111c427506bSPeter Brune   if (domainerror) {
1126188f407SPeter Brune     ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
11371b4ebd2SPeter Brune     PetscFunctionReturn(0);
11471b4ebd2SPeter Brune   }
1159bd66eb0SPeter Brune   if (linesearch->ops->vinorm) {
1169bd66eb0SPeter Brune     gnorm = fnorm;
1179bd66eb0SPeter Brune     ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1189bd66eb0SPeter Brune   } else {
11971b4ebd2SPeter Brune     ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1209bd66eb0SPeter Brune   }
1219bd66eb0SPeter Brune 
12271b4ebd2SPeter Brune   if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
12371b4ebd2SPeter Brune   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr);
12471b4ebd2SPeter Brune   if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */
12571b4ebd2SPeter Brune     if (monitor) {
126*ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
12771b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr);
128*ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
12971b4ebd2SPeter Brune     }
13071b4ebd2SPeter Brune   } else {
13171b4ebd2SPeter Brune     /* Fit points with quadratic */
13271b4ebd2SPeter Brune     lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope);
13371b4ebd2SPeter Brune     lambdaprev = lambda;
13471b4ebd2SPeter Brune     gnormprev  = gnorm;
13571b4ebd2SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
13671b4ebd2SPeter Brune     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
13771b4ebd2SPeter Brune     else                         lambda = lambdatemp;
13871b4ebd2SPeter Brune 
13971b4ebd2SPeter Brune     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
1409bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1419bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1429bd66eb0SPeter Brune     }
14371b4ebd2SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
14471b4ebd2SPeter Brune       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
14571b4ebd2SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
1466188f407SPeter Brune       ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
14771b4ebd2SPeter Brune       PetscFunctionReturn(0);
14871b4ebd2SPeter Brune     }
14971b4ebd2SPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
1506188f407SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
1516188f407SPeter Brune     if (domainerror) {
15271b4ebd2SPeter Brune       PetscFunctionReturn(0);
15371b4ebd2SPeter Brune     }
1549bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
1559bd66eb0SPeter Brune       gnorm = fnorm;
1569bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
1579bd66eb0SPeter Brune     } else {
15871b4ebd2SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
1599bd66eb0SPeter Brune     }
16071b4ebd2SPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
16171b4ebd2SPeter Brune     if (monitor) {
162*ea5d4fccSPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
16371b4ebd2SPeter Brune       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr);
164*ea5d4fccSPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
16571b4ebd2SPeter Brune     }
16671b4ebd2SPeter Brune     if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
16771b4ebd2SPeter Brune       if (monitor) {
168*ea5d4fccSPeter Brune         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
16971b4ebd2SPeter Brune         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
170*ea5d4fccSPeter Brune         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
17171b4ebd2SPeter Brune       }
17271b4ebd2SPeter Brune     } else {
17371b4ebd2SPeter Brune       /* Fit points with cubic */
17471b4ebd2SPeter Brune       for (count = 0; count < max_its; count++) {
17571b4ebd2SPeter Brune         if (lambda <= minlambda) {
17671b4ebd2SPeter Brune           if (monitor) {
177*ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
17871b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
17971b4ebd2SPeter Brune             ierr = PetscViewerASCIIPrintf(monitor,
18071b4ebd2SPeter Brune                                           "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
18171b4ebd2SPeter Brune                                           fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr);
182*ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
18371b4ebd2SPeter Brune           }
1846188f407SPeter Brune           ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
18571b4ebd2SPeter Brune           PetscFunctionReturn(0);
18671b4ebd2SPeter Brune         }
1876188f407SPeter Brune         if (bt->order == PETSCLINESEARCH_BT_CUBIC) {
18871b4ebd2SPeter Brune           t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope;
18971b4ebd2SPeter Brune           t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
19071b4ebd2SPeter Brune           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
19171b4ebd2SPeter Brune           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
19271b4ebd2SPeter Brune           d  = b*b - 3*a*initslope;
19371b4ebd2SPeter Brune           if (d < 0.0) d = 0.0;
19471b4ebd2SPeter Brune           if (a == 0.0) {
19571b4ebd2SPeter Brune             lambdatemp = -initslope/(2.0*b);
19671b4ebd2SPeter Brune           } else {
19771b4ebd2SPeter Brune             lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
19871b4ebd2SPeter Brune           }
1996188f407SPeter Brune         } else if (bt->order == PETSCLINESEARCH_BT_QUADRATIC) {
20071b4ebd2SPeter Brune           lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope);
20171b4ebd2SPeter Brune         }
20271b4ebd2SPeter Brune           lambdaprev = lambda;
20371b4ebd2SPeter Brune           gnormprev  = gnorm;
20471b4ebd2SPeter Brune         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
20571b4ebd2SPeter Brune         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
20671b4ebd2SPeter Brune         else                         lambda     = lambdatemp;
20771b4ebd2SPeter Brune         ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
20871b4ebd2SPeter Brune         if (snes->nfuncs >= snes->max_funcs) {
20971b4ebd2SPeter Brune           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
21071b4ebd2SPeter Brune           ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
21171b4ebd2SPeter Brune                             fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr);
2126188f407SPeter Brune           ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
21371b4ebd2SPeter Brune           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
21471b4ebd2SPeter Brune           PetscFunctionReturn(0);
21571b4ebd2SPeter Brune         }
21671b4ebd2SPeter Brune         ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
217c427506bSPeter Brune         ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
218c427506bSPeter Brune         if (domainerror) {
2196188f407SPeter Brune           ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
22071b4ebd2SPeter Brune           PetscFunctionReturn(0);
22171b4ebd2SPeter Brune         }
2229bd66eb0SPeter Brune         if (linesearch->ops->vinorm) {
2239bd66eb0SPeter Brune           gnorm = fnorm;
2249bd66eb0SPeter Brune           ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2259bd66eb0SPeter Brune         } else {
22671b4ebd2SPeter Brune           ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2279bd66eb0SPeter Brune         }
22871b4ebd2SPeter Brune         if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
22971b4ebd2SPeter Brune         if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
23071b4ebd2SPeter Brune           if (monitor) {
231*ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
2326188f407SPeter Brune             if (bt->order == PETSCLINESEARCH_BT_CUBIC) {
23371b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
23471b4ebd2SPeter Brune             } else {
23571b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
23671b4ebd2SPeter Brune             }
237*ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
23871b4ebd2SPeter Brune           }
23971b4ebd2SPeter Brune           break;
24071b4ebd2SPeter Brune         } else {
24171b4ebd2SPeter Brune           if (monitor) {
242*ea5d4fccSPeter Brune             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
2436188f407SPeter Brune             if (bt->order == PETSCLINESEARCH_BT_CUBIC) {
24471b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
24571b4ebd2SPeter Brune             } else {
24671b4ebd2SPeter Brune               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr);
24771b4ebd2SPeter Brune             }
248*ea5d4fccSPeter Brune             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
24971b4ebd2SPeter Brune           }
25071b4ebd2SPeter Brune         }
25171b4ebd2SPeter Brune       }
25271b4ebd2SPeter Brune     }
25371b4ebd2SPeter Brune   }
25471b4ebd2SPeter Brune 
25571b4ebd2SPeter Brune   /* postcheck */
2566188f407SPeter Brune   ierr = PetscLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr);
25771b4ebd2SPeter Brune   if (changed_y) {
25871b4ebd2SPeter Brune     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
2599bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
2609bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
2619bd66eb0SPeter Brune     }
26271b4ebd2SPeter Brune   }
263c427506bSPeter Brune   if (changed_y || changed_w) { /* recompute the function if the step has changed */
264c427506bSPeter Brune     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
26571b4ebd2SPeter Brune     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
26671b4ebd2SPeter Brune     if (domainerror) {
2676188f407SPeter Brune       ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
26871b4ebd2SPeter Brune       PetscFunctionReturn(0);
26971b4ebd2SPeter Brune     }
2709bd66eb0SPeter Brune     if (linesearch->ops->vinorm) {
2719bd66eb0SPeter Brune       gnorm = fnorm;
2729bd66eb0SPeter Brune       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
2739bd66eb0SPeter Brune     } else {
2749bd66eb0SPeter Brune       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
2759bd66eb0SPeter Brune     }
2769bd66eb0SPeter Brune     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
277c427506bSPeter Brune     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
2789bd66eb0SPeter Brune 
279c427506bSPeter Brune   }
28071b4ebd2SPeter Brune 
28171b4ebd2SPeter Brune   /* copy the solution over */
28271b4ebd2SPeter Brune   ierr = VecCopy(W, X);CHKERRQ(ierr);
283c427506bSPeter Brune   ierr = VecCopy(G, F);CHKERRQ(ierr);
284c427506bSPeter Brune   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
2856188f407SPeter Brune   ierr = PetscLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
2866188f407SPeter Brune   ierr = PetscLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
28771b4ebd2SPeter Brune   PetscFunctionReturn(0);
28871b4ebd2SPeter Brune }
28971b4ebd2SPeter Brune 
29071b4ebd2SPeter Brune #undef __FUNCT__
2916188f407SPeter Brune #define __FUNCT__ "PetscLineSearchDestroy_BT"
292731c067cSJed Brown PetscErrorCode PetscLineSearchDestroy_BT(PetscLineSearch linesearch)
293731c067cSJed Brown {
294731c067cSJed Brown   PetscErrorCode ierr;
29571b4ebd2SPeter Brune 
29671b4ebd2SPeter Brune   PetscFunctionBegin;
297731c067cSJed Brown   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
29871b4ebd2SPeter Brune   PetscFunctionReturn(0);
29971b4ebd2SPeter Brune }
30071b4ebd2SPeter Brune 
30171b4ebd2SPeter Brune 
30271b4ebd2SPeter Brune #undef __FUNCT__
3036188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetFromOptions_BT"
3046188f407SPeter Brune static PetscErrorCode PetscLineSearchSetFromOptions_BT(PetscLineSearch linesearch)
30571b4ebd2SPeter Brune {
30671b4ebd2SPeter Brune 
30771b4ebd2SPeter Brune   PetscErrorCode ierr;
3086188f407SPeter Brune   PetscLineSearch_BT    *bt;
30971b4ebd2SPeter Brune   const char       *orders[] = {"quadratic", "cubic"};
31071b4ebd2SPeter Brune   PetscInt         indx = 0;
31171b4ebd2SPeter Brune   PetscBool        flg;
31271b4ebd2SPeter Brune   PetscFunctionBegin;
31371b4ebd2SPeter Brune 
3146188f407SPeter Brune   bt = (PetscLineSearch_BT*)linesearch->data;
31571b4ebd2SPeter Brune 
3166188f407SPeter Brune   ierr = PetscOptionsHead("PetscLineSearch BT options");CHKERRQ(ierr);
3176188f407SPeter Brune   ierr = PetscOptionsReal("-linesearch_bt_alpha",   "Descent tolerance",        "PetscLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr);
3186188f407SPeter Brune   ierr = PetscOptionsEList("-linesearch_bt_order",  "Order of approximation",   "PetscLineSearchBT", orders,2,"quadratic",&indx,&flg);CHKERRQ(ierr);
31971b4ebd2SPeter Brune   if (flg) {
32071b4ebd2SPeter Brune     switch (indx) {
3216188f407SPeter Brune     case 0: bt->order = PETSCLINESEARCH_BT_QUADRATIC;
32271b4ebd2SPeter Brune       break;
3236188f407SPeter Brune     case 1: bt->order = PETSCLINESEARCH_BT_CUBIC;
32471b4ebd2SPeter Brune       break;
32571b4ebd2SPeter Brune     }
32671b4ebd2SPeter Brune   }
32771b4ebd2SPeter Brune 
32871b4ebd2SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
32971b4ebd2SPeter Brune   PetscFunctionReturn(0);
33071b4ebd2SPeter Brune }
33171b4ebd2SPeter Brune 
33271b4ebd2SPeter Brune 
33371b4ebd2SPeter Brune EXTERN_C_BEGIN
33471b4ebd2SPeter Brune #undef __FUNCT__
3356188f407SPeter Brune #define __FUNCT__ "PetscLineSearchCreate_BT"
3366188f407SPeter Brune PetscErrorCode PetscLineSearchCreate_BT(PetscLineSearch linesearch)
33771b4ebd2SPeter Brune {
33871b4ebd2SPeter Brune 
3396188f407SPeter Brune   PetscLineSearch_BT  *bt;
34071b4ebd2SPeter Brune   PetscErrorCode ierr;
34171b4ebd2SPeter Brune 
34271b4ebd2SPeter Brune   PetscFunctionBegin;
3436188f407SPeter Brune   linesearch->ops->apply          = PetscLineSearchApply_BT;
3446188f407SPeter Brune   linesearch->ops->destroy        = PetscLineSearchDestroy_BT;
3456188f407SPeter Brune   linesearch->ops->setfromoptions = PetscLineSearchSetFromOptions_BT;
34671b4ebd2SPeter Brune   linesearch->ops->reset          = PETSC_NULL;
34771b4ebd2SPeter Brune   linesearch->ops->view           = PETSC_NULL;
34871b4ebd2SPeter Brune   linesearch->ops->setup          = PETSC_NULL;
34971b4ebd2SPeter Brune 
3506188f407SPeter Brune   ierr = PetscNewLog(linesearch, PetscLineSearch_BT, &bt);CHKERRQ(ierr);
35171b4ebd2SPeter Brune   linesearch->data = (void *)bt;
35271b4ebd2SPeter Brune   linesearch->max_its = 40;
3536188f407SPeter Brune   bt->order = PETSCLINESEARCH_BT_CUBIC;
35471b4ebd2SPeter Brune   bt->alpha = 1e-4;
35571b4ebd2SPeter Brune 
35671b4ebd2SPeter Brune   PetscFunctionReturn(0);
35771b4ebd2SPeter Brune }
35871b4ebd2SPeter Brune EXTERN_C_END
359