171b4ebd2SPeter Brune #include <private/linesearchimpl.h> 271b4ebd2SPeter Brune #include <private/snesimpl.h> 371b4ebd2SPeter Brune 4f1c6b773SPeter Brune typedef enum {SNES_LINESEARCH_BT_QUADRATIC, SNES_LINESEARCH_BT_CUBIC} SNESLineSearchBTOrder; 571b4ebd2SPeter Brune 671b4ebd2SPeter Brune typedef struct { 771b4ebd2SPeter Brune PetscReal alpha; /* sufficient decrease parameter */ 8f1c6b773SPeter Brune SNESLineSearchBTOrder order; 9f1c6b773SPeter Brune } SNESLineSearch_BT; 1071b4ebd2SPeter Brune 1171b4ebd2SPeter Brune #undef __FUNCT__ 12f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_BT" 13f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) 1471b4ebd2SPeter Brune { 1571b4ebd2SPeter Brune PetscBool changed_y, changed_w; 1671b4ebd2SPeter Brune PetscErrorCode ierr; 1771b4ebd2SPeter Brune Vec X, F, Y, W, G; 1871b4ebd2SPeter Brune SNES snes; 1971b4ebd2SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm, gnormprev; 2071b4ebd2SPeter Brune PetscReal lambda, lambdatemp, lambdaprev, minlambda, maxstep, rellength, initslope, alpha; 2171b4ebd2SPeter Brune PetscReal t1, t2, a, b, d, steptol; 2271b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 2371b4ebd2SPeter Brune PetscScalar cinitslope; 2471b4ebd2SPeter Brune #endif 2571b4ebd2SPeter Brune PetscBool domainerror; 2671b4ebd2SPeter Brune PetscViewer monitor; 2771b4ebd2SPeter Brune PetscInt max_its, count; 28f1c6b773SPeter Brune SNESLineSearch_BT *bt; 2971b4ebd2SPeter Brune Mat jac; 3071b4ebd2SPeter Brune 3171b4ebd2SPeter Brune 3271b4ebd2SPeter Brune PetscFunctionBegin; 3371b4ebd2SPeter Brune 34f1c6b773SPeter Brune ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr); 35f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 36f1c6b773SPeter Brune ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 37f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 38f1c6b773SPeter Brune ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 39f1c6b773SPeter Brune ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_its); 40f1c6b773SPeter Brune bt = (SNESLineSearch_BT *)linesearch->data; 4171b4ebd2SPeter Brune 4271b4ebd2SPeter Brune alpha = bt->alpha; 4371b4ebd2SPeter Brune 4471b4ebd2SPeter Brune ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 4571b4ebd2SPeter Brune if (!jac) { 46f1c6b773SPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 4771b4ebd2SPeter Brune } 4871b4ebd2SPeter Brune /* precheck */ 49f1c6b773SPeter Brune ierr = SNESLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 50f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 5171b4ebd2SPeter Brune 5271b4ebd2SPeter Brune ierr = VecNorm(Y, NORM_2, &ynorm);CHKERRQ(ierr); 5371b4ebd2SPeter Brune if (ynorm == 0.0) { 5471b4ebd2SPeter Brune if (monitor) { 55ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 5671b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 57ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 5871b4ebd2SPeter Brune } 5971b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 6071b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 61f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 6271b4ebd2SPeter Brune PetscFunctionReturn(0); 6371b4ebd2SPeter Brune } 6471b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 6571b4ebd2SPeter Brune if (monitor) { 66ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 6771b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);CHKERRQ(ierr); 68ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 6971b4ebd2SPeter Brune } 7071b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 7171b4ebd2SPeter Brune ynorm = maxstep; 7271b4ebd2SPeter Brune } 7371b4ebd2SPeter Brune ierr = VecMaxPointwiseDivide(Y,X,&rellength);CHKERRQ(ierr); 7471b4ebd2SPeter Brune minlambda = steptol/rellength; 7571b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 7671b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 7771b4ebd2SPeter Brune ierr = VecDot(F,W,&cinitslope);CHKERRQ(ierr); 7871b4ebd2SPeter Brune initslope = PetscRealPart(cinitslope); 7971b4ebd2SPeter Brune #else 8071b4ebd2SPeter Brune ierr = VecDot(F,W,&initslope);CHKERRQ(ierr); 8171b4ebd2SPeter Brune #endif 8271b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 8371b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 8471b4ebd2SPeter Brune 85c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 869bd66eb0SPeter Brune if (linesearch->ops->viproject) { 879bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 889bd66eb0SPeter Brune } 8971b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 9071b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 9171b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 92f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 9371b4ebd2SPeter Brune PetscFunctionReturn(0); 9471b4ebd2SPeter Brune } 9571b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 96c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 97c427506bSPeter Brune if (domainerror) { 98f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 9971b4ebd2SPeter Brune PetscFunctionReturn(0); 10071b4ebd2SPeter Brune } 1019bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1029bd66eb0SPeter Brune gnorm = fnorm; 1039bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1049bd66eb0SPeter Brune } else { 10571b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1069bd66eb0SPeter Brune } 1079bd66eb0SPeter Brune 10871b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 10971b4ebd2SPeter Brune ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 11071b4ebd2SPeter Brune if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */ 11171b4ebd2SPeter Brune if (monitor) { 112ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 11371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 114ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 11571b4ebd2SPeter Brune } 11671b4ebd2SPeter Brune } else { 11771b4ebd2SPeter Brune /* Fit points with quadratic */ 11871b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope); 11971b4ebd2SPeter Brune lambdaprev = lambda; 12071b4ebd2SPeter Brune gnormprev = gnorm; 12171b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 12271b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12371b4ebd2SPeter Brune else lambda = lambdatemp; 12471b4ebd2SPeter Brune 12571b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1269bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1279bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1289bd66eb0SPeter Brune } 12971b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 13071b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 13171b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 132f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 13371b4ebd2SPeter Brune PetscFunctionReturn(0); 13471b4ebd2SPeter Brune } 13571b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 1366188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1376188f407SPeter Brune if (domainerror) { 13871b4ebd2SPeter Brune PetscFunctionReturn(0); 13971b4ebd2SPeter Brune } 1409bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1419bd66eb0SPeter Brune gnorm = fnorm; 1429bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1439bd66eb0SPeter Brune } else { 14471b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1459bd66eb0SPeter Brune } 14671b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14771b4ebd2SPeter Brune if (monitor) { 148ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 14971b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr); 150ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 15171b4ebd2SPeter Brune } 15271b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 15371b4ebd2SPeter Brune if (monitor) { 154ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 15571b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 156ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 15771b4ebd2SPeter Brune } 15871b4ebd2SPeter Brune } else { 15971b4ebd2SPeter Brune /* Fit points with cubic */ 16071b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 16171b4ebd2SPeter Brune if (lambda <= minlambda) { 16271b4ebd2SPeter Brune if (monitor) { 163ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16471b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 16571b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 16671b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 16771b4ebd2SPeter Brune fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr); 168ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16971b4ebd2SPeter Brune } 170f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 17171b4ebd2SPeter Brune PetscFunctionReturn(0); 17271b4ebd2SPeter Brune } 173f1c6b773SPeter Brune if (bt->order == SNES_LINESEARCH_BT_CUBIC) { 17471b4ebd2SPeter Brune t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope; 17571b4ebd2SPeter Brune t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 17671b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 17771b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 17871b4ebd2SPeter Brune d = b*b - 3*a*initslope; 17971b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 18071b4ebd2SPeter Brune if (a == 0.0) { 18171b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 18271b4ebd2SPeter Brune } else { 18371b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 18471b4ebd2SPeter Brune } 185f1c6b773SPeter Brune } else if (bt->order == SNES_LINESEARCH_BT_QUADRATIC) { 18671b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope); 18771b4ebd2SPeter Brune } 18871b4ebd2SPeter Brune lambdaprev = lambda; 18971b4ebd2SPeter Brune gnormprev = gnorm; 19071b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 19171b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 19271b4ebd2SPeter Brune else lambda = lambdatemp; 19371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 19471b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 19571b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 19671b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 19771b4ebd2SPeter Brune fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr); 198f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 19971b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 20071b4ebd2SPeter Brune PetscFunctionReturn(0); 20171b4ebd2SPeter Brune } 20271b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 203c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 204c427506bSPeter Brune if (domainerror) { 205f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 20671b4ebd2SPeter Brune PetscFunctionReturn(0); 20771b4ebd2SPeter Brune } 2089bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2099bd66eb0SPeter Brune gnorm = fnorm; 2109bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2119bd66eb0SPeter Brune } else { 21271b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2139bd66eb0SPeter Brune } 21471b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 21571b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 21671b4ebd2SPeter Brune if (monitor) { 217ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 218f1c6b773SPeter Brune if (bt->order == SNES_LINESEARCH_BT_CUBIC) { 21971b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 22071b4ebd2SPeter Brune } else { 22171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 22271b4ebd2SPeter Brune } 223ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 22471b4ebd2SPeter Brune } 22571b4ebd2SPeter Brune break; 22671b4ebd2SPeter Brune } else { 22771b4ebd2SPeter Brune if (monitor) { 228ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 229f1c6b773SPeter Brune if (bt->order == SNES_LINESEARCH_BT_CUBIC) { 23071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 23171b4ebd2SPeter Brune } else { 23271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 23371b4ebd2SPeter Brune } 234ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 23571b4ebd2SPeter Brune } 23671b4ebd2SPeter Brune } 23771b4ebd2SPeter Brune } 23871b4ebd2SPeter Brune } 23971b4ebd2SPeter Brune } 24071b4ebd2SPeter Brune 24171b4ebd2SPeter Brune /* postcheck */ 242f1c6b773SPeter Brune ierr = SNESLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 24371b4ebd2SPeter Brune if (changed_y) { 24471b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2459bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2469bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2479bd66eb0SPeter Brune } 24871b4ebd2SPeter Brune } 249c427506bSPeter Brune if (changed_y || changed_w) { /* recompute the function if the step has changed */ 250c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 25171b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 25271b4ebd2SPeter 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 { 2609bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2619bd66eb0SPeter Brune } 2629bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 263c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 2649bd66eb0SPeter Brune 265c427506bSPeter Brune } 26671b4ebd2SPeter Brune 26771b4ebd2SPeter Brune /* copy the solution over */ 26871b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 269c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 270c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 271f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 272f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 27371b4ebd2SPeter Brune PetscFunctionReturn(0); 27471b4ebd2SPeter Brune } 27571b4ebd2SPeter Brune 27671b4ebd2SPeter Brune #undef __FUNCT__ 277*7f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 278*7f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 279*7f1410a3SPeter Brune { 280*7f1410a3SPeter Brune PetscErrorCode ierr; 281*7f1410a3SPeter Brune PetscBool iascii; 282*7f1410a3SPeter Brune SNESLineSearch_BT *bt; 283*7f1410a3SPeter Brune PetscFunctionBegin; 284*7f1410a3SPeter Brune ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 285*7f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 286*7f1410a3SPeter Brune if (iascii) { 287*7f1410a3SPeter Brune if (bt->order == SNES_LINESEARCH_BT_CUBIC) { 288*7f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 289*7f1410a3SPeter Brune } else { 290*7f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 291*7f1410a3SPeter Brune } 292*7f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%G\n", bt->alpha);CHKERRQ(ierr); 293*7f1410a3SPeter Brune } 294*7f1410a3SPeter Brune PetscFunctionReturn(0); 295*7f1410a3SPeter Brune } 296*7f1410a3SPeter Brune 297*7f1410a3SPeter Brune 298*7f1410a3SPeter Brune #undef __FUNCT__ 299f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 300f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 301731c067cSJed Brown { 302731c067cSJed Brown PetscErrorCode ierr; 30371b4ebd2SPeter Brune 30471b4ebd2SPeter Brune PetscFunctionBegin; 305731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 30671b4ebd2SPeter Brune PetscFunctionReturn(0); 30771b4ebd2SPeter Brune } 30871b4ebd2SPeter Brune 30971b4ebd2SPeter Brune 31071b4ebd2SPeter Brune #undef __FUNCT__ 311f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 312f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 31371b4ebd2SPeter Brune { 31471b4ebd2SPeter Brune 31571b4ebd2SPeter Brune PetscErrorCode ierr; 316f1c6b773SPeter Brune SNESLineSearch_BT *bt; 31771b4ebd2SPeter Brune const char *orders[] = {"quadratic", "cubic"}; 31871b4ebd2SPeter Brune PetscInt indx = 0; 31971b4ebd2SPeter Brune PetscBool flg; 32071b4ebd2SPeter Brune PetscFunctionBegin; 32171b4ebd2SPeter Brune 322f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 32371b4ebd2SPeter Brune 324f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 3257a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 3267a35526eSPeter Brune ierr = PetscOptionsEList("-snes_linesearch_order", "Order of approximation", "SNESLineSearchBT", orders,2,"quadratic",&indx,&flg);CHKERRQ(ierr); 32771b4ebd2SPeter Brune if (flg) { 32871b4ebd2SPeter Brune switch (indx) { 329f1c6b773SPeter Brune case 0: bt->order = SNES_LINESEARCH_BT_QUADRATIC; 33071b4ebd2SPeter Brune break; 331f1c6b773SPeter Brune case 1: bt->order = SNES_LINESEARCH_BT_CUBIC; 33271b4ebd2SPeter Brune break; 33371b4ebd2SPeter Brune } 33471b4ebd2SPeter Brune } 33571b4ebd2SPeter Brune 33671b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 33771b4ebd2SPeter Brune PetscFunctionReturn(0); 33871b4ebd2SPeter Brune } 33971b4ebd2SPeter Brune 34071b4ebd2SPeter Brune 34171b4ebd2SPeter Brune #undef __FUNCT__ 342f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 343954494b2SJed Brown /*MC 344f1c6b773SPeter Brune SNES_LINESEARCH_BT - Backtracking line searches. 345954494b2SJed Brown 346954494b2SJed Brown These linesearches try a polynomial fit for the L2 norm of the error 347954494b2SJed Brown using the gradient. Failing that, they step back and try again. 348954494b2SJed Brown 349cd7522eaSPeter Brune Options Database Keys: 350cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 351cd7522eaSPeter Brune . -snes_linesearch_damping<1.0> - full steplength 352cd7522eaSPeter Brune - -snes_linesearch_order<cubic, quadratic> - order of the approximation 353cd7522eaSPeter Brune 354954494b2SJed Brown Level: advanced 355954494b2SJed Brown 356f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 357954494b2SJed Brown 358f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 359954494b2SJed Brown M*/ 360f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 36171b4ebd2SPeter Brune { 36271b4ebd2SPeter Brune 363f1c6b773SPeter Brune SNESLineSearch_BT *bt; 36471b4ebd2SPeter Brune PetscErrorCode ierr; 36571b4ebd2SPeter Brune 36671b4ebd2SPeter Brune PetscFunctionBegin; 367f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 368f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 369f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 37071b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 371*7f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 37271b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 37371b4ebd2SPeter Brune 374f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 37571b4ebd2SPeter Brune linesearch->data = (void *)bt; 37671b4ebd2SPeter Brune linesearch->max_its = 40; 377f1c6b773SPeter Brune bt->order = SNES_LINESEARCH_BT_CUBIC; 37871b4ebd2SPeter Brune bt->alpha = 1e-4; 37971b4ebd2SPeter Brune 38071b4ebd2SPeter Brune PetscFunctionReturn(0); 38171b4ebd2SPeter Brune } 382