xref: /petsc/src/snes/linesearch/impls/bt/linesearchbt.c (revision 0e9bae810fdaeb60e2713eaa8ddb89f42e079fd1)
1 #include <petsc-private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
2 #include <petsc-private/snesimpl.h>
3 
4 typedef struct {
5   PetscReal        alpha; /* sufficient decrease parameter */
6 } SNESLineSearch_BT;
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "SNESLineSearchBTSetAlpha"
10 /*@
11    SNESLineSearchBTSetAlpha - Sets the descent parameter, alpha, in the BT linesearch variant.
12 
13    Input Parameters:
14 +  linesearch - linesearch context
15 -  alpha - The descent parameter
16 
17    Level: intermediate
18 
19 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
20 @*/
21 PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha)
22 {
23   SNESLineSearch_BT  *bt;
24   PetscFunctionBegin;
25   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
26   bt = (SNESLineSearch_BT *)linesearch->data;
27   bt->alpha = alpha;
28   PetscFunctionReturn(0);
29 }
30 
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "SNESLineSearchBTGetAlpha"
34 /*@
35    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant.
36 
37    Input Parameters:
38 .  linesearch - linesearch context
39 
40    Output Parameters:
41 .  alpha - The descent parameter
42 
43    Level: intermediate
44 
45 .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
46 @*/
47 PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
48 {
49   SNESLineSearch_BT  *bt;
50   PetscFunctionBegin;
51   PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1);
52   bt = (SNESLineSearch_BT *)linesearch->data;
53   *alpha = bt->alpha;
54   PetscFunctionReturn(0);
55 }
56 
57 #undef __FUNCT__
58 #define __FUNCT__ "SNESLineSearchApply_BT"
59 static PetscErrorCode  SNESLineSearchApply_BT(SNESLineSearch linesearch)
60 {
61   PetscBool      changed_y,changed_w;
62   PetscErrorCode ierr;
63   Vec            X,F,Y,W,G;
64   SNES           snes;
65   PetscReal      fnorm, xnorm, ynorm, gnorm, gnormprev;
66   PetscReal      lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol;
67   PetscReal      t1,t2,a,b,d;
68 #if defined(PETSC_USE_COMPLEX)
69   PetscScalar    cinitslope;
70 #endif
71   PetscBool      domainerror;
72   PetscViewer    monitor;
73   PetscInt       max_its,count;
74   SNESLineSearch_BT  *bt;
75   Mat            jac;
76 
77 
78   PetscFunctionBegin;
79 
80   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
81   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
82   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
83   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
84   ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr);
85   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,PETSC_NULL,PETSC_NULL,PETSC_NULL,&max_its);
86   ierr = SNESGetTolerances(snes,PETSC_NULL,PETSC_NULL,&stol,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
87   bt = (SNESLineSearch_BT *)linesearch->data;
88 
89   alpha = bt->alpha;
90 
91   ierr = SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
92   if (!jac) SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
93 
94   /* precheck */
95   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
96   ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
97 
98   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
99   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
100   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
101   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
102 
103   if (ynorm == 0.0) {
104     if (monitor) {
105       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
106       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
107       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
108     }
109     ierr   = VecCopy(X,W);CHKERRQ(ierr);
110     ierr   = VecCopy(F,G);CHKERRQ(ierr);
111     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
112     PetscFunctionReturn(0);
113   }
114   if (ynorm > maxstep) {	/* Step too big, so scale back */
115     if (monitor) {
116       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
117       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr);
118       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
119     }
120     ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr);
121     ynorm = maxstep;
122   }
123   ierr      = MatMult(jac,Y,W);CHKERRQ(ierr);
124 #if defined(PETSC_USE_COMPLEX)
125   ierr      = VecDot(F,W,&cinitslope);CHKERRQ(ierr);
126   initslope = PetscRealPart(cinitslope);
127 #else
128   ierr      = VecDot(F,W,&initslope);CHKERRQ(ierr);
129 #endif
130   if (initslope > 0.0)  initslope = -initslope;
131   if (initslope == 0.0) initslope = -1.0;
132 
133   ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
134   if (linesearch->ops->viproject) {
135     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
136   }
137   if (snes->nfuncs >= snes->max_funcs) {
138     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
139     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
140     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
141     PetscFunctionReturn(0);
142   }
143   ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
144   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
145   if (domainerror) {
146     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
147     PetscFunctionReturn(0);
148   }
149   if (linesearch->ops->vinorm) {
150     gnorm = fnorm;
151     ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
152   } else {
153     ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
154   }
155 
156   if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
157   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
158   if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope || ynorm < stol*xnorm) { /* Sufficient reduction or step tolerance convergence */
159     if (monitor) {
160       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
161       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr);
162       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
163     }
164   } else {
165     /* Fit points with quadratic */
166     lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope);
167     lambdaprev = lambda;
168     gnormprev  = gnorm;
169     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
170     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
171     else                         lambda = lambdatemp;
172 
173     ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
174     if (linesearch->ops->viproject) {
175       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
176     }
177     if (snes->nfuncs >= snes->max_funcs) {
178       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
179       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
180       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
181       PetscFunctionReturn(0);
182     }
183     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
184     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
185     if (domainerror) {
186       PetscFunctionReturn(0);
187     }
188     if (linesearch->ops->vinorm) {
189       gnorm = fnorm;
190       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
191     } else {
192       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
193     }
194     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
195     if (monitor) {
196       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
197       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr);
198       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
199     }
200     if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
201       if (monitor) {
202         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
203         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr);
204         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
205       }
206     } else {
207       /* Fit points with cubic */
208       for (count = 0; count < max_its; count++) {
209         if (lambda <= minlambda) {
210           if (monitor) {
211             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
212             ierr = PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
213             ierr = PetscViewerASCIIPrintf(monitor,
214                                           "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
215                                           (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr);
216             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
217           }
218           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
219           PetscFunctionReturn(0);
220         }
221         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
222           t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope;
223           t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
224           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
225           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
226           d  = b*b - 3*a*initslope;
227           if (d < 0.0) d = 0.0;
228           if (a == 0.0) {
229             lambdatemp = -initslope/(2.0*b);
230           } else {
231             lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
232           }
233         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
234           lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope);
235         } else {
236           SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt");
237         }
238           lambdaprev = lambda;
239           gnormprev  = gnorm;
240         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
241         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
242         else                         lambda     = lambdatemp;
243         ierr  = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
244         if (linesearch->ops->viproject) {
245           ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr);
246         }
247         if (snes->nfuncs >= snes->max_funcs) {
248           ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
249           ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
250                             (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
251           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
252           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
253           PetscFunctionReturn(0);
254         }
255         ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
256         ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
257         if (domainerror) {
258           ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
259           PetscFunctionReturn(0);
260         }
261         if (linesearch->ops->vinorm) {
262           gnorm = fnorm;
263           ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
264         } else {
265           ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
266         }
267         if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
268         if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
269           if (monitor) {
270             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
271             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
272               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
273             } else {
274               ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr);
275             }
276             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
277           }
278           break;
279         } else {
280           if (monitor) {
281             ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
282             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
283               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);
284             } else {
285               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);
286             }
287             ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
288           }
289         }
290       }
291     }
292   }
293 
294   /* postcheck */
295   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
296   if (changed_y) {
297     ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr);
298     if (linesearch->ops->viproject) {
299       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
300     }
301   }
302   if (changed_y || changed_w) { /* recompute the function if the step has changed */
303     ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr);
304     ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
305     if (domainerror) {
306       ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
307       PetscFunctionReturn(0);
308     }
309     if (linesearch->ops->vinorm) {
310       gnorm = fnorm;
311       ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr);
312     } else {
313       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);
314     }
315     ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr);
316     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
317 
318   }
319 
320   /* copy the solution over */
321   ierr = VecCopy(W, X);CHKERRQ(ierr);
322   ierr = VecCopy(G, F);CHKERRQ(ierr);
323   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
324   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
325   ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr);
326   PetscFunctionReturn(0);
327 }
328 
329 #undef __FUNCT__
330 #define __FUNCT__ "SNESLineSearchView_BT"
331 PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
332 {
333   PetscErrorCode    ierr;
334   PetscBool         iascii;
335   SNESLineSearch_BT *bt;
336   PetscFunctionBegin;
337   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
338   bt = (SNESLineSearch_BT*)linesearch->data;
339   if (iascii) {
340     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
341     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");CHKERRQ(ierr);
342     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
343     ierr = PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");CHKERRQ(ierr);
344     }
345     ierr = PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", bt->alpha);CHKERRQ(ierr);
346   }
347   PetscFunctionReturn(0);
348 }
349 
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "SNESLineSearchDestroy_BT"
353 static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
354 {
355   PetscErrorCode ierr;
356 
357   PetscFunctionBegin;
358   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
359   PetscFunctionReturn(0);
360 }
361 
362 
363 #undef __FUNCT__
364 #define __FUNCT__ "SNESLineSearchSetFromOptions_BT"
365 static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch)
366 {
367 
368   PetscErrorCode       ierr;
369   SNESLineSearch_BT    *bt;
370   PetscFunctionBegin;
371 
372   bt = (SNESLineSearch_BT*)linesearch->data;
373 
374   ierr = PetscOptionsHead("SNESLineSearch BT options");CHKERRQ(ierr);
375   ierr = PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);CHKERRQ(ierr);
376 
377   ierr = PetscOptionsTail();CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 
382 #undef __FUNCT__
383 #define __FUNCT__ "SNESLineSearchCreate_BT"
384 /*MC
385    SNESLINESEARCHBT - Backtracking line search.
386 
387    This line search finds the minimum of a polynomial fitting of the L2 norm of the
388    function. If this fit does not satisfy the conditions for progress, the interval shrinks
389    and the fit is reattempted at most max_it times or until lambda is below minlambda.
390 
391    Options Database Keys:
392 +  -snes_linesearch_alpha<1e-4> - slope descent parameter
393 .  -snes_linesearch_damping<1.0> - initial step length
394 .  -snes_linesearch_max_it<40> - maximum number of shrinking step
395 .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
396 -  -snes_linesearch_order<cubic,quadratic> - order of the approximation
397 
398    Level: advanced
399 
400    Notes:
401    This line search is taken from "Numerical Methods for Unconstrained
402    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
403 
404 .keywords: SNES, SNESLineSearch, damping
405 
406 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
407 M*/
408 PETSC_EXTERN_C PetscErrorCode SNESLineSearchCreate_BT(SNESLineSearch linesearch)
409 {
410 
411   SNESLineSearch_BT  *bt;
412   PetscErrorCode ierr;
413 
414   PetscFunctionBegin;
415   linesearch->ops->apply          = SNESLineSearchApply_BT;
416   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
417   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
418   linesearch->ops->reset          = PETSC_NULL;
419   linesearch->ops->view           = SNESLineSearchView_BT;
420   linesearch->ops->setup          = PETSC_NULL;
421 
422   ierr = PetscNewLog(linesearch, SNESLineSearch_BT, &bt);CHKERRQ(ierr);
423   linesearch->data = (void *)bt;
424   linesearch->max_its = 40;
425   linesearch->order = SNES_LINESEARCH_ORDER_CUBIC;
426   bt->alpha = 1e-4;
427 
428   PetscFunctionReturn(0);
429 }
430