171b4ebd2SPeter Brune #include <private/linesearchimpl.h> 271b4ebd2SPeter Brune #include <private/snesimpl.h> 371b4ebd2SPeter Brune 4*6188f407SPeter Brune typedef enum {PETSCLINESEARCH_BT_QUADRATIC, PETSCLINESEARCH_BT_CUBIC} PetscLineSearchBTOrder; 571b4ebd2SPeter Brune 671b4ebd2SPeter Brune typedef struct { 771b4ebd2SPeter Brune PetscReal alpha; /* sufficient decrease parameter */ 8*6188f407SPeter Brune PetscLineSearchBTOrder order; 9*6188f407SPeter Brune } PetscLineSearch_BT; 1071b4ebd2SPeter Brune 1171b4ebd2SPeter Brune /*MC 12*6188f407SPeter Brune PetscLineSearchBT - Backtracking line searches. 1371b4ebd2SPeter Brune 1471b4ebd2SPeter Brune These linesearches try a polynomial fit for the L2 norm of the error 1571b4ebd2SPeter Brune using the gradient. Failing that, they step back and try again. 1671b4ebd2SPeter Brune 1771b4ebd2SPeter Brune Level: advanced 1871b4ebd2SPeter Brune 19*6188f407SPeter Brune .keywords: SNES, PetscLineSearch, damping 2071b4ebd2SPeter Brune 21*6188f407SPeter Brune .seealso: PetscLineSearchCreate(), PetscLineSearchSetType() 2271b4ebd2SPeter Brune M*/ 2371b4ebd2SPeter Brune 2471b4ebd2SPeter Brune #undef __FUNCT__ 25*6188f407SPeter Brune #define __FUNCT__ "PetscLineSearchApply_BT" 2671b4ebd2SPeter Brune 27*6188f407SPeter Brune PetscErrorCode PetscLineSearchApply_BT(PetscLineSearch linesearch) 2871b4ebd2SPeter Brune { 2971b4ebd2SPeter Brune PetscBool changed_y, changed_w; 3071b4ebd2SPeter Brune PetscErrorCode ierr; 3171b4ebd2SPeter Brune Vec X, F, Y, W, G; 3271b4ebd2SPeter Brune SNES snes; 3371b4ebd2SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm, gnormprev; 3471b4ebd2SPeter Brune PetscReal lambda, lambdatemp, lambdaprev, minlambda, maxstep, rellength, initslope, alpha; 3571b4ebd2SPeter Brune PetscReal t1, t2, a, b, d, steptol; 3671b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 3771b4ebd2SPeter Brune PetscScalar cinitslope; 3871b4ebd2SPeter Brune #endif 3971b4ebd2SPeter Brune PetscBool domainerror; 4071b4ebd2SPeter Brune PetscViewer monitor; 4171b4ebd2SPeter Brune PetscInt max_its, count; 42*6188f407SPeter Brune PetscLineSearch_BT *bt; 4371b4ebd2SPeter Brune Mat jac; 4471b4ebd2SPeter Brune 4571b4ebd2SPeter Brune 4671b4ebd2SPeter Brune PetscFunctionBegin; 4771b4ebd2SPeter Brune 48*6188f407SPeter Brune ierr = PetscLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr); 49*6188f407SPeter Brune ierr = PetscLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 50*6188f407SPeter Brune ierr = PetscLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 51*6188f407SPeter Brune ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 52*6188f407SPeter Brune ierr = PetscLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 53*6188f407SPeter Brune ierr = PetscLineSearchGetTolerances(linesearch, &steptol, &maxstep, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_its); 54*6188f407SPeter Brune bt = (PetscLineSearch_BT *)linesearch->data; 5571b4ebd2SPeter Brune 5671b4ebd2SPeter Brune alpha = bt->alpha; 5771b4ebd2SPeter Brune 5871b4ebd2SPeter Brune ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 5971b4ebd2SPeter Brune if (!jac) { 60*6188f407SPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "PetscLineSearchBT requires a Jacobian matrix"); 6171b4ebd2SPeter Brune } 6271b4ebd2SPeter Brune /* precheck */ 63*6188f407SPeter Brune ierr = PetscLineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 64*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 6571b4ebd2SPeter Brune 6671b4ebd2SPeter Brune ierr = VecNorm(Y, NORM_2, &ynorm);CHKERRQ(ierr); 6771b4ebd2SPeter Brune if (ynorm == 0.0) { 6871b4ebd2SPeter Brune if (monitor) { 6971b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 7071b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 7171b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 7271b4ebd2SPeter Brune } 7371b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 7471b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 75*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 7671b4ebd2SPeter Brune PetscFunctionReturn(0); 7771b4ebd2SPeter Brune } 7871b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 7971b4ebd2SPeter Brune if (monitor) { 8071b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 8171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);CHKERRQ(ierr); 8271b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 8371b4ebd2SPeter Brune } 8471b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 8571b4ebd2SPeter Brune ynorm = maxstep; 8671b4ebd2SPeter Brune } 8771b4ebd2SPeter Brune ierr = VecMaxPointwiseDivide(Y,X,&rellength);CHKERRQ(ierr); 8871b4ebd2SPeter Brune minlambda = steptol/rellength; 8971b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 9071b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 9171b4ebd2SPeter Brune ierr = VecDot(F,W,&cinitslope);CHKERRQ(ierr); 9271b4ebd2SPeter Brune initslope = PetscRealPart(cinitslope); 9371b4ebd2SPeter Brune #else 9471b4ebd2SPeter Brune ierr = VecDot(F,W,&initslope);CHKERRQ(ierr); 9571b4ebd2SPeter Brune #endif 9671b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 9771b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 9871b4ebd2SPeter Brune 99c427506bSPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 10071b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 10171b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 10271b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 103*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 10471b4ebd2SPeter Brune PetscFunctionReturn(0); 10571b4ebd2SPeter Brune } 10671b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 107c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 108c427506bSPeter Brune if (domainerror) { 109*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 11071b4ebd2SPeter Brune PetscFunctionReturn(0); 11171b4ebd2SPeter Brune } 11271b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 11371b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 11471b4ebd2SPeter Brune ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 11571b4ebd2SPeter Brune if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */ 11671b4ebd2SPeter Brune if (monitor) { 11771b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 11871b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 11971b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 12071b4ebd2SPeter Brune } 12171b4ebd2SPeter Brune } else { 12271b4ebd2SPeter Brune /* Fit points with quadratic */ 12371b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope); 12471b4ebd2SPeter Brune lambdaprev = lambda; 12571b4ebd2SPeter Brune gnormprev = gnorm; 12671b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 12771b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 12871b4ebd2SPeter Brune else lambda = lambdatemp; 12971b4ebd2SPeter Brune 13071b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 13171b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 13271b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 13371b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 134*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 13571b4ebd2SPeter Brune PetscFunctionReturn(0); 13671b4ebd2SPeter Brune } 13771b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 138*6188f407SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 139*6188f407SPeter Brune if (domainerror) { 14071b4ebd2SPeter Brune PetscFunctionReturn(0); 14171b4ebd2SPeter Brune } 14271b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 14371b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 14471b4ebd2SPeter Brune if (monitor) { 14571b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 14671b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr); 14771b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 14871b4ebd2SPeter Brune } 14971b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 15071b4ebd2SPeter Brune if (monitor) { 15171b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 15271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 15371b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 15471b4ebd2SPeter Brune } 15571b4ebd2SPeter Brune } else { 15671b4ebd2SPeter Brune /* Fit points with cubic */ 15771b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 15871b4ebd2SPeter Brune if (lambda <= minlambda) { 15971b4ebd2SPeter Brune if (monitor) { 16071b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 16171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 16271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 16371b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 16471b4ebd2SPeter Brune fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr); 16571b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 16671b4ebd2SPeter Brune } 167*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 16871b4ebd2SPeter Brune PetscFunctionReturn(0); 16971b4ebd2SPeter Brune } 170*6188f407SPeter Brune if (bt->order == PETSCLINESEARCH_BT_CUBIC) { 17171b4ebd2SPeter Brune t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope; 17271b4ebd2SPeter Brune t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 17371b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 17471b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 17571b4ebd2SPeter Brune d = b*b - 3*a*initslope; 17671b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 17771b4ebd2SPeter Brune if (a == 0.0) { 17871b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 17971b4ebd2SPeter Brune } else { 18071b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 18171b4ebd2SPeter Brune } 182*6188f407SPeter Brune } else if (bt->order == PETSCLINESEARCH_BT_QUADRATIC) { 18371b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope); 18471b4ebd2SPeter Brune } 18571b4ebd2SPeter Brune lambdaprev = lambda; 18671b4ebd2SPeter Brune gnormprev = gnorm; 18771b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 18871b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 18971b4ebd2SPeter Brune else lambda = lambdatemp; 19071b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 19171b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 19271b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 19371b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 19471b4ebd2SPeter Brune fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr); 195*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 19671b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 19771b4ebd2SPeter Brune PetscFunctionReturn(0); 19871b4ebd2SPeter Brune } 19971b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 200c427506bSPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 201c427506bSPeter Brune if (domainerror) { 202*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 20371b4ebd2SPeter Brune PetscFunctionReturn(0); 20471b4ebd2SPeter Brune } 20571b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 20671b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 20771b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 20871b4ebd2SPeter Brune if (monitor) { 20971b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 210*6188f407SPeter Brune if (bt->order == PETSCLINESEARCH_BT_CUBIC) { 21171b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 21271b4ebd2SPeter Brune } else { 21371b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 21471b4ebd2SPeter Brune } 21571b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 21671b4ebd2SPeter Brune } 21771b4ebd2SPeter Brune break; 21871b4ebd2SPeter Brune } else { 21971b4ebd2SPeter Brune if (monitor) { 22071b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 221*6188f407SPeter Brune if (bt->order == PETSCLINESEARCH_BT_CUBIC) { 22271b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 22371b4ebd2SPeter Brune } else { 22471b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 22571b4ebd2SPeter Brune } 22671b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 22771b4ebd2SPeter Brune } 22871b4ebd2SPeter Brune } 22971b4ebd2SPeter Brune } 23071b4ebd2SPeter Brune } 23171b4ebd2SPeter Brune } 23271b4ebd2SPeter Brune 23371b4ebd2SPeter Brune /* postcheck */ 234*6188f407SPeter Brune ierr = PetscLineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 23571b4ebd2SPeter Brune if (changed_y) { 23671b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 23771b4ebd2SPeter Brune } 238c427506bSPeter Brune if (changed_y || changed_w) { /* recompute the function if the step has changed */ 239c427506bSPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 24071b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 24171b4ebd2SPeter Brune if (domainerror) { 242*6188f407SPeter Brune ierr = PetscLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 24371b4ebd2SPeter Brune PetscFunctionReturn(0); 24471b4ebd2SPeter Brune } 245c427506bSPeter Brune ierr = VecNormBegin(G,NORM_2,&gnorm);CHKERRQ(ierr); 246c427506bSPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 247c427506bSPeter Brune ierr = VecNormBegin(Y,NORM_2,&ynorm);CHKERRQ(ierr); 248c427506bSPeter Brune ierr = VecNormEnd(G,NORM_2,&gnorm);CHKERRQ(ierr); 249c427506bSPeter Brune ierr = VecNormEnd(Y,NORM_2,&ynorm);CHKERRQ(ierr); 250c427506bSPeter Brune } 25171b4ebd2SPeter Brune 25271b4ebd2SPeter Brune /* copy the solution over */ 25371b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 254c427506bSPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 255c427506bSPeter Brune ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); 256*6188f407SPeter Brune ierr = PetscLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 257*6188f407SPeter Brune ierr = PetscLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); 25871b4ebd2SPeter Brune PetscFunctionReturn(0); 25971b4ebd2SPeter Brune } 26071b4ebd2SPeter Brune 26171b4ebd2SPeter Brune #undef __FUNCT__ 262*6188f407SPeter Brune #define __FUNCT__ "PetscLineSearchDestroy_BT" 263*6188f407SPeter Brune PetscErrorCode PetscLineSearchDestroy_BT(PetscLineSearch linesearch) { 26471b4ebd2SPeter Brune 26571b4ebd2SPeter Brune PetscFunctionBegin; 26671b4ebd2SPeter Brune PetscFunctionReturn(0); 26771b4ebd2SPeter Brune 26871b4ebd2SPeter Brune } 26971b4ebd2SPeter Brune 27071b4ebd2SPeter Brune 27171b4ebd2SPeter Brune #undef __FUNCT__ 272*6188f407SPeter Brune #define __FUNCT__ "PetscLineSearchSetFromOptions_BT" 273*6188f407SPeter Brune static PetscErrorCode PetscLineSearchSetFromOptions_BT(PetscLineSearch linesearch) 27471b4ebd2SPeter Brune { 27571b4ebd2SPeter Brune 27671b4ebd2SPeter Brune PetscErrorCode ierr; 277*6188f407SPeter Brune PetscLineSearch_BT *bt; 27871b4ebd2SPeter Brune const char *orders[] = {"quadratic", "cubic"}; 27971b4ebd2SPeter Brune PetscInt indx = 0; 28071b4ebd2SPeter Brune PetscBool flg; 28171b4ebd2SPeter Brune PetscFunctionBegin; 28271b4ebd2SPeter Brune 283*6188f407SPeter Brune bt = (PetscLineSearch_BT*)linesearch->data; 28471b4ebd2SPeter Brune 285*6188f407SPeter Brune ierr = PetscOptionsHead("PetscLineSearch BT options");CHKERRQ(ierr); 286*6188f407SPeter Brune ierr = PetscOptionsReal("-linesearch_bt_alpha", "Descent tolerance", "PetscLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 287*6188f407SPeter Brune ierr = PetscOptionsEList("-linesearch_bt_order", "Order of approximation", "PetscLineSearchBT", orders,2,"quadratic",&indx,&flg);CHKERRQ(ierr); 28871b4ebd2SPeter Brune if (flg) { 28971b4ebd2SPeter Brune switch (indx) { 290*6188f407SPeter Brune case 0: bt->order = PETSCLINESEARCH_BT_QUADRATIC; 29171b4ebd2SPeter Brune break; 292*6188f407SPeter Brune case 1: bt->order = PETSCLINESEARCH_BT_CUBIC; 29371b4ebd2SPeter Brune break; 29471b4ebd2SPeter Brune } 29571b4ebd2SPeter Brune } 29671b4ebd2SPeter Brune 29771b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 29871b4ebd2SPeter Brune PetscFunctionReturn(0); 29971b4ebd2SPeter Brune } 30071b4ebd2SPeter Brune 30171b4ebd2SPeter Brune 30271b4ebd2SPeter Brune EXTERN_C_BEGIN 30371b4ebd2SPeter Brune #undef __FUNCT__ 304*6188f407SPeter Brune #define __FUNCT__ "PetscLineSearchCreate_BT" 305*6188f407SPeter Brune PetscErrorCode PetscLineSearchCreate_BT(PetscLineSearch linesearch) 30671b4ebd2SPeter Brune { 30771b4ebd2SPeter Brune 308*6188f407SPeter Brune PetscLineSearch_BT *bt; 30971b4ebd2SPeter Brune PetscErrorCode ierr; 31071b4ebd2SPeter Brune 31171b4ebd2SPeter Brune PetscFunctionBegin; 312*6188f407SPeter Brune linesearch->ops->apply = PetscLineSearchApply_BT; 313*6188f407SPeter Brune linesearch->ops->destroy = PetscLineSearchDestroy_BT; 314*6188f407SPeter Brune linesearch->ops->setfromoptions = PetscLineSearchSetFromOptions_BT; 31571b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 31671b4ebd2SPeter Brune linesearch->ops->view = PETSC_NULL; 31771b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 31871b4ebd2SPeter Brune 319*6188f407SPeter Brune ierr = PetscNewLog(linesearch, PetscLineSearch_BT, &bt);CHKERRQ(ierr); 32071b4ebd2SPeter Brune linesearch->data = (void *)bt; 32171b4ebd2SPeter Brune linesearch->max_its = 40; 322*6188f407SPeter Brune bt->order = PETSCLINESEARCH_BT_CUBIC; 32371b4ebd2SPeter Brune bt->alpha = 1e-4; 32471b4ebd2SPeter Brune 32571b4ebd2SPeter Brune PetscFunctionReturn(0); 32671b4ebd2SPeter Brune } 32771b4ebd2SPeter Brune EXTERN_C_END 328