1*b45d2f2cSJed Brown #include <petsc-private/linesearchimpl.h> 2*b45d2f2cSJed 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__ 9f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_BT" 10f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) 1171b4ebd2SPeter Brune { 1271b4ebd2SPeter Brune PetscBool changed_y, changed_w; 1371b4ebd2SPeter Brune PetscErrorCode ierr; 1471b4ebd2SPeter Brune Vec X, F, Y, W, G; 1571b4ebd2SPeter Brune SNES snes; 1671b4ebd2SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm, gnormprev; 1771b4ebd2SPeter Brune PetscReal lambda, lambdatemp, lambdaprev, minlambda, maxstep, rellength, initslope, alpha; 1871b4ebd2SPeter Brune PetscReal t1, t2, a, b, d, steptol; 1971b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 2071b4ebd2SPeter Brune PetscScalar cinitslope; 2171b4ebd2SPeter Brune #endif 2271b4ebd2SPeter Brune PetscBool domainerror; 2371b4ebd2SPeter Brune PetscViewer monitor; 2471b4ebd2SPeter Brune PetscInt max_its, count; 25f1c6b773SPeter Brune SNESLineSearch_BT *bt; 2671b4ebd2SPeter Brune Mat jac; 2771b4ebd2SPeter Brune 2871b4ebd2SPeter Brune 2971b4ebd2SPeter Brune PetscFunctionBegin; 3071b4ebd2SPeter Brune 31f1c6b773SPeter Brune ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr); 32f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 33f1c6b773SPeter Brune ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 34f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 35f1c6b773SPeter Brune ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 36f1c6b773SPeter Brune ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_its); 37f1c6b773SPeter Brune bt = (SNESLineSearch_BT *)linesearch->data; 3871b4ebd2SPeter Brune 3971b4ebd2SPeter Brune alpha = bt->alpha; 4071b4ebd2SPeter Brune 4171b4ebd2SPeter Brune ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 4271b4ebd2SPeter Brune if (!jac) { 43f1c6b773SPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 4471b4ebd2SPeter Brune } 4571b4ebd2SPeter Brune /* precheck */ 46f1c6b773SPeter Brune ierr = SNESLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 47f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 4871b4ebd2SPeter Brune 4971b4ebd2SPeter Brune ierr = VecNorm(Y, NORM_2, &ynorm);CHKERRQ(ierr); 5071b4ebd2SPeter Brune if (ynorm == 0.0) { 5171b4ebd2SPeter Brune if (monitor) { 52ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 5371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 54ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 5571b4ebd2SPeter Brune } 5671b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 5771b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 58f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 5971b4ebd2SPeter Brune PetscFunctionReturn(0); 6071b4ebd2SPeter Brune } 6171b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 6271b4ebd2SPeter Brune if (monitor) { 63ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 6471b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);CHKERRQ(ierr); 65ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 6671b4ebd2SPeter Brune } 6771b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 6871b4ebd2SPeter Brune ynorm = maxstep; 6971b4ebd2SPeter Brune } 7071b4ebd2SPeter Brune ierr = VecMaxPointwiseDivide(Y,X,&rellength);CHKERRQ(ierr); 7171b4ebd2SPeter Brune minlambda = steptol/rellength; 7271b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 7371b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 7471b4ebd2SPeter Brune ierr = VecDot(F,W,&cinitslope);CHKERRQ(ierr); 7571b4ebd2SPeter Brune initslope = PetscRealPart(cinitslope); 7671b4ebd2SPeter Brune #else 7771b4ebd2SPeter Brune ierr = VecDot(F,W,&initslope);CHKERRQ(ierr); 7871b4ebd2SPeter Brune #endif 7971b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 8071b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 8171b4ebd2SPeter Brune 82c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 839bd66eb0SPeter Brune if (linesearch->ops->viproject) { 849bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 859bd66eb0SPeter Brune } 8671b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 8771b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 8871b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 89f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 9071b4ebd2SPeter Brune PetscFunctionReturn(0); 9171b4ebd2SPeter Brune } 9271b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 93c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 94c427506bSPeter Brune if (domainerror) { 95f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 9671b4ebd2SPeter Brune PetscFunctionReturn(0); 9771b4ebd2SPeter Brune } 989bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 999bd66eb0SPeter Brune gnorm = fnorm; 1009bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1019bd66eb0SPeter Brune } else { 10271b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1039bd66eb0SPeter Brune } 1049bd66eb0SPeter Brune 10571b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 10671b4ebd2SPeter Brune ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 10771b4ebd2SPeter Brune if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */ 10871b4ebd2SPeter Brune if (monitor) { 109ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 11071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 111ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 11271b4ebd2SPeter Brune } 11371b4ebd2SPeter Brune } else { 11471b4ebd2SPeter Brune /* Fit points with quadratic */ 11571b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope); 11671b4ebd2SPeter Brune lambdaprev = lambda; 11771b4ebd2SPeter Brune gnormprev = gnorm; 11871b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 11971b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12071b4ebd2SPeter Brune else lambda = lambdatemp; 12171b4ebd2SPeter Brune 12271b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1239bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1249bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1259bd66eb0SPeter Brune } 12671b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 12771b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 12871b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 129f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 13071b4ebd2SPeter Brune PetscFunctionReturn(0); 13171b4ebd2SPeter Brune } 13271b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 1336188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1346188f407SPeter Brune if (domainerror) { 13571b4ebd2SPeter Brune PetscFunctionReturn(0); 13671b4ebd2SPeter Brune } 1379bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1389bd66eb0SPeter Brune gnorm = fnorm; 1399bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1409bd66eb0SPeter Brune } else { 14171b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1429bd66eb0SPeter Brune } 14371b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14471b4ebd2SPeter Brune if (monitor) { 145ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 14671b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr); 147ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 14871b4ebd2SPeter Brune } 14971b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 15071b4ebd2SPeter Brune if (monitor) { 151ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 15271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 153ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 15471b4ebd2SPeter Brune } 15571b4ebd2SPeter Brune } else { 15671b4ebd2SPeter Brune /* Fit points with cubic */ 15771b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 15871b4ebd2SPeter Brune if (lambda <= minlambda) { 15971b4ebd2SPeter Brune if (monitor) { 160ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 16271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 16371b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 16471b4ebd2SPeter Brune fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr); 165ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16671b4ebd2SPeter Brune } 167f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 16871b4ebd2SPeter Brune PetscFunctionReturn(0); 16971b4ebd2SPeter Brune } 17059405d5eSPeter Brune if (linesearch->order == SNES_LINESEARCH_CUBIC) { 17171b4ebd2SPeter Brune t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope; 17271b4ebd2SPeter Brune t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 17371b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 17471b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 17571b4ebd2SPeter Brune d = b*b - 3*a*initslope; 17671b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 17771b4ebd2SPeter Brune if (a == 0.0) { 17871b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 17971b4ebd2SPeter Brune } else { 18071b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 18171b4ebd2SPeter Brune } 18259405d5eSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_QUADRATIC) { 18371b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope); 18459405d5eSPeter Brune } else { 18559405d5eSPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt"); 18671b4ebd2SPeter Brune } 18771b4ebd2SPeter Brune lambdaprev = lambda; 18871b4ebd2SPeter Brune gnormprev = gnorm; 18971b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 19071b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 19171b4ebd2SPeter Brune else lambda = lambdatemp; 19271b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 19371b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 19471b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 19571b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 19671b4ebd2SPeter Brune fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr); 197f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 19871b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 19971b4ebd2SPeter Brune PetscFunctionReturn(0); 20071b4ebd2SPeter Brune } 20171b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 202c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 203c427506bSPeter Brune if (domainerror) { 204f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 20571b4ebd2SPeter Brune PetscFunctionReturn(0); 20671b4ebd2SPeter Brune } 2079bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2089bd66eb0SPeter Brune gnorm = fnorm; 2099bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2109bd66eb0SPeter Brune } else { 21171b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2129bd66eb0SPeter Brune } 21371b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 21471b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 21571b4ebd2SPeter Brune if (monitor) { 216ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 21759405d5eSPeter Brune if (linesearch->order == SNES_LINESEARCH_CUBIC) { 21871b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 21971b4ebd2SPeter Brune } else { 22071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 22171b4ebd2SPeter Brune } 222ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 22371b4ebd2SPeter Brune } 22471b4ebd2SPeter Brune break; 22571b4ebd2SPeter Brune } else { 22671b4ebd2SPeter Brune if (monitor) { 227ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 22859405d5eSPeter Brune if (linesearch->order == SNES_LINESEARCH_CUBIC) { 22971b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 23071b4ebd2SPeter Brune } else { 23171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 23271b4ebd2SPeter Brune } 233ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 23471b4ebd2SPeter Brune } 23571b4ebd2SPeter Brune } 23671b4ebd2SPeter Brune } 23771b4ebd2SPeter Brune } 23871b4ebd2SPeter Brune } 23971b4ebd2SPeter Brune 24071b4ebd2SPeter Brune /* postcheck */ 241f1c6b773SPeter Brune ierr = SNESLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 24271b4ebd2SPeter Brune if (changed_y) { 24371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2449bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2459bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2469bd66eb0SPeter Brune } 24771b4ebd2SPeter Brune } 248c427506bSPeter Brune if (changed_y || changed_w) { /* recompute the function if the step has changed */ 249c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 25071b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 25171b4ebd2SPeter Brune if (domainerror) { 252f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 25371b4ebd2SPeter Brune PetscFunctionReturn(0); 25471b4ebd2SPeter Brune } 2559bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2569bd66eb0SPeter Brune gnorm = fnorm; 2579bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2589bd66eb0SPeter Brune } else { 2599bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2609bd66eb0SPeter Brune } 2619bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 262c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 2639bd66eb0SPeter Brune 264c427506bSPeter Brune } 26571b4ebd2SPeter Brune 26671b4ebd2SPeter Brune /* copy the solution over */ 26771b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 268c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 269c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 270f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 271f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 27271b4ebd2SPeter Brune PetscFunctionReturn(0); 27371b4ebd2SPeter Brune } 27471b4ebd2SPeter Brune 27571b4ebd2SPeter Brune #undef __FUNCT__ 2767f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 2777f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 2787f1410a3SPeter Brune { 2797f1410a3SPeter Brune PetscErrorCode ierr; 2807f1410a3SPeter Brune PetscBool iascii; 2817f1410a3SPeter Brune SNESLineSearch_BT *bt; 2827f1410a3SPeter Brune PetscFunctionBegin; 2837f1410a3SPeter Brune ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 2847f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 2857f1410a3SPeter Brune if (iascii) { 28659405d5eSPeter Brune if (linesearch->order == SNES_LINESEARCH_CUBIC) { 2877f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 28859405d5eSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_QUADRATIC) { 2897f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 2907f1410a3SPeter Brune } 2917f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%G\n", bt->alpha);CHKERRQ(ierr); 2927f1410a3SPeter Brune } 2937f1410a3SPeter Brune PetscFunctionReturn(0); 2947f1410a3SPeter Brune } 2957f1410a3SPeter Brune 2967f1410a3SPeter Brune 2977f1410a3SPeter Brune #undef __FUNCT__ 298f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 299f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 300731c067cSJed Brown { 301731c067cSJed Brown PetscErrorCode ierr; 30271b4ebd2SPeter Brune 30371b4ebd2SPeter Brune PetscFunctionBegin; 304731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 30571b4ebd2SPeter Brune PetscFunctionReturn(0); 30671b4ebd2SPeter Brune } 30771b4ebd2SPeter Brune 30871b4ebd2SPeter Brune 30971b4ebd2SPeter Brune #undef __FUNCT__ 310f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 311f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 31271b4ebd2SPeter Brune { 31371b4ebd2SPeter Brune 31471b4ebd2SPeter Brune PetscErrorCode ierr; 315f1c6b773SPeter Brune SNESLineSearch_BT *bt; 31671b4ebd2SPeter Brune PetscFunctionBegin; 31771b4ebd2SPeter Brune 318f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 31971b4ebd2SPeter Brune 320f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 3217a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 32271b4ebd2SPeter Brune 32371b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 32471b4ebd2SPeter Brune PetscFunctionReturn(0); 32571b4ebd2SPeter Brune } 32671b4ebd2SPeter Brune 32771b4ebd2SPeter Brune 32871b4ebd2SPeter Brune #undef __FUNCT__ 329f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 330954494b2SJed Brown /*MC 331f1c6b773SPeter Brune SNES_LINESEARCH_BT - Backtracking line searches. 332954494b2SJed Brown 333954494b2SJed Brown These linesearches try a polynomial fit for the L2 norm of the error 334954494b2SJed Brown using the gradient. Failing that, they step back and try again. 335954494b2SJed Brown 336cd7522eaSPeter Brune Options Database Keys: 337cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 338cd7522eaSPeter Brune . -snes_linesearch_damping<1.0> - full steplength 339cd7522eaSPeter Brune - -snes_linesearch_order<cubic, quadratic> - order of the approximation 340cd7522eaSPeter Brune 341954494b2SJed Brown Level: advanced 342954494b2SJed Brown 343f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 344954494b2SJed Brown 345f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 346954494b2SJed Brown M*/ 347f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 34871b4ebd2SPeter Brune { 34971b4ebd2SPeter Brune 350f1c6b773SPeter Brune SNESLineSearch_BT *bt; 35171b4ebd2SPeter Brune PetscErrorCode ierr; 35271b4ebd2SPeter Brune 35371b4ebd2SPeter Brune PetscFunctionBegin; 354f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 355f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 356f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 35771b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 3587f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 35971b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 36071b4ebd2SPeter Brune 361f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 36271b4ebd2SPeter Brune linesearch->data = (void *)bt; 36371b4ebd2SPeter Brune linesearch->max_its = 40; 36459405d5eSPeter Brune linesearch->order = SNES_LINESEARCH_CUBIC; 36571b4ebd2SPeter Brune bt->alpha = 1e-4; 36671b4ebd2SPeter Brune 36771b4ebd2SPeter Brune PetscFunctionReturn(0); 36871b4ebd2SPeter Brune } 369