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); 8522d28d08SBarry Smith ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,PETSC_NULL,PETSC_NULL,PETSC_NULL,&max_its);CHKERRQ(ierr); 863b288129SPeter Brune ierr = SNESGetTolerances(snes,PETSC_NULL,PETSC_NULL,&stol,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 87bd4a8a71SBarry Smith ierr = SNESGetObjective(snes,&objective,PETSC_NULL);CHKERRQ(ierr); 88f1c6b773SPeter Brune bt = (SNESLineSearch_BT *)linesearch->data; 8971b4ebd2SPeter Brune 9071b4ebd2SPeter Brune alpha = bt->alpha; 9171b4ebd2SPeter Brune 9271b4ebd2SPeter Brune ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 938a430afbSPeter Brune 94bd4a8a71SBarry Smith if (!jac && !objective) SETERRQ(((PetscObject)linesearch)->comm, 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 1748a430afbSPeter Brune if (PetscIsInfOrNanReal(g)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 175bd4a8a71SBarry Smith if (!objective) { 176c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1778a430afbSPeter Brune } 1788a430afbSPeter Brune if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */ 17971b4ebd2SPeter Brune if (monitor) { 180ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 181bd4a8a71SBarry Smith if (!objective) { 182c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1838a430afbSPeter Brune } else { 1848a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr); 1858a430afbSPeter Brune } 186ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 18771b4ebd2SPeter Brune } 18871b4ebd2SPeter Brune } else { 189c21ba15cSPeter Brune /* Since the full step didn't work and the step is tiny, quit */ 190c21ba15cSPeter Brune if (stol*xnorm > ynorm) { 191c21ba15cSPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 192c21ba15cSPeter Brune ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr); 193c21ba15cSPeter Brune PetscFunctionReturn(0); 194c21ba15cSPeter Brune } 19571b4ebd2SPeter Brune /* Fit points with quadratic */ 1968a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*lambda*initslope); 19771b4ebd2SPeter Brune lambdaprev = lambda; 1988a430afbSPeter Brune gprev = g; 19971b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 20071b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 20171b4ebd2SPeter Brune else lambda = lambdatemp; 20271b4ebd2SPeter Brune 20371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2049bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2059bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2069bd66eb0SPeter Brune } 20771b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 20871b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 20971b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 210f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 21171b4ebd2SPeter Brune PetscFunctionReturn(0); 21271b4ebd2SPeter Brune } 213bd4a8a71SBarry Smith if (objective) { 2148a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 2158a430afbSPeter Brune } else { 21671b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 2176188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 2186188f407SPeter Brune if (domainerror) { 21971b4ebd2SPeter Brune PetscFunctionReturn(0); 22071b4ebd2SPeter Brune } 2219bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2229bd66eb0SPeter Brune gnorm = fnorm; 2239bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2249bd66eb0SPeter Brune } else { 22571b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2268a430afbSPeter Brune g = gnorm*gnorm; 2279bd66eb0SPeter Brune } 2288a430afbSPeter Brune } 2298a430afbSPeter Brune if (PetscIsInfOrNanReal(g)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 23071b4ebd2SPeter Brune if (monitor) { 231ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 232bd4a8a71SBarry Smith if (!objective) { 233c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); 2348a430afbSPeter Brune } else { 2358a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr); 2368a430afbSPeter Brune } 237ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 23871b4ebd2SPeter Brune } 2398a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */ 24071b4ebd2SPeter Brune if (monitor) { 241ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 243ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24471b4ebd2SPeter Brune } 24571b4ebd2SPeter Brune } else { 24671b4ebd2SPeter Brune /* Fit points with cubic */ 24771b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 24871b4ebd2SPeter Brune if (lambda <= minlambda) { 24971b4ebd2SPeter Brune if (monitor) { 250ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 252bd4a8a71SBarry Smith if (!objective) { 25371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 25471b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 255c69d1a72SBarry Smith (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2568a430afbSPeter Brune } else { 2578a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 2588a430afbSPeter Brune " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 2598a430afbSPeter Brune (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2608a430afbSPeter Brune } 261ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 26271b4ebd2SPeter Brune } 263f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 26471b4ebd2SPeter Brune PetscFunctionReturn(0); 26571b4ebd2SPeter Brune } 266b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 2678a430afbSPeter Brune t1 = .5*(g - f) - lambda*initslope; 2688a430afbSPeter Brune t2 = .5*(gprev - f) - lambdaprev*initslope; 26971b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27071b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27171b4ebd2SPeter Brune d = b*b - 3*a*initslope; 27271b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 27371b4ebd2SPeter Brune if (a == 0.0) { 27471b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 27571b4ebd2SPeter Brune } else { 27671b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 27771b4ebd2SPeter Brune } 278b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 2798a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*initslope); 280*f23aa3ddSBarry Smith } else SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt"); 28171b4ebd2SPeter Brune lambdaprev = lambda; 2828a430afbSPeter Brune gprev = g; 28371b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 28471b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 28571b4ebd2SPeter Brune else lambda = lambdatemp; 28671b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 287b7b2e573SPeter Brune if (linesearch->ops->viproject) { 288b7b2e573SPeter Brune ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); 289b7b2e573SPeter Brune } 29071b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 29171b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 292bd4a8a71SBarry Smith if (!objective) { 29371b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 294c69d1a72SBarry Smith (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 2958a430afbSPeter Brune } 296f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 29771b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 29871b4ebd2SPeter Brune PetscFunctionReturn(0); 29971b4ebd2SPeter Brune } 300bd4a8a71SBarry Smith if (objective) { 3018a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 3028a430afbSPeter Brune } else { 30371b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 304c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 305c427506bSPeter Brune if (domainerror) { 306f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 30771b4ebd2SPeter Brune PetscFunctionReturn(0); 30871b4ebd2SPeter Brune } 3099bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3109bd66eb0SPeter Brune gnorm = fnorm; 3119bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3129bd66eb0SPeter Brune } else { 31371b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3148a430afbSPeter Brune g = gnorm*gnorm; 3158a430afbSPeter Brune } 3169bd66eb0SPeter Brune } 31771b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3188a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */ 31971b4ebd2SPeter Brune if (monitor) { 320ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 321bd4a8a71SBarry Smith if (!objective) { 322b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 323c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 32471b4ebd2SPeter Brune } else { 325c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 32671b4ebd2SPeter Brune } 327ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3288a430afbSPeter Brune } else { 3298a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3308a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3318a430afbSPeter Brune } else { 3328a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3338a430afbSPeter Brune } 3348a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3358a430afbSPeter Brune } 33671b4ebd2SPeter Brune } 33771b4ebd2SPeter Brune break; 33871b4ebd2SPeter Brune } else { 33971b4ebd2SPeter Brune if (monitor) { 340ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 341bd4a8a71SBarry Smith if (!objective) { 342b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 343c69d1a72SBarry 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); 34471b4ebd2SPeter Brune } else { 345c69d1a72SBarry 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); 34671b4ebd2SPeter Brune } 347ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3488a430afbSPeter Brune } else { 3498a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3508a430afbSPeter 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); 3518a430afbSPeter Brune } else { 3528a430afbSPeter 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); 3538a430afbSPeter Brune } 3548a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3558a430afbSPeter Brune } 35671b4ebd2SPeter Brune } 35771b4ebd2SPeter Brune } 35871b4ebd2SPeter Brune } 35971b4ebd2SPeter Brune } 36071b4ebd2SPeter Brune } 36171b4ebd2SPeter Brune 36271b4ebd2SPeter Brune /* postcheck */ 3637b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 36471b4ebd2SPeter Brune if (changed_y) { 36571b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 3669bd66eb0SPeter Brune if (linesearch->ops->viproject) { 3679bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 3689bd66eb0SPeter Brune } 36971b4ebd2SPeter Brune } 370bd4a8a71SBarry Smith if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ 371c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 37271b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 37371b4ebd2SPeter Brune if (domainerror) { 374f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 37571b4ebd2SPeter Brune PetscFunctionReturn(0); 37671b4ebd2SPeter Brune } 3779bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3789bd66eb0SPeter Brune gnorm = fnorm; 3799bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3809bd66eb0SPeter Brune } else { 3819bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3829bd66eb0SPeter Brune } 3839bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 384c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3859bd66eb0SPeter Brune 386c427506bSPeter Brune } 38771b4ebd2SPeter Brune 38871b4ebd2SPeter Brune /* copy the solution over */ 38971b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 390c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 391c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 392f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 393f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 39471b4ebd2SPeter Brune PetscFunctionReturn(0); 39571b4ebd2SPeter Brune } 39671b4ebd2SPeter Brune 39771b4ebd2SPeter Brune #undef __FUNCT__ 3987f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 3997f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 4007f1410a3SPeter Brune { 4017f1410a3SPeter Brune PetscErrorCode ierr; 4027f1410a3SPeter Brune PetscBool iascii; 4037f1410a3SPeter Brune SNESLineSearch_BT *bt; 404075cc632SBarry Smith 4057f1410a3SPeter Brune PetscFunctionBegin; 406251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4077f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 4087f1410a3SPeter Brune if (iascii) { 409b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 4107f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 411b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 4127f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 4137f1410a3SPeter Brune } 414007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 4157f1410a3SPeter Brune } 4167f1410a3SPeter Brune PetscFunctionReturn(0); 4177f1410a3SPeter Brune } 4187f1410a3SPeter Brune 4197f1410a3SPeter Brune 4207f1410a3SPeter Brune #undef __FUNCT__ 421f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 422f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 423731c067cSJed Brown { 424731c067cSJed Brown PetscErrorCode ierr; 42571b4ebd2SPeter Brune 42671b4ebd2SPeter Brune PetscFunctionBegin; 427731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 42871b4ebd2SPeter Brune PetscFunctionReturn(0); 42971b4ebd2SPeter Brune } 43071b4ebd2SPeter Brune 43171b4ebd2SPeter Brune 43271b4ebd2SPeter Brune #undef __FUNCT__ 433f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 434f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 43571b4ebd2SPeter Brune { 43671b4ebd2SPeter Brune 43771b4ebd2SPeter Brune PetscErrorCode ierr; 438f1c6b773SPeter Brune SNESLineSearch_BT *bt; 43971b4ebd2SPeter Brune 4406e111a19SKarl Rupp PetscFunctionBegin; 441f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 44271b4ebd2SPeter Brune 443f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 4447a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 44571b4ebd2SPeter Brune 44671b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 44771b4ebd2SPeter Brune PetscFunctionReturn(0); 44871b4ebd2SPeter Brune } 44971b4ebd2SPeter Brune 45071b4ebd2SPeter Brune 45171b4ebd2SPeter Brune #undef __FUNCT__ 452f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 453954494b2SJed Brown /*MC 454db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 455954494b2SJed Brown 456db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 457db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 458db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 459954494b2SJed Brown 460cd7522eaSPeter Brune Options Database Keys: 461cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 462db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 463db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 464db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 465cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 466cd7522eaSPeter Brune 467954494b2SJed Brown Level: advanced 468954494b2SJed Brown 469db609ea7SPeter Brune Notes: 470db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 471db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 472db609ea7SPeter Brune 473f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 474954494b2SJed Brown 475f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 476954494b2SJed Brown M*/ 477f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 47871b4ebd2SPeter Brune { 47971b4ebd2SPeter Brune 480f1c6b773SPeter Brune SNESLineSearch_BT *bt; 48171b4ebd2SPeter Brune PetscErrorCode ierr; 48271b4ebd2SPeter Brune 48371b4ebd2SPeter Brune PetscFunctionBegin; 484f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 485f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 486f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 48771b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 4887f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 48971b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 49071b4ebd2SPeter Brune 491f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 49271b4ebd2SPeter Brune linesearch->data = (void *)bt; 49371b4ebd2SPeter Brune linesearch->max_its = 40; 494b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 49571b4ebd2SPeter Brune bt->alpha = 1e-4; 49671b4ebd2SPeter Brune 49771b4ebd2SPeter Brune PetscFunctionReturn(0); 49871b4ebd2SPeter Brune } 499