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 { 16071b4ebd2SPeter Brune ierr = SNESComputeFunction(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) { 196*ced04f9dSPeter Brune ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr); 197c21ba15cSPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 198c21ba15cSPeter Brune ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr); 199c21ba15cSPeter Brune PetscFunctionReturn(0); 200c21ba15cSPeter Brune } 20171b4ebd2SPeter Brune /* Fit points with quadratic */ 2028a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*lambda*initslope); 20371b4ebd2SPeter Brune lambdaprev = lambda; 2048a430afbSPeter Brune gprev = g; 20571b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 20671b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 20771b4ebd2SPeter Brune else lambda = lambdatemp; 20871b4ebd2SPeter Brune 20971b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2109bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2119bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2129bd66eb0SPeter Brune } 21371b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 21471b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 21571b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 216f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 21771b4ebd2SPeter Brune PetscFunctionReturn(0); 21871b4ebd2SPeter Brune } 219bd4a8a71SBarry Smith if (objective) { 2208a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 2218a430afbSPeter Brune } else { 22271b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 2236188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 224f5af7f23SKarl Rupp if (domainerror) PetscFunctionReturn(0); 225f5af7f23SKarl Rupp 2269bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2279bd66eb0SPeter Brune gnorm = fnorm; 2289bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2299bd66eb0SPeter Brune } else { 23071b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2319bd66eb0SPeter Brune } 232f1ab69e3SPeter Brune g = PetscSqr(gnorm); 2338a430afbSPeter Brune } 234c98378a5SBarry Smith if (PetscIsInfOrNanReal(g)) { 235c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 236c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 237c98378a5SBarry Smith PetscFunctionReturn(0); 238c98378a5SBarry Smith } 23971b4ebd2SPeter Brune if (monitor) { 240ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 241bd4a8a71SBarry Smith if (!objective) { 242c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); 2438a430afbSPeter Brune } else { 2448a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr); 2458a430afbSPeter Brune } 246ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24771b4ebd2SPeter Brune } 2488a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */ 24971b4ebd2SPeter Brune if (monitor) { 250ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 252ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25371b4ebd2SPeter Brune } 25471b4ebd2SPeter Brune } else { 25571b4ebd2SPeter Brune /* Fit points with cubic */ 25671b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 25771b4ebd2SPeter Brune if (lambda <= minlambda) { 25871b4ebd2SPeter Brune if (monitor) { 259ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 26071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 261bd4a8a71SBarry Smith if (!objective) { 26271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 26371b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 264c69d1a72SBarry Smith (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2658a430afbSPeter Brune } else { 2668a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 2678a430afbSPeter Brune " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 2688a430afbSPeter Brune (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2698a430afbSPeter Brune } 270ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 27171b4ebd2SPeter Brune } 272f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 27371b4ebd2SPeter Brune PetscFunctionReturn(0); 27471b4ebd2SPeter Brune } 275b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 2768a430afbSPeter Brune t1 = .5*(g - f) - lambda*initslope; 2778a430afbSPeter Brune t2 = .5*(gprev - f) - lambdaprev*initslope; 27871b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27971b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 28071b4ebd2SPeter Brune d = b*b - 3*a*initslope; 28171b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 282f5af7f23SKarl Rupp if (a == 0.0) lambdatemp = -initslope/(2.0*b); 283f5af7f23SKarl Rupp else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 284f5af7f23SKarl Rupp 285b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 2868a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*initslope); 287ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt"); 28871b4ebd2SPeter Brune lambdaprev = lambda; 2898a430afbSPeter Brune gprev = g; 29071b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 29171b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 29271b4ebd2SPeter Brune else lambda = lambdatemp; 29371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 294b7b2e573SPeter Brune if (linesearch->ops->viproject) { 295b7b2e573SPeter Brune ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); 296b7b2e573SPeter Brune } 29771b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 29871b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 299bd4a8a71SBarry Smith if (!objective) { 30071b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 301c69d1a72SBarry Smith (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 3028a430afbSPeter Brune } 303f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 30471b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 30571b4ebd2SPeter Brune PetscFunctionReturn(0); 30671b4ebd2SPeter Brune } 307bd4a8a71SBarry Smith if (objective) { 3088a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 3098a430afbSPeter Brune } else { 31071b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 311c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 312c427506bSPeter Brune if (domainerror) { 313f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 31471b4ebd2SPeter Brune PetscFunctionReturn(0); 31571b4ebd2SPeter Brune } 3169bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3179bd66eb0SPeter Brune gnorm = fnorm; 3189bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3199bd66eb0SPeter Brune } else { 32071b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3218a430afbSPeter Brune } 322f1ab69e3SPeter Brune g = PetscSqr(gnorm); 3239bd66eb0SPeter Brune } 324c98378a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) { 325c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 326c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 327c98378a5SBarry Smith PetscFunctionReturn(0); 328c98378a5SBarry Smith } 3298a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */ 33071b4ebd2SPeter Brune if (monitor) { 331ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 332bd4a8a71SBarry Smith if (!objective) { 333b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 334c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 33571b4ebd2SPeter Brune } else { 336c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 33771b4ebd2SPeter Brune } 338ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3398a430afbSPeter Brune } else { 3408a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3418a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3428a430afbSPeter Brune } else { 3438a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3448a430afbSPeter Brune } 3458a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3468a430afbSPeter Brune } 34771b4ebd2SPeter Brune } 34871b4ebd2SPeter Brune break; 349f5af7f23SKarl Rupp } else if (monitor) { 350ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 351bd4a8a71SBarry Smith if (!objective) { 352b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 353c69d1a72SBarry 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); 35471b4ebd2SPeter Brune } else { 355c69d1a72SBarry 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); 35671b4ebd2SPeter Brune } 357ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3588a430afbSPeter Brune } else { 3598a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3608a430afbSPeter 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); 3618a430afbSPeter Brune } else { 3628a430afbSPeter 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); 3638a430afbSPeter Brune } 3648a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3658a430afbSPeter Brune } 36671b4ebd2SPeter Brune } 36771b4ebd2SPeter Brune } 36871b4ebd2SPeter Brune } 36971b4ebd2SPeter Brune } 37071b4ebd2SPeter Brune 37171b4ebd2SPeter Brune /* postcheck */ 3727b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 37371b4ebd2SPeter Brune if (changed_y) { 37471b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 3759bd66eb0SPeter Brune if (linesearch->ops->viproject) { 3769bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 3779bd66eb0SPeter Brune } 37871b4ebd2SPeter Brune } 379bd4a8a71SBarry Smith if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ 380c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 38171b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 38271b4ebd2SPeter Brune if (domainerror) { 383f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 38471b4ebd2SPeter Brune PetscFunctionReturn(0); 38571b4ebd2SPeter Brune } 3869bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3879bd66eb0SPeter Brune gnorm = fnorm; 3889bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3899bd66eb0SPeter Brune } else { 3909bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3919bd66eb0SPeter Brune } 3929bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 393c98378a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) { 394c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 395c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 396c98378a5SBarry Smith PetscFunctionReturn(0); 397c98378a5SBarry Smith } 398c427506bSPeter Brune } 39971b4ebd2SPeter Brune 40071b4ebd2SPeter Brune /* copy the solution over */ 40171b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 402c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 403c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 404f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 405f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 40671b4ebd2SPeter Brune PetscFunctionReturn(0); 40771b4ebd2SPeter Brune } 40871b4ebd2SPeter Brune 40971b4ebd2SPeter Brune #undef __FUNCT__ 4107f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 4117f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 4127f1410a3SPeter Brune { 4137f1410a3SPeter Brune PetscErrorCode ierr; 4147f1410a3SPeter Brune PetscBool iascii; 4157f1410a3SPeter Brune SNESLineSearch_BT *bt; 416075cc632SBarry Smith 4177f1410a3SPeter Brune PetscFunctionBegin; 418251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4197f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 4207f1410a3SPeter Brune if (iascii) { 421b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 4227f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 423b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 4247f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 4257f1410a3SPeter Brune } 426007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 4277f1410a3SPeter Brune } 4287f1410a3SPeter Brune PetscFunctionReturn(0); 4297f1410a3SPeter Brune } 4307f1410a3SPeter Brune 4317f1410a3SPeter Brune 4327f1410a3SPeter Brune #undef __FUNCT__ 433f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 434f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 435731c067cSJed Brown { 436731c067cSJed Brown PetscErrorCode ierr; 43771b4ebd2SPeter Brune 43871b4ebd2SPeter Brune PetscFunctionBegin; 439731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 44071b4ebd2SPeter Brune PetscFunctionReturn(0); 44171b4ebd2SPeter Brune } 44271b4ebd2SPeter Brune 44371b4ebd2SPeter Brune 44471b4ebd2SPeter Brune #undef __FUNCT__ 445f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 446f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 44771b4ebd2SPeter Brune { 44871b4ebd2SPeter Brune 44971b4ebd2SPeter Brune PetscErrorCode ierr; 450f1c6b773SPeter Brune SNESLineSearch_BT *bt; 45171b4ebd2SPeter Brune 4526e111a19SKarl Rupp PetscFunctionBegin; 453f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 45471b4ebd2SPeter Brune 455f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 4560298fd71SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr); 45771b4ebd2SPeter Brune 45871b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 45971b4ebd2SPeter Brune PetscFunctionReturn(0); 46071b4ebd2SPeter Brune } 46171b4ebd2SPeter Brune 46271b4ebd2SPeter Brune 46371b4ebd2SPeter Brune #undef __FUNCT__ 464f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 465954494b2SJed Brown /*MC 466db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 467954494b2SJed Brown 468db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 469db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 470db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 471954494b2SJed Brown 472cd7522eaSPeter Brune Options Database Keys: 473cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 474db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 475db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 476db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 477cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 478cd7522eaSPeter Brune 479954494b2SJed Brown Level: advanced 480954494b2SJed Brown 481db609ea7SPeter Brune Notes: 482db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 483db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 484db609ea7SPeter Brune 485f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 486954494b2SJed Brown 487f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 488954494b2SJed Brown M*/ 4898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 49071b4ebd2SPeter Brune { 49171b4ebd2SPeter Brune 492f1c6b773SPeter Brune SNESLineSearch_BT *bt; 49371b4ebd2SPeter Brune PetscErrorCode ierr; 49471b4ebd2SPeter Brune 49571b4ebd2SPeter Brune PetscFunctionBegin; 496f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 497f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 498f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 4990298fd71SBarry Smith linesearch->ops->reset = NULL; 5007f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 5010298fd71SBarry Smith linesearch->ops->setup = NULL; 50271b4ebd2SPeter Brune 503f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 504f5af7f23SKarl Rupp 50571b4ebd2SPeter Brune linesearch->data = (void*)bt; 50671b4ebd2SPeter Brune linesearch->max_its = 40; 507b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 50871b4ebd2SPeter Brune bt->alpha = 1e-4; 50971b4ebd2SPeter Brune PetscFunctionReturn(0); 51071b4ebd2SPeter Brune } 511