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; 246e111a19SKarl Rupp 252f4102e2SPeter Brune PetscFunctionBegin; 262f4102e2SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 272f4102e2SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 282f4102e2SPeter Brune bt->alpha = alpha; 292f4102e2SPeter Brune PetscFunctionReturn(0); 302f4102e2SPeter Brune } 312f4102e2SPeter Brune 322f4102e2SPeter Brune 332f4102e2SPeter Brune #undef __FUNCT__ 342f4102e2SPeter Brune #define __FUNCT__ "SNESLineSearchBTGetAlpha" 352f4102e2SPeter Brune /*@ 362f4102e2SPeter Brune SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant. 372f4102e2SPeter Brune 382f4102e2SPeter Brune Input Parameters: 392f4102e2SPeter Brune . linesearch - linesearch context 408e557f58SPeter Brune 418e557f58SPeter Brune Output Parameters: 422f4102e2SPeter Brune . alpha - The descent parameter 432f4102e2SPeter Brune 442f4102e2SPeter Brune Level: intermediate 452f4102e2SPeter Brune 461a4f838cSPeter Brune .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT 472f4102e2SPeter Brune @*/ 482f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha) 492f4102e2SPeter Brune { 502f4102e2SPeter Brune SNESLineSearch_BT *bt; 516e111a19SKarl Rupp 522f4102e2SPeter Brune PetscFunctionBegin; 532f4102e2SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 542f4102e2SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 552f4102e2SPeter Brune *alpha = bt->alpha; 562f4102e2SPeter Brune PetscFunctionReturn(0); 572f4102e2SPeter Brune } 582f4102e2SPeter Brune 592f4102e2SPeter Brune #undef __FUNCT__ 60f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_BT" 61f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) 6271b4ebd2SPeter Brune { 6371b4ebd2SPeter Brune PetscBool changed_y,changed_w; 6471b4ebd2SPeter Brune PetscErrorCode ierr; 6571b4ebd2SPeter Brune Vec X,F,Y,W,G; 6671b4ebd2SPeter Brune SNES snes; 678a430afbSPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm; 683b288129SPeter Brune PetscReal lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol; 69dfcd3828SPeter Brune PetscReal t1,t2,a,b,d; 708a430afbSPeter Brune PetscReal f; 718a430afbSPeter Brune PetscReal g,gprev; 7271b4ebd2SPeter Brune PetscBool domainerror; 7371b4ebd2SPeter Brune PetscViewer monitor; 7471b4ebd2SPeter Brune PetscInt max_its,count; 75f1c6b773SPeter Brune SNESLineSearch_BT *bt; 7671b4ebd2SPeter Brune Mat jac; 77bd4a8a71SBarry Smith PetscErrorCode (*objective)(SNES,Vec,PetscReal*,void*); 7871b4ebd2SPeter Brune 7971b4ebd2SPeter Brune PetscFunctionBegin; 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); 850298fd71SBarry Smith ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,NULL,NULL,NULL,&max_its);CHKERRQ(ierr); 860298fd71SBarry Smith ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr); 870298fd71SBarry Smith ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr); 88f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 8971b4ebd2SPeter Brune 9071b4ebd2SPeter Brune alpha = bt->alpha; 9171b4ebd2SPeter Brune 920298fd71SBarry Smith ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr); 938a430afbSPeter Brune 94ce94432eSBarry Smith if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 95c69d1a72SBarry Smith 9671b4ebd2SPeter Brune /* precheck */ 977b1df9c1SPeter Brune ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); 98f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 9971b4ebd2SPeter Brune 1003b288129SPeter Brune ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr); 1013b288129SPeter Brune ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 1023b288129SPeter Brune ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr); 1033b288129SPeter Brune ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 1043b288129SPeter Brune 10571b4ebd2SPeter Brune if (ynorm == 0.0) { 10671b4ebd2SPeter Brune if (monitor) { 107ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 10871b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 109ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 11071b4ebd2SPeter Brune } 11171b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 11271b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 11371cee68bSPeter Brune ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr); 114f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 11571b4ebd2SPeter Brune PetscFunctionReturn(0); 11671b4ebd2SPeter Brune } 11771b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 11871b4ebd2SPeter Brune if (monitor) { 119ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 120c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr); 121ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 12271b4ebd2SPeter Brune } 12371b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 12471b4ebd2SPeter Brune ynorm = maxstep; 12571b4ebd2SPeter Brune } 1268a430afbSPeter Brune 1278a430afbSPeter Brune /* if the SNES has an objective set, use that instead of the function value */ 128bd4a8a71SBarry Smith if (objective) { 1298a430afbSPeter Brune ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr); 1308a430afbSPeter Brune } else { 1318a430afbSPeter Brune f = fnorm*fnorm; 1328a430afbSPeter Brune } 1338a430afbSPeter Brune 1348a430afbSPeter Brune /* compute the initial slope */ 135bd4a8a71SBarry Smith if (objective) { 1368a430afbSPeter Brune /* slope comes from the function (assumed to be the gradient of the objective */ 13767392de3SBarry Smith ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr); 1388a430afbSPeter Brune } else { 1398a430afbSPeter Brune /* slope comes from the normal equations */ 14071b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 14167392de3SBarry Smith ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr); 14271b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 14371b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 1448a430afbSPeter Brune } 14571b4ebd2SPeter Brune 146c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1479bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1489bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1499bd66eb0SPeter Brune } 15071b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 15171b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 15271b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 153f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 15471b4ebd2SPeter Brune PetscFunctionReturn(0); 15571b4ebd2SPeter Brune } 1568a430afbSPeter Brune 157bd4a8a71SBarry Smith if (objective) { 1588a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 1598a430afbSPeter Brune } else { 160ce5a860aSPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); 161c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 162c427506bSPeter Brune if (domainerror) { 163f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 16471b4ebd2SPeter Brune PetscFunctionReturn(0); 16571b4ebd2SPeter Brune } 1669bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1679bd66eb0SPeter Brune gnorm = fnorm; 1689bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1699bd66eb0SPeter Brune } else { 17071b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1718a430afbSPeter Brune } 17279e7e81bSJed Brown g = PetscSqr(gnorm); 1739bd66eb0SPeter Brune } 1749bd66eb0SPeter Brune 175c98378a5SBarry Smith if (PetscIsInfOrNanReal(g)) { 176c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 177c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 178c98378a5SBarry Smith PetscFunctionReturn(0); 179c98378a5SBarry Smith } 180bd4a8a71SBarry Smith if (!objective) { 181c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1828a430afbSPeter Brune } 1838a430afbSPeter Brune if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */ 18471b4ebd2SPeter Brune if (monitor) { 185ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 186bd4a8a71SBarry Smith if (!objective) { 187c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1888a430afbSPeter Brune } else { 1898a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr); 1908a430afbSPeter Brune } 191ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 19271b4ebd2SPeter Brune } 19371b4ebd2SPeter Brune } else { 194c21ba15cSPeter Brune /* Since the full step didn't work and the step is tiny, quit */ 195c21ba15cSPeter Brune if (stol*xnorm > ynorm) { 196c21ba15cSPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 197c21ba15cSPeter Brune ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr); 198c21ba15cSPeter Brune PetscFunctionReturn(0); 199c21ba15cSPeter Brune } 20071b4ebd2SPeter Brune /* Fit points with quadratic */ 2018a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*lambda*initslope); 20271b4ebd2SPeter Brune lambdaprev = lambda; 2038a430afbSPeter Brune gprev = g; 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 20871b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2099bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2109bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2119bd66eb0SPeter Brune } 21271b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 21371b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 21471b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 215f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 21671b4ebd2SPeter Brune PetscFunctionReturn(0); 21771b4ebd2SPeter Brune } 218bd4a8a71SBarry Smith if (objective) { 2198a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 2208a430afbSPeter Brune } else { 221ce5a860aSPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); 2226188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 223f5af7f23SKarl Rupp if (domainerror) PetscFunctionReturn(0); 224f5af7f23SKarl Rupp 2259bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2269bd66eb0SPeter Brune gnorm = fnorm; 2279bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2289bd66eb0SPeter Brune } else { 22971b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2309bd66eb0SPeter Brune } 231f1ab69e3SPeter Brune g = PetscSqr(gnorm); 2328a430afbSPeter Brune } 233c98378a5SBarry Smith if (PetscIsInfOrNanReal(g)) { 234c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 235c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 236c98378a5SBarry Smith PetscFunctionReturn(0); 237c98378a5SBarry Smith } 23871b4ebd2SPeter Brune if (monitor) { 239ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 240bd4a8a71SBarry Smith if (!objective) { 241c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); 2428a430afbSPeter Brune } else { 2438a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr); 2448a430afbSPeter Brune } 245ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24671b4ebd2SPeter Brune } 2478a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */ 24871b4ebd2SPeter Brune if (monitor) { 249ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 251ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25271b4ebd2SPeter Brune } 25371b4ebd2SPeter Brune } else { 25471b4ebd2SPeter Brune /* Fit points with cubic */ 25571b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 25671b4ebd2SPeter Brune if (lambda <= minlambda) { 25771b4ebd2SPeter Brune if (monitor) { 258ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25971b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 260bd4a8a71SBarry Smith if (!objective) { 26171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 26271b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 263c69d1a72SBarry Smith (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2648a430afbSPeter Brune } else { 2658a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 2668a430afbSPeter Brune " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 2678a430afbSPeter Brune (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2688a430afbSPeter Brune } 269ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 27071b4ebd2SPeter Brune } 271f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 27271b4ebd2SPeter Brune PetscFunctionReturn(0); 27371b4ebd2SPeter Brune } 274b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 2758a430afbSPeter Brune t1 = .5*(g - f) - lambda*initslope; 2768a430afbSPeter Brune t2 = .5*(gprev - f) - lambdaprev*initslope; 27771b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27871b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27971b4ebd2SPeter Brune d = b*b - 3*a*initslope; 28071b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 281f5af7f23SKarl Rupp if (a == 0.0) lambdatemp = -initslope/(2.0*b); 282f5af7f23SKarl Rupp else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 283f5af7f23SKarl Rupp 284b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 2858a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*initslope); 286ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt"); 28771b4ebd2SPeter Brune lambdaprev = lambda; 2888a430afbSPeter Brune gprev = g; 28971b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 29071b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 29171b4ebd2SPeter Brune else lambda = lambdatemp; 29271b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 293b7b2e573SPeter Brune if (linesearch->ops->viproject) { 294b7b2e573SPeter Brune ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); 295b7b2e573SPeter Brune } 29671b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 29771b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 298bd4a8a71SBarry Smith if (!objective) { 29971b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 300c69d1a72SBarry Smith (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 3018a430afbSPeter Brune } 302f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 30371b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 30471b4ebd2SPeter Brune PetscFunctionReturn(0); 30571b4ebd2SPeter Brune } 306bd4a8a71SBarry Smith if (objective) { 3078a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 3088a430afbSPeter Brune } else { 309ce5a860aSPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); 310c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 311c427506bSPeter Brune if (domainerror) { 312f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 31371b4ebd2SPeter Brune PetscFunctionReturn(0); 31471b4ebd2SPeter Brune } 3159bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3169bd66eb0SPeter Brune gnorm = fnorm; 3179bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3189bd66eb0SPeter Brune } else { 31971b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3208a430afbSPeter Brune } 321f1ab69e3SPeter Brune g = PetscSqr(gnorm); 3229bd66eb0SPeter Brune } 323c98378a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) { 324c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 325c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 326c98378a5SBarry Smith PetscFunctionReturn(0); 327c98378a5SBarry Smith } 3288a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */ 32971b4ebd2SPeter Brune if (monitor) { 330ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 331bd4a8a71SBarry Smith if (!objective) { 332b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 333c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 33471b4ebd2SPeter Brune } else { 335c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 33671b4ebd2SPeter Brune } 337ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3388a430afbSPeter Brune } else { 3398a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3408a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3418a430afbSPeter Brune } else { 3428a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3438a430afbSPeter Brune } 3448a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3458a430afbSPeter Brune } 34671b4ebd2SPeter Brune } 34771b4ebd2SPeter Brune break; 348f5af7f23SKarl Rupp } else if (monitor) { 349ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 350bd4a8a71SBarry Smith if (!objective) { 351b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 352c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 35371b4ebd2SPeter Brune } else { 354c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 35571b4ebd2SPeter Brune } 356ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3578a430afbSPeter Brune } else { 3588a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3598a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3608a430afbSPeter Brune } else { 3618a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3628a430afbSPeter Brune } 3638a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3648a430afbSPeter Brune } 36571b4ebd2SPeter Brune } 36671b4ebd2SPeter Brune } 36771b4ebd2SPeter Brune } 36871b4ebd2SPeter Brune } 36971b4ebd2SPeter Brune 37071b4ebd2SPeter Brune /* postcheck */ 3717b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 37271b4ebd2SPeter Brune if (changed_y) { 37371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 3749bd66eb0SPeter Brune if (linesearch->ops->viproject) { 3759bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 3769bd66eb0SPeter Brune } 37771b4ebd2SPeter Brune } 378bd4a8a71SBarry Smith if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ 379ce5a860aSPeter Brune ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); 38071b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 38171b4ebd2SPeter Brune if (domainerror) { 382f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 38371b4ebd2SPeter Brune PetscFunctionReturn(0); 38471b4ebd2SPeter Brune } 3859bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3869bd66eb0SPeter Brune gnorm = fnorm; 3879bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3889bd66eb0SPeter Brune } else { 3899bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3909bd66eb0SPeter Brune } 3919bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 392c98378a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) { 393c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 394c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 395c98378a5SBarry Smith PetscFunctionReturn(0); 396c98378a5SBarry Smith } 397c427506bSPeter Brune } 39871b4ebd2SPeter Brune 39971b4ebd2SPeter Brune /* copy the solution over */ 40071b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 401c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 402c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 403f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 404f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 40571b4ebd2SPeter Brune PetscFunctionReturn(0); 40671b4ebd2SPeter Brune } 40771b4ebd2SPeter Brune 40871b4ebd2SPeter Brune #undef __FUNCT__ 4097f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 4107f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 4117f1410a3SPeter Brune { 4127f1410a3SPeter Brune PetscErrorCode ierr; 4137f1410a3SPeter Brune PetscBool iascii; 4147f1410a3SPeter Brune SNESLineSearch_BT *bt; 415075cc632SBarry Smith 4167f1410a3SPeter Brune PetscFunctionBegin; 417251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4187f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 4197f1410a3SPeter Brune if (iascii) { 420b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 4217f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 422b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 4237f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 4247f1410a3SPeter Brune } 425007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 4267f1410a3SPeter Brune } 4277f1410a3SPeter Brune PetscFunctionReturn(0); 4287f1410a3SPeter Brune } 4297f1410a3SPeter Brune 4307f1410a3SPeter Brune 4317f1410a3SPeter Brune #undef __FUNCT__ 432f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 433f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 434731c067cSJed Brown { 435731c067cSJed Brown PetscErrorCode ierr; 43671b4ebd2SPeter Brune 43771b4ebd2SPeter Brune PetscFunctionBegin; 438731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 43971b4ebd2SPeter Brune PetscFunctionReturn(0); 44071b4ebd2SPeter Brune } 44171b4ebd2SPeter Brune 44271b4ebd2SPeter Brune 44371b4ebd2SPeter Brune #undef __FUNCT__ 444f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 445f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 44671b4ebd2SPeter Brune { 44771b4ebd2SPeter Brune 44871b4ebd2SPeter Brune PetscErrorCode ierr; 449f1c6b773SPeter Brune SNESLineSearch_BT *bt; 45071b4ebd2SPeter Brune 4516e111a19SKarl Rupp PetscFunctionBegin; 452f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 45371b4ebd2SPeter Brune 454f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 4550298fd71SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr); 45671b4ebd2SPeter Brune 45771b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 45871b4ebd2SPeter Brune PetscFunctionReturn(0); 45971b4ebd2SPeter Brune } 46071b4ebd2SPeter Brune 46171b4ebd2SPeter Brune 46271b4ebd2SPeter Brune #undef __FUNCT__ 463f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 464954494b2SJed Brown /*MC 465db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 466954494b2SJed Brown 467db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 468db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 469db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 470954494b2SJed Brown 471cd7522eaSPeter Brune Options Database Keys: 472cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 473db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 474db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 475db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 476cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 477cd7522eaSPeter Brune 478954494b2SJed Brown Level: advanced 479954494b2SJed Brown 480db609ea7SPeter Brune Notes: 481db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 482db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 483db609ea7SPeter Brune 484f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 485954494b2SJed Brown 486f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 487954494b2SJed Brown M*/ 4888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 48971b4ebd2SPeter Brune { 49071b4ebd2SPeter Brune 491f1c6b773SPeter Brune SNESLineSearch_BT *bt; 49271b4ebd2SPeter Brune PetscErrorCode ierr; 49371b4ebd2SPeter Brune 49471b4ebd2SPeter Brune PetscFunctionBegin; 495f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 496f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 497f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 4980298fd71SBarry Smith linesearch->ops->reset = NULL; 4997f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 5000298fd71SBarry Smith linesearch->ops->setup = NULL; 50171b4ebd2SPeter Brune 502*b00a9115SJed Brown ierr = PetscNewLog(linesearch,&bt);CHKERRQ(ierr); 503f5af7f23SKarl Rupp 50471b4ebd2SPeter Brune linesearch->data = (void*)bt; 50571b4ebd2SPeter Brune linesearch->max_its = 40; 506b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 50771b4ebd2SPeter Brune bt->alpha = 1e-4; 50871b4ebd2SPeter Brune PetscFunctionReturn(0); 50971b4ebd2SPeter Brune } 510