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); 113f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 11471b4ebd2SPeter Brune PetscFunctionReturn(0); 11571b4ebd2SPeter Brune } 11671b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 11771b4ebd2SPeter Brune if (monitor) { 118ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 119c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr); 120ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 12171b4ebd2SPeter Brune } 12271b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 12371b4ebd2SPeter Brune ynorm = maxstep; 12471b4ebd2SPeter Brune } 1258a430afbSPeter Brune 1268a430afbSPeter Brune /* if the SNES has an objective set, use that instead of the function value */ 127bd4a8a71SBarry Smith if (objective) { 1288a430afbSPeter Brune ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr); 1298a430afbSPeter Brune } else { 1308a430afbSPeter Brune f = fnorm*fnorm; 1318a430afbSPeter Brune } 1328a430afbSPeter Brune 1338a430afbSPeter Brune /* compute the initial slope */ 134bd4a8a71SBarry Smith if (objective) { 1358a430afbSPeter Brune /* slope comes from the function (assumed to be the gradient of the objective */ 13667392de3SBarry Smith ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr); 1378a430afbSPeter Brune } else { 1388a430afbSPeter Brune /* slope comes from the normal equations */ 13971b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 14067392de3SBarry Smith ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr); 14171b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 14271b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 1438a430afbSPeter Brune } 14471b4ebd2SPeter Brune 145c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1469bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1479bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1489bd66eb0SPeter Brune } 14971b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 15071b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 15171b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 152f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 15371b4ebd2SPeter Brune PetscFunctionReturn(0); 15471b4ebd2SPeter Brune } 1558a430afbSPeter Brune 156bd4a8a71SBarry Smith if (objective) { 1578a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 1588a430afbSPeter Brune } else { 15971b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 160c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 161c427506bSPeter Brune if (domainerror) { 162f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 16371b4ebd2SPeter Brune PetscFunctionReturn(0); 16471b4ebd2SPeter Brune } 1659bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1669bd66eb0SPeter Brune gnorm = fnorm; 1679bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1689bd66eb0SPeter Brune } else { 16971b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1708a430afbSPeter Brune } 17179e7e81bSJed Brown g = PetscSqr(gnorm); 1729bd66eb0SPeter Brune } 1739bd66eb0SPeter Brune 174c98378a5SBarry Smith if (PetscIsInfOrNanReal(g)) { 175c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 176c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 177c98378a5SBarry Smith PetscFunctionReturn(0); 178c98378a5SBarry Smith } 179bd4a8a71SBarry Smith if (!objective) { 180c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1818a430afbSPeter Brune } 1828a430afbSPeter Brune if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */ 18371b4ebd2SPeter Brune if (monitor) { 184ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 185bd4a8a71SBarry Smith if (!objective) { 186c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1878a430afbSPeter Brune } else { 1888a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr); 1898a430afbSPeter Brune } 190ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 19171b4ebd2SPeter Brune } 19271b4ebd2SPeter Brune } else { 193c21ba15cSPeter Brune /* Since the full step didn't work and the step is tiny, quit */ 194c21ba15cSPeter Brune if (stol*xnorm > ynorm) { 195c21ba15cSPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 196c21ba15cSPeter Brune ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr); 197c21ba15cSPeter Brune PetscFunctionReturn(0); 198c21ba15cSPeter Brune } 19971b4ebd2SPeter Brune /* Fit points with quadratic */ 2008a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*lambda*initslope); 20171b4ebd2SPeter Brune lambdaprev = lambda; 2028a430afbSPeter Brune gprev = g; 20371b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 20471b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 20571b4ebd2SPeter Brune else lambda = lambdatemp; 20671b4ebd2SPeter Brune 20771b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2089bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2099bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2109bd66eb0SPeter Brune } 21171b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 21271b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 21371b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 214f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 21571b4ebd2SPeter Brune PetscFunctionReturn(0); 21671b4ebd2SPeter Brune } 217bd4a8a71SBarry Smith if (objective) { 2188a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 2198a430afbSPeter Brune } else { 22071b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 2216188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 222f5af7f23SKarl Rupp if (domainerror) PetscFunctionReturn(0); 223f5af7f23SKarl Rupp 2249bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2259bd66eb0SPeter Brune gnorm = fnorm; 2269bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2279bd66eb0SPeter Brune } else { 22871b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2298a430afbSPeter Brune g = gnorm*gnorm; 2309bd66eb0SPeter Brune } 2318a430afbSPeter Brune } 232c98378a5SBarry Smith if (PetscIsInfOrNanReal(g)) { 233c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 234c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 235c98378a5SBarry Smith PetscFunctionReturn(0); 236c98378a5SBarry Smith } 23771b4ebd2SPeter Brune if (monitor) { 238ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 239bd4a8a71SBarry Smith if (!objective) { 240c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); 2418a430afbSPeter Brune } else { 2428a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr); 2438a430afbSPeter Brune } 244ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24571b4ebd2SPeter Brune } 2468a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */ 24771b4ebd2SPeter Brune if (monitor) { 248ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24971b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 250ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25171b4ebd2SPeter Brune } 25271b4ebd2SPeter Brune } else { 25371b4ebd2SPeter Brune /* Fit points with cubic */ 25471b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 25571b4ebd2SPeter Brune if (lambda <= minlambda) { 25671b4ebd2SPeter Brune if (monitor) { 257ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25871b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 259bd4a8a71SBarry Smith if (!objective) { 26071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 26171b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 262c69d1a72SBarry Smith (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2638a430afbSPeter Brune } else { 2648a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 2658a430afbSPeter Brune " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 2668a430afbSPeter Brune (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2678a430afbSPeter Brune } 268ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 26971b4ebd2SPeter Brune } 270f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 27171b4ebd2SPeter Brune PetscFunctionReturn(0); 27271b4ebd2SPeter Brune } 273b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 2748a430afbSPeter Brune t1 = .5*(g - f) - lambda*initslope; 2758a430afbSPeter Brune t2 = .5*(gprev - f) - lambdaprev*initslope; 27671b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27771b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27871b4ebd2SPeter Brune d = b*b - 3*a*initslope; 27971b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 280f5af7f23SKarl Rupp if (a == 0.0) lambdatemp = -initslope/(2.0*b); 281f5af7f23SKarl Rupp else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 282f5af7f23SKarl Rupp 283b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 2848a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*initslope); 285ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt"); 28671b4ebd2SPeter Brune lambdaprev = lambda; 2878a430afbSPeter Brune gprev = g; 28871b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 28971b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 29071b4ebd2SPeter Brune else lambda = lambdatemp; 29171b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 292b7b2e573SPeter Brune if (linesearch->ops->viproject) { 293b7b2e573SPeter Brune ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); 294b7b2e573SPeter Brune } 29571b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 29671b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 297bd4a8a71SBarry Smith if (!objective) { 29871b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 299c69d1a72SBarry Smith (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 3008a430afbSPeter Brune } 301f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 30271b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 30371b4ebd2SPeter Brune PetscFunctionReturn(0); 30471b4ebd2SPeter Brune } 305bd4a8a71SBarry Smith if (objective) { 3068a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 3078a430afbSPeter Brune } else { 30871b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 309c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 310c427506bSPeter Brune if (domainerror) { 311f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 31271b4ebd2SPeter Brune PetscFunctionReturn(0); 31371b4ebd2SPeter Brune } 3149bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3159bd66eb0SPeter Brune gnorm = fnorm; 3169bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3179bd66eb0SPeter Brune } else { 31871b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3198a430afbSPeter Brune g = gnorm*gnorm; 3208a430afbSPeter Brune } 3219bd66eb0SPeter Brune } 322c98378a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) { 323c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 324c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 325c98378a5SBarry Smith PetscFunctionReturn(0); 326c98378a5SBarry Smith } 3278a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */ 32871b4ebd2SPeter Brune if (monitor) { 329ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 330bd4a8a71SBarry Smith if (!objective) { 331b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 332c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 33371b4ebd2SPeter Brune } else { 334c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 33571b4ebd2SPeter Brune } 336ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3378a430afbSPeter Brune } else { 3388a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3398a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3408a430afbSPeter Brune } else { 3418a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3428a430afbSPeter Brune } 3438a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3448a430afbSPeter Brune } 34571b4ebd2SPeter Brune } 34671b4ebd2SPeter Brune break; 347f5af7f23SKarl Rupp } else if (monitor) { 348ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 349bd4a8a71SBarry Smith if (!objective) { 350b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 351c69d1a72SBarry 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); 35271b4ebd2SPeter Brune } else { 353c69d1a72SBarry 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); 35471b4ebd2SPeter Brune } 355ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3568a430afbSPeter Brune } else { 3578a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3588a430afbSPeter 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); 3598a430afbSPeter Brune } else { 3608a430afbSPeter 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); 3618a430afbSPeter Brune } 3628a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3638a430afbSPeter Brune } 36471b4ebd2SPeter Brune } 36571b4ebd2SPeter Brune } 36671b4ebd2SPeter Brune } 36771b4ebd2SPeter Brune } 36871b4ebd2SPeter Brune 36971b4ebd2SPeter Brune /* postcheck */ 3707b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 37171b4ebd2SPeter Brune if (changed_y) { 37271b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 3739bd66eb0SPeter Brune if (linesearch->ops->viproject) { 3749bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 3759bd66eb0SPeter Brune } 37671b4ebd2SPeter Brune } 377bd4a8a71SBarry Smith if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ 378c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 37971b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 38071b4ebd2SPeter Brune if (domainerror) { 381f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 38271b4ebd2SPeter Brune PetscFunctionReturn(0); 38371b4ebd2SPeter Brune } 3849bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3859bd66eb0SPeter Brune gnorm = fnorm; 3869bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3879bd66eb0SPeter Brune } else { 3889bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3899bd66eb0SPeter Brune } 3909bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 391c98378a5SBarry Smith if (PetscIsInfOrNanReal(gnorm)) { 392c98378a5SBarry Smith ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 393c98378a5SBarry Smith ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); 394c98378a5SBarry Smith PetscFunctionReturn(0); 395c98378a5SBarry Smith } 396c427506bSPeter Brune } 39771b4ebd2SPeter Brune 39871b4ebd2SPeter Brune /* copy the solution over */ 39971b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 400c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 401c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 402f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 403f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 40471b4ebd2SPeter Brune PetscFunctionReturn(0); 40571b4ebd2SPeter Brune } 40671b4ebd2SPeter Brune 40771b4ebd2SPeter Brune #undef __FUNCT__ 4087f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 4097f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 4107f1410a3SPeter Brune { 4117f1410a3SPeter Brune PetscErrorCode ierr; 4127f1410a3SPeter Brune PetscBool iascii; 4137f1410a3SPeter Brune SNESLineSearch_BT *bt; 414075cc632SBarry Smith 4157f1410a3SPeter Brune PetscFunctionBegin; 416251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4177f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 4187f1410a3SPeter Brune if (iascii) { 419b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 4207f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 421b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 4227f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 4237f1410a3SPeter Brune } 424007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 4257f1410a3SPeter Brune } 4267f1410a3SPeter Brune PetscFunctionReturn(0); 4277f1410a3SPeter Brune } 4287f1410a3SPeter Brune 4297f1410a3SPeter Brune 4307f1410a3SPeter Brune #undef __FUNCT__ 431f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 432f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 433731c067cSJed Brown { 434731c067cSJed Brown PetscErrorCode ierr; 43571b4ebd2SPeter Brune 43671b4ebd2SPeter Brune PetscFunctionBegin; 437731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 43871b4ebd2SPeter Brune PetscFunctionReturn(0); 43971b4ebd2SPeter Brune } 44071b4ebd2SPeter Brune 44171b4ebd2SPeter Brune 44271b4ebd2SPeter Brune #undef __FUNCT__ 443f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 444f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 44571b4ebd2SPeter Brune { 44671b4ebd2SPeter Brune 44771b4ebd2SPeter Brune PetscErrorCode ierr; 448f1c6b773SPeter Brune SNESLineSearch_BT *bt; 44971b4ebd2SPeter Brune 4506e111a19SKarl Rupp PetscFunctionBegin; 451f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 45271b4ebd2SPeter Brune 453f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 4540298fd71SBarry Smith ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, NULL);CHKERRQ(ierr); 45571b4ebd2SPeter Brune 45671b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 45771b4ebd2SPeter Brune PetscFunctionReturn(0); 45871b4ebd2SPeter Brune } 45971b4ebd2SPeter Brune 46071b4ebd2SPeter Brune 46171b4ebd2SPeter Brune #undef __FUNCT__ 462f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 463954494b2SJed Brown /*MC 464db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 465954494b2SJed Brown 466db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 467db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 468db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 469954494b2SJed Brown 470cd7522eaSPeter Brune Options Database Keys: 471cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 472db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 473db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 474db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 475cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 476cd7522eaSPeter Brune 477954494b2SJed Brown Level: advanced 478954494b2SJed Brown 479db609ea7SPeter Brune Notes: 480db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 481db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 482db609ea7SPeter Brune 483f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 484954494b2SJed Brown 485f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 486954494b2SJed Brown M*/ 487*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 48871b4ebd2SPeter Brune { 48971b4ebd2SPeter Brune 490f1c6b773SPeter Brune SNESLineSearch_BT *bt; 49171b4ebd2SPeter Brune PetscErrorCode ierr; 49271b4ebd2SPeter Brune 49371b4ebd2SPeter Brune PetscFunctionBegin; 494f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 495f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 496f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 4970298fd71SBarry Smith linesearch->ops->reset = NULL; 4987f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 4990298fd71SBarry Smith linesearch->ops->setup = NULL; 50071b4ebd2SPeter Brune 501f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 502f5af7f23SKarl Rupp 50371b4ebd2SPeter Brune linesearch->data = (void*)bt; 50471b4ebd2SPeter Brune linesearch->max_its = 40; 505b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 50671b4ebd2SPeter Brune bt->alpha = 1e-4; 50771b4ebd2SPeter Brune PetscFunctionReturn(0); 50871b4ebd2SPeter Brune } 509