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; 242f4102e2SPeter Brune PetscFunctionBegin; 252f4102e2SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 262f4102e2SPeter Brune bt = (SNESLineSearch_BT *)linesearch->data; 272f4102e2SPeter Brune bt->alpha = alpha; 282f4102e2SPeter Brune PetscFunctionReturn(0); 292f4102e2SPeter Brune } 302f4102e2SPeter Brune 312f4102e2SPeter Brune 322f4102e2SPeter Brune #undef __FUNCT__ 332f4102e2SPeter Brune #define __FUNCT__ "SNESLineSearchBTGetAlpha" 342f4102e2SPeter Brune /*@ 352f4102e2SPeter Brune SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant. 362f4102e2SPeter Brune 372f4102e2SPeter Brune Input Parameters: 382f4102e2SPeter Brune . linesearch - linesearch context 398e557f58SPeter Brune 408e557f58SPeter Brune Output Parameters: 412f4102e2SPeter Brune . alpha - The descent parameter 422f4102e2SPeter Brune 432f4102e2SPeter Brune Level: intermediate 442f4102e2SPeter Brune 451a4f838cSPeter Brune .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT 462f4102e2SPeter Brune @*/ 472f4102e2SPeter Brune PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha) 482f4102e2SPeter Brune { 492f4102e2SPeter Brune SNESLineSearch_BT *bt; 502f4102e2SPeter Brune PetscFunctionBegin; 512f4102e2SPeter Brune PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 522f4102e2SPeter Brune bt = (SNESLineSearch_BT *)linesearch->data; 532f4102e2SPeter Brune *alpha = bt->alpha; 542f4102e2SPeter Brune PetscFunctionReturn(0); 552f4102e2SPeter Brune } 562f4102e2SPeter Brune 572f4102e2SPeter Brune #undef __FUNCT__ 58f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_BT" 59f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) 6071b4ebd2SPeter Brune { 6171b4ebd2SPeter Brune PetscBool changed_y,changed_w; 6271b4ebd2SPeter Brune PetscErrorCode ierr; 6371b4ebd2SPeter Brune Vec X,F,Y,W,G; 6471b4ebd2SPeter Brune SNES snes; 658a430afbSPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm; 663b288129SPeter Brune PetscReal lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol; 67dfcd3828SPeter Brune PetscReal t1,t2,a,b,d; 688a430afbSPeter Brune PetscReal f; 698a430afbSPeter Brune PetscReal g,gprev; 7071b4ebd2SPeter Brune PetscBool domainerror; 7171b4ebd2SPeter Brune PetscViewer monitor; 7271b4ebd2SPeter Brune PetscInt max_its,count; 73f1c6b773SPeter Brune SNESLineSearch_BT *bt; 7471b4ebd2SPeter Brune Mat jac; 75*bd4a8a71SBarry Smith PetscErrorCode (*objective)(SNES,Vec,PetscReal *,void*); 7671b4ebd2SPeter Brune 7771b4ebd2SPeter Brune PetscFunctionBegin; 7871b4ebd2SPeter Brune 79f1c6b773SPeter Brune ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr); 80f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 81f1c6b773SPeter Brune ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 82f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 83f1c6b773SPeter Brune ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 8422d28d08SBarry Smith ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,PETSC_NULL,PETSC_NULL,PETSC_NULL,&max_its);CHKERRQ(ierr); 853b288129SPeter Brune ierr = SNESGetTolerances(snes,PETSC_NULL,PETSC_NULL,&stol,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 86*bd4a8a71SBarry Smith ierr = SNESGetObjective(snes,&objective,PETSC_NULL);CHKERRQ(ierr); 87f1c6b773SPeter Brune bt = (SNESLineSearch_BT *)linesearch->data; 8871b4ebd2SPeter Brune 8971b4ebd2SPeter Brune alpha = bt->alpha; 9071b4ebd2SPeter Brune 9171b4ebd2SPeter Brune ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 928a430afbSPeter Brune 93*bd4a8a71SBarry Smith if (!jac && !objective) SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 94c69d1a72SBarry Smith 9571b4ebd2SPeter Brune /* precheck */ 967b1df9c1SPeter Brune ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); 97f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 9871b4ebd2SPeter Brune 993b288129SPeter Brune ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr); 1003b288129SPeter Brune ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 1013b288129SPeter Brune ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr); 1023b288129SPeter Brune ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 1033b288129SPeter Brune 10471b4ebd2SPeter Brune if (ynorm == 0.0) { 10571b4ebd2SPeter Brune if (monitor) { 106ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 10771b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 108ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 10971b4ebd2SPeter Brune } 11071b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 11171b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 112f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 11371b4ebd2SPeter Brune PetscFunctionReturn(0); 11471b4ebd2SPeter Brune } 11571b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 11671b4ebd2SPeter Brune if (monitor) { 117ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 118c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr); 119ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 12071b4ebd2SPeter Brune } 12171b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 12271b4ebd2SPeter Brune ynorm = maxstep; 12371b4ebd2SPeter Brune } 1248a430afbSPeter Brune 1258a430afbSPeter Brune /* if the SNES has an objective set, use that instead of the function value */ 126*bd4a8a71SBarry Smith if (objective) { 1278a430afbSPeter Brune ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr); 1288a430afbSPeter Brune } else { 1298a430afbSPeter Brune f = fnorm*fnorm; 1308a430afbSPeter Brune } 1318a430afbSPeter Brune 1328a430afbSPeter Brune /* compute the initial slope */ 133*bd4a8a71SBarry Smith if (objective) { 1348a430afbSPeter Brune /* slope comes from the function (assumed to be the gradient of the objective */ 13567392de3SBarry Smith ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr); 1368a430afbSPeter Brune } else { 1378a430afbSPeter Brune /* slope comes from the normal equations */ 13871b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 13967392de3SBarry Smith ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr); 14071b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 14171b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 1428a430afbSPeter Brune } 14371b4ebd2SPeter Brune 144c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1459bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1469bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1479bd66eb0SPeter Brune } 14871b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 14971b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 15071b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 151f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 15271b4ebd2SPeter Brune PetscFunctionReturn(0); 15371b4ebd2SPeter Brune } 1548a430afbSPeter Brune 155*bd4a8a71SBarry Smith if (objective) { 1568a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 1578a430afbSPeter Brune } else { 15871b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 159c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 160c427506bSPeter Brune if (domainerror) { 161f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 16271b4ebd2SPeter Brune PetscFunctionReturn(0); 16371b4ebd2SPeter Brune } 1649bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1659bd66eb0SPeter Brune gnorm = fnorm; 1669bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1679bd66eb0SPeter Brune } else { 16871b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1698a430afbSPeter Brune } 17079e7e81bSJed Brown g = PetscSqr(gnorm); 1719bd66eb0SPeter Brune } 1729bd66eb0SPeter Brune 1738a430afbSPeter Brune if (PetscIsInfOrNanReal(g)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 174*bd4a8a71SBarry Smith if (!objective) { 175c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1768a430afbSPeter Brune } 1778a430afbSPeter Brune if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */ 17871b4ebd2SPeter Brune if (monitor) { 179ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 180*bd4a8a71SBarry Smith if (!objective) { 181c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1828a430afbSPeter Brune } else { 1838a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr); 1848a430afbSPeter Brune } 185ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 18671b4ebd2SPeter Brune } 18771b4ebd2SPeter Brune } else { 188c21ba15cSPeter Brune /* Since the full step didn't work and the step is tiny, quit */ 189c21ba15cSPeter Brune if (stol*xnorm > ynorm) { 190c21ba15cSPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 191c21ba15cSPeter Brune ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr); 192c21ba15cSPeter Brune PetscFunctionReturn(0); 193c21ba15cSPeter Brune } 19471b4ebd2SPeter Brune /* Fit points with quadratic */ 1958a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*lambda*initslope); 19671b4ebd2SPeter Brune lambdaprev = lambda; 1978a430afbSPeter Brune gprev = g; 19871b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 19971b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 20071b4ebd2SPeter Brune else lambda = lambdatemp; 20171b4ebd2SPeter Brune 20271b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2039bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2049bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2059bd66eb0SPeter Brune } 20671b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 20771b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 20871b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 209f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 21071b4ebd2SPeter Brune PetscFunctionReturn(0); 21171b4ebd2SPeter Brune } 212*bd4a8a71SBarry Smith if (objective) { 2138a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 2148a430afbSPeter Brune } else { 21571b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 2166188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 2176188f407SPeter Brune if (domainerror) { 21871b4ebd2SPeter Brune PetscFunctionReturn(0); 21971b4ebd2SPeter Brune } 2209bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2219bd66eb0SPeter Brune gnorm = fnorm; 2229bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2239bd66eb0SPeter Brune } else { 22471b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2258a430afbSPeter Brune g = gnorm*gnorm; 2269bd66eb0SPeter Brune } 2278a430afbSPeter Brune } 2288a430afbSPeter Brune if (PetscIsInfOrNanReal(g)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 22971b4ebd2SPeter Brune if (monitor) { 230ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 231*bd4a8a71SBarry Smith if (!objective) { 232c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); 2338a430afbSPeter Brune } else { 2348a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr); 2358a430afbSPeter Brune } 236ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 23771b4ebd2SPeter Brune } 2388a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */ 23971b4ebd2SPeter Brune if (monitor) { 240ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 242ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 24371b4ebd2SPeter Brune } 24471b4ebd2SPeter Brune } else { 24571b4ebd2SPeter Brune /* Fit points with cubic */ 24671b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 24771b4ebd2SPeter Brune if (lambda <= minlambda) { 24871b4ebd2SPeter Brune if (monitor) { 249ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 25071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 251*bd4a8a71SBarry Smith if (!objective) { 25271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 25371b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 254c69d1a72SBarry Smith (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2558a430afbSPeter Brune } else { 2568a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 2578a430afbSPeter Brune " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 2588a430afbSPeter Brune (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 2598a430afbSPeter Brune } 260ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 26171b4ebd2SPeter Brune } 262f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 26371b4ebd2SPeter Brune PetscFunctionReturn(0); 26471b4ebd2SPeter Brune } 265b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 2668a430afbSPeter Brune t1 = .5*(g - f) - lambda*initslope; 2678a430afbSPeter Brune t2 = .5*(gprev - f) - lambdaprev*initslope; 26871b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 26971b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 27071b4ebd2SPeter Brune d = b*b - 3*a*initslope; 27171b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 27271b4ebd2SPeter Brune if (a == 0.0) { 27371b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 27471b4ebd2SPeter Brune } else { 27571b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 27671b4ebd2SPeter Brune } 277b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 2788a430afbSPeter Brune lambdatemp = -initslope/(g - f - 2.0*initslope); 27959405d5eSPeter Brune } else { 28059405d5eSPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt"); 28171b4ebd2SPeter Brune } 28271b4ebd2SPeter Brune lambdaprev = lambda; 2838a430afbSPeter Brune gprev = g; 28471b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 28571b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 28671b4ebd2SPeter Brune else lambda = lambdatemp; 28771b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 288b7b2e573SPeter Brune if (linesearch->ops->viproject) { 289b7b2e573SPeter Brune ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); 290b7b2e573SPeter Brune } 29171b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 29271b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 293*bd4a8a71SBarry Smith if (!objective) { 29471b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 295c69d1a72SBarry Smith (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 2968a430afbSPeter Brune } 297f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 29871b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 29971b4ebd2SPeter Brune PetscFunctionReturn(0); 30071b4ebd2SPeter Brune } 301*bd4a8a71SBarry Smith if (objective) { 3028a430afbSPeter Brune ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); 3038a430afbSPeter Brune } else { 30471b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 305c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 306c427506bSPeter Brune if (domainerror) { 307f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 30871b4ebd2SPeter Brune PetscFunctionReturn(0); 30971b4ebd2SPeter Brune } 3109bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3119bd66eb0SPeter Brune gnorm = fnorm; 3129bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3139bd66eb0SPeter Brune } else { 31471b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3158a430afbSPeter Brune g = gnorm*gnorm; 3168a430afbSPeter Brune } 3179bd66eb0SPeter Brune } 31871b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3198a430afbSPeter Brune if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */ 32071b4ebd2SPeter Brune if (monitor) { 321ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 322*bd4a8a71SBarry Smith if (!objective) { 323b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 324c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 32571b4ebd2SPeter Brune } else { 326c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 32771b4ebd2SPeter Brune } 328ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3298a430afbSPeter Brune } else { 3308a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3318a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3328a430afbSPeter Brune } else { 3338a430afbSPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); 3348a430afbSPeter Brune } 3358a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3368a430afbSPeter Brune } 33771b4ebd2SPeter Brune } 33871b4ebd2SPeter Brune break; 33971b4ebd2SPeter Brune } else { 34071b4ebd2SPeter Brune if (monitor) { 341ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 342*bd4a8a71SBarry Smith if (!objective) { 343b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 344c69d1a72SBarry 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); 34571b4ebd2SPeter Brune } else { 346c69d1a72SBarry 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); 34771b4ebd2SPeter Brune } 348ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3498a430afbSPeter Brune } else { 3508a430afbSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3518a430afbSPeter 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); 3528a430afbSPeter Brune } else { 3538a430afbSPeter 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); 3548a430afbSPeter Brune } 3558a430afbSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 3568a430afbSPeter Brune } 35771b4ebd2SPeter Brune } 35871b4ebd2SPeter Brune } 35971b4ebd2SPeter Brune } 36071b4ebd2SPeter Brune } 36171b4ebd2SPeter Brune } 36271b4ebd2SPeter Brune 36371b4ebd2SPeter Brune /* postcheck */ 3647b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 36571b4ebd2SPeter Brune if (changed_y) { 36671b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 3679bd66eb0SPeter Brune if (linesearch->ops->viproject) { 3689bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 3699bd66eb0SPeter Brune } 37071b4ebd2SPeter Brune } 371*bd4a8a71SBarry Smith if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ 372c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 37371b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 37471b4ebd2SPeter Brune if (domainerror) { 375f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 37671b4ebd2SPeter Brune PetscFunctionReturn(0); 37771b4ebd2SPeter Brune } 3789bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3799bd66eb0SPeter Brune gnorm = fnorm; 3809bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3819bd66eb0SPeter Brune } else { 3829bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3839bd66eb0SPeter Brune } 3849bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 385c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3869bd66eb0SPeter Brune 387c427506bSPeter Brune } 38871b4ebd2SPeter Brune 38971b4ebd2SPeter Brune /* copy the solution over */ 39071b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 391c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 392c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 393f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 394f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 39571b4ebd2SPeter Brune PetscFunctionReturn(0); 39671b4ebd2SPeter Brune } 39771b4ebd2SPeter Brune 39871b4ebd2SPeter Brune #undef __FUNCT__ 3997f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 4007f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 4017f1410a3SPeter Brune { 4027f1410a3SPeter Brune PetscErrorCode ierr; 4037f1410a3SPeter Brune PetscBool iascii; 4047f1410a3SPeter Brune SNESLineSearch_BT *bt; 405075cc632SBarry Smith 4067f1410a3SPeter Brune PetscFunctionBegin; 407251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4087f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 4097f1410a3SPeter Brune if (iascii) { 410b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 4117f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 412b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 4137f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 4147f1410a3SPeter Brune } 415007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 4167f1410a3SPeter Brune } 4177f1410a3SPeter Brune PetscFunctionReturn(0); 4187f1410a3SPeter Brune } 4197f1410a3SPeter Brune 4207f1410a3SPeter Brune 4217f1410a3SPeter Brune #undef __FUNCT__ 422f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 423f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 424731c067cSJed Brown { 425731c067cSJed Brown PetscErrorCode ierr; 42671b4ebd2SPeter Brune 42771b4ebd2SPeter Brune PetscFunctionBegin; 428731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 42971b4ebd2SPeter Brune PetscFunctionReturn(0); 43071b4ebd2SPeter Brune } 43171b4ebd2SPeter Brune 43271b4ebd2SPeter Brune 43371b4ebd2SPeter Brune #undef __FUNCT__ 434f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 435f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 43671b4ebd2SPeter Brune { 43771b4ebd2SPeter Brune 43871b4ebd2SPeter Brune PetscErrorCode ierr; 439f1c6b773SPeter Brune SNESLineSearch_BT *bt; 44071b4ebd2SPeter Brune PetscFunctionBegin; 44171b4ebd2SPeter Brune 442f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 44371b4ebd2SPeter Brune 444f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 4457a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 44671b4ebd2SPeter Brune 44771b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 44871b4ebd2SPeter Brune PetscFunctionReturn(0); 44971b4ebd2SPeter Brune } 45071b4ebd2SPeter Brune 45171b4ebd2SPeter Brune 45271b4ebd2SPeter Brune #undef __FUNCT__ 453f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 454954494b2SJed Brown /*MC 455db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 456954494b2SJed Brown 457db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 458db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 459db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 460954494b2SJed Brown 461cd7522eaSPeter Brune Options Database Keys: 462cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 463db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 464db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 465db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 466cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 467cd7522eaSPeter Brune 468954494b2SJed Brown Level: advanced 469954494b2SJed Brown 470db609ea7SPeter Brune Notes: 471db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 472db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 473db609ea7SPeter Brune 474f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 475954494b2SJed Brown 476f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 477954494b2SJed Brown M*/ 478f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 47971b4ebd2SPeter Brune { 48071b4ebd2SPeter Brune 481f1c6b773SPeter Brune SNESLineSearch_BT *bt; 48271b4ebd2SPeter Brune PetscErrorCode ierr; 48371b4ebd2SPeter Brune 48471b4ebd2SPeter Brune PetscFunctionBegin; 485f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 486f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 487f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 48871b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 4897f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 49071b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 49171b4ebd2SPeter Brune 492f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 49371b4ebd2SPeter Brune linesearch->data = (void *)bt; 49471b4ebd2SPeter Brune linesearch->max_its = 40; 495b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 49671b4ebd2SPeter Brune bt->alpha = 1e-4; 49771b4ebd2SPeter Brune 49871b4ebd2SPeter Brune PetscFunctionReturn(0); 49971b4ebd2SPeter Brune } 500