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); 9271b4ebd2SPeter Brune if (!jac) { 93f1c6b773SPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); 9471b4ebd2SPeter Brune } 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); 11871b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),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 } 12471b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 12571b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 12671b4ebd2SPeter Brune ierr = VecDot(F,W,&cinitslope);CHKERRQ(ierr); 12771b4ebd2SPeter Brune initslope = PetscRealPart(cinitslope); 12871b4ebd2SPeter Brune #else 12971b4ebd2SPeter Brune ierr = VecDot(F,W,&initslope);CHKERRQ(ierr); 13071b4ebd2SPeter Brune #endif 13171b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 13271b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 13371b4ebd2SPeter Brune 134c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1359bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1369bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1379bd66eb0SPeter Brune } 13871b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 13971b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 14071b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 141f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 14271b4ebd2SPeter Brune PetscFunctionReturn(0); 14371b4ebd2SPeter Brune } 14471b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 145c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 146c427506bSPeter Brune if (domainerror) { 147f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 14871b4ebd2SPeter Brune PetscFunctionReturn(0); 14971b4ebd2SPeter Brune } 1509bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1519bd66eb0SPeter Brune gnorm = fnorm; 1529bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1539bd66eb0SPeter Brune } else { 15471b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1559bd66eb0SPeter Brune } 1569bd66eb0SPeter Brune 15771b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 15871b4ebd2SPeter Brune ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 1593b288129SPeter Brune if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope || ynorm < stol*xnorm) { /* Sufficient reduction or step tolerance convergence */ 16071b4ebd2SPeter Brune if (monitor) { 161ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 163ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 16471b4ebd2SPeter Brune } 16571b4ebd2SPeter Brune } else { 16671b4ebd2SPeter Brune /* Fit points with quadratic */ 16771b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope); 16871b4ebd2SPeter Brune lambdaprev = lambda; 16971b4ebd2SPeter Brune gnormprev = gnorm; 17071b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 17171b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 17271b4ebd2SPeter Brune else lambda = lambdatemp; 17371b4ebd2SPeter Brune 17471b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 1759bd66eb0SPeter Brune if (linesearch->ops->viproject) { 1769bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 1779bd66eb0SPeter Brune } 17871b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 17971b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 18071b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 181f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 18271b4ebd2SPeter Brune PetscFunctionReturn(0); 18371b4ebd2SPeter Brune } 18471b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 1856188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 1866188f407SPeter Brune if (domainerror) { 18771b4ebd2SPeter Brune PetscFunctionReturn(0); 18871b4ebd2SPeter Brune } 1899bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 1909bd66eb0SPeter Brune gnorm = fnorm; 1919bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 1929bd66eb0SPeter Brune } else { 19371b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 1949bd66eb0SPeter Brune } 19571b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 19671b4ebd2SPeter Brune if (monitor) { 197ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 19871b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr); 199ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 20071b4ebd2SPeter Brune } 20171b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 20271b4ebd2SPeter Brune if (monitor) { 203ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 20471b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 205ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 20671b4ebd2SPeter Brune } 20771b4ebd2SPeter Brune } else { 20871b4ebd2SPeter Brune /* Fit points with cubic */ 20971b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 21071b4ebd2SPeter Brune if (lambda <= minlambda) { 21171b4ebd2SPeter Brune if (monitor) { 212ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 21371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 21471b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 21571b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 21671b4ebd2SPeter Brune fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr); 217ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 21871b4ebd2SPeter Brune } 219f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 22071b4ebd2SPeter Brune PetscFunctionReturn(0); 22171b4ebd2SPeter Brune } 222b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 22371b4ebd2SPeter Brune t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope; 22471b4ebd2SPeter Brune t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 22571b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 22671b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 22771b4ebd2SPeter Brune d = b*b - 3*a*initslope; 22871b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 22971b4ebd2SPeter Brune if (a == 0.0) { 23071b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 23171b4ebd2SPeter Brune } else { 23271b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 23371b4ebd2SPeter Brune } 234b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 23571b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope); 23659405d5eSPeter Brune } else { 23759405d5eSPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt"); 23871b4ebd2SPeter Brune } 23971b4ebd2SPeter Brune lambdaprev = lambda; 24071b4ebd2SPeter Brune gnormprev = gnorm; 24171b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 24271b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 24371b4ebd2SPeter Brune else lambda = lambdatemp; 24471b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 245*b7b2e573SPeter Brune if (linesearch->ops->viproject) { 246*b7b2e573SPeter Brune ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); 247*b7b2e573SPeter Brune } 24871b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 24971b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 25071b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 25171b4ebd2SPeter Brune fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr); 252f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 25371b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 25471b4ebd2SPeter Brune PetscFunctionReturn(0); 25571b4ebd2SPeter Brune } 25671b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 257c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 258c427506bSPeter Brune if (domainerror) { 259f1c6b773SPeter Brune ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 26071b4ebd2SPeter Brune PetscFunctionReturn(0); 26171b4ebd2SPeter Brune } 2629bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 2639bd66eb0SPeter Brune gnorm = fnorm; 2649bd66eb0SPeter Brune ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); 2659bd66eb0SPeter Brune } else { 26671b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 2679bd66eb0SPeter Brune } 26871b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 26971b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 27071b4ebd2SPeter Brune if (monitor) { 271ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 272b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 27371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 27471b4ebd2SPeter Brune } else { 27571b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 27671b4ebd2SPeter Brune } 277ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 27871b4ebd2SPeter Brune } 27971b4ebd2SPeter Brune break; 28071b4ebd2SPeter Brune } else { 28171b4ebd2SPeter Brune if (monitor) { 282ea5d4fccSPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 283b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 28471b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 28571b4ebd2SPeter Brune } else { 28671b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 28771b4ebd2SPeter Brune } 288ea5d4fccSPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); 28971b4ebd2SPeter Brune } 29071b4ebd2SPeter Brune } 29171b4ebd2SPeter Brune } 29271b4ebd2SPeter Brune } 29371b4ebd2SPeter Brune } 29471b4ebd2SPeter Brune 29571b4ebd2SPeter Brune /* postcheck */ 2967b1df9c1SPeter Brune ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 29771b4ebd2SPeter Brune if (changed_y) { 29871b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 2999bd66eb0SPeter Brune if (linesearch->ops->viproject) { 3009bd66eb0SPeter Brune ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); 3019bd66eb0SPeter Brune } 30271b4ebd2SPeter Brune } 303c427506bSPeter Brune if (changed_y || changed_w) { /* recompute the function if the step has changed */ 304c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 30571b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 30671b4ebd2SPeter 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 { 3149bd66eb0SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 3159bd66eb0SPeter Brune } 3169bd66eb0SPeter Brune ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 317c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 3189bd66eb0SPeter Brune 319c427506bSPeter Brune } 32071b4ebd2SPeter Brune 32171b4ebd2SPeter Brune /* copy the solution over */ 32271b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 323c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 324c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 325f1c6b773SPeter Brune ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 326f1c6b773SPeter Brune ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 32771b4ebd2SPeter Brune PetscFunctionReturn(0); 32871b4ebd2SPeter Brune } 32971b4ebd2SPeter Brune 33071b4ebd2SPeter Brune #undef __FUNCT__ 3317f1410a3SPeter Brune #define __FUNCT__ "SNESLineSearchView_BT" 3327f1410a3SPeter Brune PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer) 3337f1410a3SPeter Brune { 3347f1410a3SPeter Brune PetscErrorCode ierr; 3357f1410a3SPeter Brune PetscBool iascii; 3367f1410a3SPeter Brune SNESLineSearch_BT *bt; 3377f1410a3SPeter Brune PetscFunctionBegin; 338251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3397f1410a3SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 3407f1410a3SPeter Brune if (iascii) { 341b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { 3427f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: cubic\n");CHKERRQ(ierr); 343b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { 3447f1410a3SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " interpolation: quadratic\n");CHKERRQ(ierr); 3457f1410a3SPeter Brune } 346007b6e36SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " alpha=%e\n", bt->alpha);CHKERRQ(ierr); 3477f1410a3SPeter Brune } 3487f1410a3SPeter Brune PetscFunctionReturn(0); 3497f1410a3SPeter Brune } 3507f1410a3SPeter Brune 3517f1410a3SPeter Brune 3527f1410a3SPeter Brune #undef __FUNCT__ 353f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchDestroy_BT" 354f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch) 355731c067cSJed Brown { 356731c067cSJed Brown PetscErrorCode ierr; 35771b4ebd2SPeter Brune 35871b4ebd2SPeter Brune PetscFunctionBegin; 359731c067cSJed Brown ierr = PetscFree(linesearch->data);CHKERRQ(ierr); 36071b4ebd2SPeter Brune PetscFunctionReturn(0); 36171b4ebd2SPeter Brune } 36271b4ebd2SPeter Brune 36371b4ebd2SPeter Brune 36471b4ebd2SPeter Brune #undef __FUNCT__ 365f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchSetFromOptions_BT" 366f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch) 36771b4ebd2SPeter Brune { 36871b4ebd2SPeter Brune 36971b4ebd2SPeter Brune PetscErrorCode ierr; 370f1c6b773SPeter Brune SNESLineSearch_BT *bt; 37171b4ebd2SPeter Brune PetscFunctionBegin; 37271b4ebd2SPeter Brune 373f1c6b773SPeter Brune bt = (SNESLineSearch_BT*)linesearch->data; 37471b4ebd2SPeter Brune 375f1c6b773SPeter Brune ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr); 3767a35526eSPeter Brune ierr = PetscOptionsReal("-snes_linesearch_alpha", "Descent tolerance", "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 37771b4ebd2SPeter Brune 37871b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 37971b4ebd2SPeter Brune PetscFunctionReturn(0); 38071b4ebd2SPeter Brune } 38171b4ebd2SPeter Brune 38271b4ebd2SPeter Brune 38371b4ebd2SPeter Brune #undef __FUNCT__ 384f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_BT" 385954494b2SJed Brown /*MC 386db609ea7SPeter Brune SNESLINESEARCHBT - Backtracking line search. 387954494b2SJed Brown 388db609ea7SPeter Brune This line search finds the minimum of a polynomial fitting of the L2 norm of the 389db609ea7SPeter Brune function. If this fit does not satisfy the conditions for progress, the interval shrinks 390db609ea7SPeter Brune and the fit is reattempted at most max_it times or until lambda is below minlambda. 391954494b2SJed Brown 392cd7522eaSPeter Brune Options Database Keys: 393cd7522eaSPeter Brune + -snes_linesearch_alpha<1e-4> - slope descent parameter 394db609ea7SPeter Brune . -snes_linesearch_damping<1.0> - initial step length 395db609ea7SPeter Brune . -snes_linesearch_max_it<40> - maximum number of shrinking step 396db609ea7SPeter Brune . -snes_linesearch_minlambda<1e-12> - minimum step length allowed 397cd7522eaSPeter Brune - -snes_linesearch_order<cubic,quadratic> - order of the approximation 398cd7522eaSPeter Brune 399954494b2SJed Brown Level: advanced 400954494b2SJed Brown 401db609ea7SPeter Brune Notes: 402db609ea7SPeter Brune This line search is taken from "Numerical Methods for Unconstrained 403db609ea7SPeter Brune Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 404db609ea7SPeter Brune 405f1c6b773SPeter Brune .keywords: SNES, SNESLineSearch, damping 406954494b2SJed Brown 407f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType() 408954494b2SJed Brown M*/ 409f1c6b773SPeter Brune PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch) 41071b4ebd2SPeter Brune { 41171b4ebd2SPeter Brune 412f1c6b773SPeter Brune SNESLineSearch_BT *bt; 41371b4ebd2SPeter Brune PetscErrorCode ierr; 41471b4ebd2SPeter Brune 41571b4ebd2SPeter Brune PetscFunctionBegin; 416f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_BT; 417f1c6b773SPeter Brune linesearch->ops->destroy = SNESLineSearchDestroy_BT; 418f1c6b773SPeter Brune linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT; 41971b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 4207f1410a3SPeter Brune linesearch->ops->view = SNESLineSearchView_BT; 42171b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 42271b4ebd2SPeter Brune 423f1c6b773SPeter Brune ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr); 42471b4ebd2SPeter Brune linesearch->data = (void *)bt; 42571b4ebd2SPeter Brune linesearch->max_its = 40; 426b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_CUBIC; 42771b4ebd2SPeter Brune bt->alpha = 1e-4; 42871b4ebd2SPeter Brune 42971b4ebd2SPeter Brune PetscFunctionReturn(0); 43071b4ebd2SPeter Brune } 431