1*71b4ebd2SPeter Brune #include <private/linesearchimpl.h> 2*71b4ebd2SPeter Brune #include <private/snesimpl.h> 3*71b4ebd2SPeter Brune 4*71b4ebd2SPeter Brune typedef enum {LINESEARCH_BT_QUADRATIC, LINESEARCH_BT_CUBIC} LineSearchBTOrder; 5*71b4ebd2SPeter Brune 6*71b4ebd2SPeter Brune typedef struct { 7*71b4ebd2SPeter Brune PetscReal alpha; /* sufficient decrease parameter */ 8*71b4ebd2SPeter Brune LineSearchBTOrder order; 9*71b4ebd2SPeter Brune } LineSearch_BT; 10*71b4ebd2SPeter Brune 11*71b4ebd2SPeter Brune /*MC 12*71b4ebd2SPeter Brune LineSearchBT - Backtracking line searches. 13*71b4ebd2SPeter Brune 14*71b4ebd2SPeter Brune These linesearches try a polynomial fit for the L2 norm of the error 15*71b4ebd2SPeter Brune using the gradient. Failing that, they step back and try again. 16*71b4ebd2SPeter Brune 17*71b4ebd2SPeter Brune Level: advanced 18*71b4ebd2SPeter Brune 19*71b4ebd2SPeter Brune .keywords: SNES, LineSearch, damping 20*71b4ebd2SPeter Brune 21*71b4ebd2SPeter Brune .seealso: LineSearchCreate(), LineSearchSetType() 22*71b4ebd2SPeter Brune M*/ 23*71b4ebd2SPeter Brune 24*71b4ebd2SPeter Brune #undef __FUNCT__ 25*71b4ebd2SPeter Brune #define __FUNCT__ "LineSearchApply_BT" 26*71b4ebd2SPeter Brune 27*71b4ebd2SPeter Brune PetscErrorCode LineSearchApply_BT(LineSearch linesearch) 28*71b4ebd2SPeter Brune { 29*71b4ebd2SPeter Brune PetscBool changed_y, changed_w; 30*71b4ebd2SPeter Brune PetscErrorCode ierr; 31*71b4ebd2SPeter Brune Vec X, F, Y, W, G; 32*71b4ebd2SPeter Brune SNES snes; 33*71b4ebd2SPeter Brune PetscReal fnorm, xnorm, ynorm, gnorm, gnormprev; 34*71b4ebd2SPeter Brune PetscReal lambda, lambdatemp, lambdaprev, minlambda, maxstep, rellength, initslope, alpha; 35*71b4ebd2SPeter Brune PetscReal t1, t2, a, b, d, steptol; 36*71b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 37*71b4ebd2SPeter Brune PetscScalar cinitslope; 38*71b4ebd2SPeter Brune #endif 39*71b4ebd2SPeter Brune PetscBool domainerror; 40*71b4ebd2SPeter Brune PetscViewer monitor; 41*71b4ebd2SPeter Brune PetscInt max_its, count; 42*71b4ebd2SPeter Brune LineSearch_BT *bt; 43*71b4ebd2SPeter Brune Mat jac; 44*71b4ebd2SPeter Brune 45*71b4ebd2SPeter Brune 46*71b4ebd2SPeter Brune PetscFunctionBegin; 47*71b4ebd2SPeter Brune 48*71b4ebd2SPeter Brune ierr = LineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr); 49*71b4ebd2SPeter Brune ierr = LineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 50*71b4ebd2SPeter Brune ierr = LineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); 51*71b4ebd2SPeter Brune ierr = LineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 52*71b4ebd2SPeter Brune ierr = LineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); 53*71b4ebd2SPeter Brune ierr = LineSearchGetTolerances(linesearch, &steptol, &maxstep, PETSC_NULL, PETSC_NULL, PETSC_NULL, &max_its); 54*71b4ebd2SPeter Brune bt = (LineSearch_BT *)linesearch->data; 55*71b4ebd2SPeter Brune 56*71b4ebd2SPeter Brune alpha = bt->alpha; 57*71b4ebd2SPeter Brune 58*71b4ebd2SPeter Brune ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 59*71b4ebd2SPeter Brune if (!jac) { 60*71b4ebd2SPeter Brune SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "LineSearchBT requires a Jacobian matrix"); 61*71b4ebd2SPeter Brune } 62*71b4ebd2SPeter Brune 63*71b4ebd2SPeter Brune ierr = PetscLogEventBegin(LineSearch_Apply,snes,X,F,G);CHKERRQ(ierr); 64*71b4ebd2SPeter Brune 65*71b4ebd2SPeter Brune /* precheck */ 66*71b4ebd2SPeter Brune ierr = LineSearchPreCheck(linesearch, &changed_y);CHKERRQ(ierr); 67*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); 68*71b4ebd2SPeter Brune 69*71b4ebd2SPeter Brune ierr = VecNorm(Y, NORM_2, &ynorm);CHKERRQ(ierr); 70*71b4ebd2SPeter Brune if (ynorm == 0.0) { 71*71b4ebd2SPeter Brune if (monitor) { 72*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 73*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 74*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 75*71b4ebd2SPeter Brune } 76*71b4ebd2SPeter Brune ierr = VecCopy(X,W);CHKERRQ(ierr); 77*71b4ebd2SPeter Brune ierr = VecCopy(F,G);CHKERRQ(ierr); 78*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 79*71b4ebd2SPeter Brune PetscFunctionReturn(0); 80*71b4ebd2SPeter Brune } 81*71b4ebd2SPeter Brune if (ynorm > maxstep) { /* Step too big, so scale back */ 82*71b4ebd2SPeter Brune if (monitor) { 83*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 84*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);CHKERRQ(ierr); 85*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 86*71b4ebd2SPeter Brune } 87*71b4ebd2SPeter Brune ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); 88*71b4ebd2SPeter Brune ynorm = maxstep; 89*71b4ebd2SPeter Brune } 90*71b4ebd2SPeter Brune ierr = VecMaxPointwiseDivide(Y,X,&rellength);CHKERRQ(ierr); 91*71b4ebd2SPeter Brune minlambda = steptol/rellength; 92*71b4ebd2SPeter Brune ierr = MatMult(jac,Y,W);CHKERRQ(ierr); 93*71b4ebd2SPeter Brune #if defined(PETSC_USE_COMPLEX) 94*71b4ebd2SPeter Brune ierr = VecDot(F,W,&cinitslope);CHKERRQ(ierr); 95*71b4ebd2SPeter Brune initslope = PetscRealPart(cinitslope); 96*71b4ebd2SPeter Brune #else 97*71b4ebd2SPeter Brune ierr = VecDot(F,W,&initslope);CHKERRQ(ierr); 98*71b4ebd2SPeter Brune #endif 99*71b4ebd2SPeter Brune if (initslope > 0.0) initslope = -initslope; 100*71b4ebd2SPeter Brune if (initslope == 0.0) initslope = -1.0; 101*71b4ebd2SPeter Brune 102*71b4ebd2SPeter Brune ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr); 103*71b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 104*71b4ebd2SPeter Brune ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 105*71b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 106*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 107*71b4ebd2SPeter Brune PetscFunctionReturn(0); 108*71b4ebd2SPeter Brune } 109*71b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 110*71b4ebd2SPeter Brune if (snes->domainerror) { 111*71b4ebd2SPeter Brune ierr = PetscLogEventEnd(SNES_LineSearch,snes,X,F,G);CHKERRQ(ierr); 112*71b4ebd2SPeter Brune PetscFunctionReturn(0); 113*71b4ebd2SPeter Brune } 114*71b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 115*71b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 116*71b4ebd2SPeter Brune ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 117*71b4ebd2SPeter Brune if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* Sufficient reduction */ 118*71b4ebd2SPeter Brune if (monitor) { 119*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 120*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);CHKERRQ(ierr); 121*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 122*71b4ebd2SPeter Brune } 123*71b4ebd2SPeter Brune } else { 124*71b4ebd2SPeter Brune /* Fit points with quadratic */ 125*71b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope); 126*71b4ebd2SPeter Brune lambdaprev = lambda; 127*71b4ebd2SPeter Brune gnormprev = gnorm; 128*71b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 129*71b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 130*71b4ebd2SPeter Brune else lambda = lambdatemp; 131*71b4ebd2SPeter Brune 132*71b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 133*71b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 134*71b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 135*71b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 136*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 137*71b4ebd2SPeter Brune PetscFunctionReturn(0); 138*71b4ebd2SPeter Brune } 139*71b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 140*71b4ebd2SPeter Brune if (snes->domainerror) { 141*71b4ebd2SPeter Brune ierr = PetscLogEventEnd(SNES_LineSearch,snes,X,F,G);CHKERRQ(ierr); 142*71b4ebd2SPeter Brune PetscFunctionReturn(0); 143*71b4ebd2SPeter Brune } 144*71b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 145*71b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 146*71b4ebd2SPeter Brune if (monitor) { 147*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 148*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",gnorm);CHKERRQ(ierr); 149*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 150*71b4ebd2SPeter Brune } 151*71b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 152*71b4ebd2SPeter Brune if (monitor) { 153*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 154*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); 155*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 156*71b4ebd2SPeter Brune } 157*71b4ebd2SPeter Brune } else { 158*71b4ebd2SPeter Brune /* Fit points with cubic */ 159*71b4ebd2SPeter Brune for (count = 0; count < max_its; count++) { 160*71b4ebd2SPeter Brune if (lambda <= minlambda) { 161*71b4ebd2SPeter Brune if (monitor) { 162*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 163*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 164*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor, 165*71b4ebd2SPeter Brune " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 166*71b4ebd2SPeter Brune fnorm, gnorm, ynorm, minlambda, lambda, initslope);CHKERRQ(ierr); 167*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 168*71b4ebd2SPeter Brune } 169*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 170*71b4ebd2SPeter Brune PetscFunctionReturn(0); 171*71b4ebd2SPeter Brune } 172*71b4ebd2SPeter Brune if (bt->order == LINESEARCH_BT_CUBIC) { 173*71b4ebd2SPeter Brune t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope; 174*71b4ebd2SPeter Brune t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 175*71b4ebd2SPeter Brune a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 176*71b4ebd2SPeter Brune b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 177*71b4ebd2SPeter Brune d = b*b - 3*a*initslope; 178*71b4ebd2SPeter Brune if (d < 0.0) d = 0.0; 179*71b4ebd2SPeter Brune if (a == 0.0) { 180*71b4ebd2SPeter Brune lambdatemp = -initslope/(2.0*b); 181*71b4ebd2SPeter Brune } else { 182*71b4ebd2SPeter Brune lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); 183*71b4ebd2SPeter Brune } 184*71b4ebd2SPeter Brune } else if (bt->order == LINESEARCH_BT_QUADRATIC) { 185*71b4ebd2SPeter Brune lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope); 186*71b4ebd2SPeter Brune } 187*71b4ebd2SPeter Brune lambdaprev = lambda; 188*71b4ebd2SPeter Brune gnormprev = gnorm; 189*71b4ebd2SPeter Brune if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 190*71b4ebd2SPeter Brune if (lambdatemp <= .1*lambda) lambda = .1*lambda; 191*71b4ebd2SPeter Brune else lambda = lambdatemp; 192*71b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 193*71b4ebd2SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 194*71b4ebd2SPeter Brune ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 195*71b4ebd2SPeter Brune ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", 196*71b4ebd2SPeter Brune fnorm,gnorm,ynorm,lambda,initslope);CHKERRQ(ierr); 197*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); 198*71b4ebd2SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 199*71b4ebd2SPeter Brune PetscFunctionReturn(0); 200*71b4ebd2SPeter Brune } 201*71b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); 202*71b4ebd2SPeter Brune if (snes->domainerror) { 203*71b4ebd2SPeter Brune ierr = PetscLogEventEnd(SNES_LineSearch,snes,X,F,G);CHKERRQ(ierr); 204*71b4ebd2SPeter Brune PetscFunctionReturn(0); 205*71b4ebd2SPeter Brune } 206*71b4ebd2SPeter Brune ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); 207*71b4ebd2SPeter Brune if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 208*71b4ebd2SPeter Brune if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 209*71b4ebd2SPeter Brune if (monitor) { 210*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 211*71b4ebd2SPeter Brune if (bt->order == LINESEARCH_BT_CUBIC) { 212*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 213*71b4ebd2SPeter Brune } else { 214*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 215*71b4ebd2SPeter Brune } 216*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 217*71b4ebd2SPeter Brune } 218*71b4ebd2SPeter Brune break; 219*71b4ebd2SPeter Brune } else { 220*71b4ebd2SPeter Brune if (monitor) { 221*71b4ebd2SPeter Brune ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 222*71b4ebd2SPeter Brune if (bt->order == LINESEARCH_BT_CUBIC) { 223*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 224*71b4ebd2SPeter Brune } else { 225*71b4ebd2SPeter Brune ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);CHKERRQ(ierr); 226*71b4ebd2SPeter Brune } 227*71b4ebd2SPeter Brune ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 228*71b4ebd2SPeter Brune } 229*71b4ebd2SPeter Brune } 230*71b4ebd2SPeter Brune } 231*71b4ebd2SPeter Brune } 232*71b4ebd2SPeter Brune } 233*71b4ebd2SPeter Brune 234*71b4ebd2SPeter Brune /* postcheck */ 235*71b4ebd2SPeter Brune ierr = LineSearchPostCheck(linesearch, &changed_y, &changed_w);CHKERRQ(ierr); 236*71b4ebd2SPeter Brune if (changed_y) { 237*71b4ebd2SPeter Brune ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); 238*71b4ebd2SPeter Brune } 239*71b4ebd2SPeter Brune ierr = SNESComputeFunction(snes,W,F);CHKERRQ(ierr); 240*71b4ebd2SPeter Brune ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); 241*71b4ebd2SPeter Brune if (domainerror) { 242*71b4ebd2SPeter Brune ierr = LineSearchSetSuccess(linesearch, PETSC_FALSE); 243*71b4ebd2SPeter Brune PetscFunctionReturn(0); 244*71b4ebd2SPeter Brune } 245*71b4ebd2SPeter Brune 246*71b4ebd2SPeter Brune /* copy the solution over */ 247*71b4ebd2SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 248*71b4ebd2SPeter Brune ierr = LineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); 249*71b4ebd2SPeter Brune ierr = LineSearchComputeNorms(linesearch);CHKERRQ(ierr); 250*71b4ebd2SPeter Brune PetscFunctionReturn(0); 251*71b4ebd2SPeter Brune } 252*71b4ebd2SPeter Brune 253*71b4ebd2SPeter Brune #undef __FUNCT__ 254*71b4ebd2SPeter Brune #define __FUNCT__ "LineSearchDestroy_BT" 255*71b4ebd2SPeter Brune PetscErrorCode LineSearchDestroy_BT(LineSearch linesearch) { 256*71b4ebd2SPeter Brune 257*71b4ebd2SPeter Brune PetscFunctionBegin; 258*71b4ebd2SPeter Brune PetscFunctionReturn(0); 259*71b4ebd2SPeter Brune 260*71b4ebd2SPeter Brune } 261*71b4ebd2SPeter Brune 262*71b4ebd2SPeter Brune 263*71b4ebd2SPeter Brune #undef __FUNCT__ 264*71b4ebd2SPeter Brune #define __FUNCT__ "LineSearchSetFromOptions_BT" 265*71b4ebd2SPeter Brune static PetscErrorCode LineSearchSetFromOptions_BT(LineSearch linesearch) 266*71b4ebd2SPeter Brune { 267*71b4ebd2SPeter Brune 268*71b4ebd2SPeter Brune PetscErrorCode ierr; 269*71b4ebd2SPeter Brune LineSearch_BT *bt; 270*71b4ebd2SPeter Brune const char *orders[] = {"quadratic", "cubic"}; 271*71b4ebd2SPeter Brune PetscInt indx = 0; 272*71b4ebd2SPeter Brune PetscBool flg; 273*71b4ebd2SPeter Brune PetscFunctionBegin; 274*71b4ebd2SPeter Brune 275*71b4ebd2SPeter Brune bt = (LineSearch_BT*)linesearch->data; 276*71b4ebd2SPeter Brune 277*71b4ebd2SPeter Brune ierr = PetscOptionsHead("LineSearch BT options");CHKERRQ(ierr); 278*71b4ebd2SPeter Brune ierr = PetscOptionsReal("-linesearch_bt_alpha", "Descent tolerance", "LineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr); 279*71b4ebd2SPeter Brune ierr = PetscOptionsEList("-linesearch_bt_order", "Order of approximation", "LineSearchBT", orders,2,"quadratic",&indx,&flg);CHKERRQ(ierr); 280*71b4ebd2SPeter Brune if (flg) { 281*71b4ebd2SPeter Brune switch (indx) { 282*71b4ebd2SPeter Brune case 0: bt->order = LINESEARCH_BT_QUADRATIC; 283*71b4ebd2SPeter Brune break; 284*71b4ebd2SPeter Brune case 1: bt->order = LINESEARCH_BT_CUBIC; 285*71b4ebd2SPeter Brune break; 286*71b4ebd2SPeter Brune } 287*71b4ebd2SPeter Brune } 288*71b4ebd2SPeter Brune 289*71b4ebd2SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 290*71b4ebd2SPeter Brune PetscFunctionReturn(0); 291*71b4ebd2SPeter Brune } 292*71b4ebd2SPeter Brune 293*71b4ebd2SPeter Brune 294*71b4ebd2SPeter Brune EXTERN_C_BEGIN 295*71b4ebd2SPeter Brune #undef __FUNCT__ 296*71b4ebd2SPeter Brune #define __FUNCT__ "LineSearchCreate_BT" 297*71b4ebd2SPeter Brune PetscErrorCode LineSearchCreate_BT(LineSearch linesearch) 298*71b4ebd2SPeter Brune { 299*71b4ebd2SPeter Brune 300*71b4ebd2SPeter Brune LineSearch_BT *bt; 301*71b4ebd2SPeter Brune PetscErrorCode ierr; 302*71b4ebd2SPeter Brune 303*71b4ebd2SPeter Brune PetscFunctionBegin; 304*71b4ebd2SPeter Brune linesearch->ops->apply = LineSearchApply_BT; 305*71b4ebd2SPeter Brune linesearch->ops->destroy = LineSearchDestroy_BT; 306*71b4ebd2SPeter Brune linesearch->ops->setfromoptions = LineSearchSetFromOptions_BT; 307*71b4ebd2SPeter Brune linesearch->ops->reset = PETSC_NULL; 308*71b4ebd2SPeter Brune linesearch->ops->view = PETSC_NULL; 309*71b4ebd2SPeter Brune linesearch->ops->setup = PETSC_NULL; 310*71b4ebd2SPeter Brune 311*71b4ebd2SPeter Brune ierr = PetscNewLog(linesearch, LineSearch_BT, &bt);CHKERRQ(ierr); 312*71b4ebd2SPeter Brune linesearch->data = (void *)bt; 313*71b4ebd2SPeter Brune linesearch->max_its = 40; 314*71b4ebd2SPeter Brune bt->order = LINESEARCH_BT_CUBIC; 315*71b4ebd2SPeter Brune bt->alpha = 1e-4; 316*71b4ebd2SPeter Brune 317*71b4ebd2SPeter Brune PetscFunctionReturn(0); 318*71b4ebd2SPeter Brune } 319*71b4ebd2SPeter Brune EXTERN_C_END 320