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; 6571b4ebd2SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm, gnormprev; 663b288129SPeter Brune PetscReal lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol; 67dfcd3828SPeter Brune PetscReal t1,t2,a,b,d; 6871b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 6971b4ebd2SPeter Brune PetscScalar cinitslope; 7071b4ebd2SPeter Brune #endif 7171b4ebd2SPeter Brune PetscBool domainerror; 7271b4ebd2SPeter Brune PetscViewer monitor; 7371b4ebd2SPeter Brune PetscInt max_its,count; 74f1c6b773SPeter Brune SNESLineSearch_BT *bt; 7571b4ebd2SPeter Brune Mat jac; 7671b4ebd2SPeter Brune 7771b4ebd2SPeter Brune 7871b4ebd2SPeter Brune PetscFunctionBegin; 7971b4ebd2SPeter Brune 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); 85dfcd3828SPeter Brune ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,PETSC_NULL,PETSC_NULL,PETSC_NULL,&max_its); 863b288129SPeter Brune ierr = SNESGetTolerances(snes,PETSC_NULL,PETSC_NULL,&stol,PETSC_NULL,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); 92*c69d1a72SBarry Smith if (!jac) SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 93*c69d1a72SBarry Smith 9471b4ebd2SPeter Brune /* precheck */ 957b1df9c1SPeter Brune ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); 96f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 9771b4ebd2SPeter Brune 983b288129SPeter Brune ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr); 993b288129SPeter Brune ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 1003b288129SPeter Brune ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr); 1013b288129SPeter Brune ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 1023b288129SPeter Brune 10371b4ebd2SPeter Brune if (ynorm == 0.0) { 10471b4ebd2SPeter Brune if (monitor) { 105ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 10671b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 107ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 10871b4ebd2SPeter Brune } 10971b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 11071b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 111f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 11271b4ebd2SPeter Brune PetscFunctionReturn(0); 11371b4ebd2SPeter Brune } 11471b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 11571b4ebd2SPeter Brune if (monitor) { 116ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 117*c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr); 118ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 11971b4ebd2SPeter Brune } 12071b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 12171b4ebd2SPeter Brune ynorm = maxstep; 12271b4ebd2SPeter Brune } 12371b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 12471b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 12571b4ebd2SPeter Brune ierr = VecDot(F,W,&cinitslope);CHKERRQ(ierr); 12671b4ebd2SPeter Brune initslope = PetscRealPart(cinitslope); 12771b4ebd2SPeter Brune #else 12871b4ebd2SPeter Brune ierr = VecDot(F,W,&initslope);CHKERRQ(ierr); 12971b4ebd2SPeter Brune #endif 13071b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 13171b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 13271b4ebd2SPeter Brune 133c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1349bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1359bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1369bd66eb0SPeter Brune } 13771b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 13871b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 13971b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 140f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 14171b4ebd2SPeter Brune PetscFunctionReturn(0); 14271b4ebd2SPeter Brune } 14371b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 144c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 145c427506bSPeter Brune if (domainerror) { 146f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 14771b4ebd2SPeter Brune PetscFunctionReturn(0); 14871b4ebd2SPeter Brune } 1499bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1509bd66eb0SPeter Brune gnorm = fnorm; 1519bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1529bd66eb0SPeter Brune } else { 15371b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1549bd66eb0SPeter Brune } 1559bd66eb0SPeter Brune 15671b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 157*c69d1a72SBarry Smith ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 1583b288129SPeter Brune if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope || ynorm < stol*xnorm) { /* Sufficient reduction or step tolerance convergence */ 15971b4ebd2SPeter Brune if (monitor) { 160ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 161*c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); 162ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16371b4ebd2SPeter Brune } 16471b4ebd2SPeter Brune } else { 16571b4ebd2SPeter Brune /* Fit points with quadratic */ 16671b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope); 16771b4ebd2SPeter Brune lambdaprev = lambda; 16871b4ebd2SPeter Brune gnormprev = gnorm; 16971b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 17071b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 17171b4ebd2SPeter Brune else lambda = lambdatemp; 17271b4ebd2SPeter Brune 17371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1749bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1759bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1769bd66eb0SPeter Brune } 17771b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 17871b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 17971b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 180f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 18171b4ebd2SPeter Brune PetscFunctionReturn(0); 18271b4ebd2SPeter Brune } 18371b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 1846188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1856188f407SPeter Brune if (domainerror) { 18671b4ebd2SPeter Brune PetscFunctionReturn(0); 18771b4ebd2SPeter Brune } 1889bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1899bd66eb0SPeter Brune gnorm = fnorm; 1909bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1919bd66eb0SPeter Brune } else { 19271b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1939bd66eb0SPeter Brune } 19471b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 19571b4ebd2SPeter Brune if (monitor) { 196ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 197*c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); 198ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 19971b4ebd2SPeter Brune } 20071b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 20171b4ebd2SPeter Brune if (monitor) { 202ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 20371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 204ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 20571b4ebd2SPeter Brune } 20671b4ebd2SPeter Brune } else { 20771b4ebd2SPeter Brune /* Fit points with cubic */ 20871b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 20971b4ebd2SPeter Brune if (lambda <= minlambda) { 21071b4ebd2SPeter Brune if (monitor) { 211ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 21271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 21371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 21471b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 215*c69d1a72SBarry Smith (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); 216ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 21771b4ebd2SPeter Brune } 218f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 21971b4ebd2SPeter Brune PetscFunctionReturn(0); 22071b4ebd2SPeter Brune } 221b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 22271b4ebd2SPeter Brune t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope; 22371b4ebd2SPeter Brune t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 22471b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 22571b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 22671b4ebd2SPeter Brune d = b*b - 3*a*initslope; 22771b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 22871b4ebd2SPeter Brune if (a == 0.0) { 22971b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 23071b4ebd2SPeter Brune } else { 23171b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 23271b4ebd2SPeter Brune } 233b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 23471b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope); 23559405d5eSPeter Brune } else { 23659405d5eSPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt"); 23771b4ebd2SPeter Brune } 23871b4ebd2SPeter Brune lambdaprev = lambda; 23971b4ebd2SPeter Brune gnormprev = gnorm; 24071b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 24171b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 24271b4ebd2SPeter Brune else lambda = lambdatemp; 24371b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 24471b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 24571b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 24671b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 247*c69d1a72SBarry Smith (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 248f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 24971b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 25071b4ebd2SPeter Brune PetscFunctionReturn(0); 25171b4ebd2SPeter Brune } 25271b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 253c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 254c427506bSPeter Brune if (domainerror) { 255f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 25671b4ebd2SPeter Brune PetscFunctionReturn(0); 25771b4ebd2SPeter Brune } 2589bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2599bd66eb0SPeter Brune gnorm = fnorm; 2609bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2619bd66eb0SPeter Brune } else { 26271b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2639bd66eb0SPeter Brune } 26471b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 26571b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 26671b4ebd2SPeter Brune if (monitor) { 267ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 268b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 269*c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 27071b4ebd2SPeter Brune } else { 271*c69d1a72SBarry Smith ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); 27271b4ebd2SPeter Brune } 273ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 27471b4ebd2SPeter Brune } 27571b4ebd2SPeter Brune break; 27671b4ebd2SPeter Brune } else { 27771b4ebd2SPeter Brune if (monitor) { 278ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 279b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 280*c69d1a72SBarry 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); 28171b4ebd2SPeter Brune } else { 282*c69d1a72SBarry 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); 28371b4ebd2SPeter Brune } 284ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 28571b4ebd2SPeter Brune } 28671b4ebd2SPeter Brune } 28771b4ebd2SPeter Brune } 28871b4ebd2SPeter Brune } 28971b4ebd2SPeter Brune } 29071b4ebd2SPeter Brune 29171b4ebd2SPeter Brune /* postcheck */ 2927b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 29371b4ebd2SPeter Brune if (changed_y) { 29471b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2959bd66eb0SPeter Brune if (linesearch->ops->viproject) { 2969bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 2979bd66eb0SPeter Brune } 29871b4ebd2SPeter Brune } 299c427506bSPeter Brune if (changed_y || changed_w) { /* recompute the function if the step has changed */ 300c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 30171b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 30271b4ebd2SPeter Brune if (domainerror) { 303f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 30471b4ebd2SPeter Brune PetscFunctionReturn(0); 30571b4ebd2SPeter Brune } 3069bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 3079bd66eb0SPeter Brune gnorm = fnorm; 3089bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 3099bd66eb0SPeter Brune } else { 3109bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3119bd66eb0SPeter Brune } 3129bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 313c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3149bd66eb0SPeter Brune 315c427506bSPeter Brune } 31671b4ebd2SPeter Brune 31771b4ebd2SPeter Brune /* copy the solution over */ 31871b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 319c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 320c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 321f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 322f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 32371b4ebd2SPeter Brune PetscFunctionReturn(0); 32471b4ebd2SPeter Brune } 32571b4ebd2SPeter Brune 32671b4ebd2SPeter Brune #undef __FUNCT__ 3277f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 3287f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 3297f1410a3SPeter Brune { 3307f1410a3SPeter Brune PetscErrorCode ierr; 3317f1410a3SPeter Brune PetscBool iascii; 3327f1410a3SPeter Brune SNESLineSearch_BT *bt; 3337f1410a3SPeter Brune PetscFunctionBegin; 334251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3357f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 3367f1410a3SPeter Brune if (iascii) { 337b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3387f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 339b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 3407f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 3417f1410a3SPeter Brune } 342007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 3437f1410a3SPeter Brune } 3447f1410a3SPeter Brune PetscFunctionReturn(0); 3457f1410a3SPeter Brune } 3467f1410a3SPeter Brune 3477f1410a3SPeter Brune 3487f1410a3SPeter Brune #undef __FUNCT__ 349f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 350f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 351731c067cSJed Brown { 352731c067cSJed Brown PetscErrorCode ierr; 35371b4ebd2SPeter Brune 35471b4ebd2SPeter Brune PetscFunctionBegin; 355731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 35671b4ebd2SPeter Brune PetscFunctionReturn(0); 35771b4ebd2SPeter Brune } 35871b4ebd2SPeter Brune 35971b4ebd2SPeter Brune 36071b4ebd2SPeter Brune #undef __FUNCT__ 361f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 362f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 36371b4ebd2SPeter Brune { 36471b4ebd2SPeter Brune 36571b4ebd2SPeter Brune PetscErrorCode ierr; 366f1c6b773SPeter Brune SNESLineSearch_BT *bt; 36771b4ebd2SPeter Brune PetscFunctionBegin; 36871b4ebd2SPeter Brune 369f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 37071b4ebd2SPeter Brune 371f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 3727a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 37371b4ebd2SPeter Brune 37471b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 37571b4ebd2SPeter Brune PetscFunctionReturn(0); 37671b4ebd2SPeter Brune } 37771b4ebd2SPeter Brune 37871b4ebd2SPeter Brune 37971b4ebd2SPeter Brune #undef __FUNCT__ 380f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 381954494b2SJed Brown /*MC 382db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 383954494b2SJed Brown 384db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 385db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 386db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 387954494b2SJed Brown 388cd7522eaSPeter Brune Options Database Keys: 389cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 390db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 391db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 392db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 393cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 394cd7522eaSPeter Brune 395954494b2SJed Brown Level: advanced 396954494b2SJed Brown 397db609ea7SPeter Brune Notes: 398db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 399db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 400db609ea7SPeter Brune 401f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 402954494b2SJed Brown 403f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 404954494b2SJed Brown M*/ 405f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 40671b4ebd2SPeter Brune { 40771b4ebd2SPeter Brune 408f1c6b773SPeter Brune SNESLineSearch_BT *bt; 40971b4ebd2SPeter Brune PetscErrorCode ierr; 41071b4ebd2SPeter Brune 41171b4ebd2SPeter Brune PetscFunctionBegin; 412f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 413f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 414f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 41571b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 4167f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 41771b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 41871b4ebd2SPeter Brune 419f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 42071b4ebd2SPeter Brune linesearch->data = (void *)bt; 42171b4ebd2SPeter Brune linesearch->max_its = 40; 422b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 42371b4ebd2SPeter Brune bt->alpha = 1e-4; 42471b4ebd2SPeter Brune 42571b4ebd2SPeter Brune PetscFunctionReturn(0); 42671b4ebd2SPeter Brune } 427